Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PyrGeom.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PyrGeom.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Pyramidic geometry information.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <SpatialDomains/PyrGeom.h>
37 
40 #include <StdRegions/StdPyrExp.h>
41 #include <SpatialDomains/SegGeom.h>
42 #include <SpatialDomains/TriGeom.h>
44 
45 
46 namespace Nektar
47 {
48  namespace SpatialDomains
49  {
51  {
53  }
54 
56  Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
57  {
59 
60  /// Copy the face shared pointers
61  m_faces.insert(m_faces.begin(), faces, faces+PyrGeom::kNfaces);
62 
63  /// Set up orientation vectors with correct amount of elements.
64  m_eorient.resize(kNedges);
65  m_forient.resize(kNfaces);
66 
71  SetUpXmap();
72  SetUpCoeffs(m_xmap->GetNcoeffs());
73  }
74 
76  {
77 
78  }
79 
81  {
83  {
84  int i;
85  GeomType Gtype = eRegular;
86 
87  v_FillGeom();
88 
89  // check to see if expansions are linear
90  for(i = 0; i < m_coordim; ++i)
91  {
92  if (m_xmap->GetBasisNumModes(0) != 2 ||
93  m_xmap->GetBasisNumModes(1) != 2 ||
94  m_xmap->GetBasisNumModes(2) != 2 )
95  {
96  Gtype = eDeformed;
97  }
98  }
99 
100  // check to see if all quadrilateral faces are parallelograms
101  if(Gtype == eRegular)
102  {
103  // Ensure each face is a parallelogram? Check this.
104  for (i = 0; i < m_coordim; i++)
105  {
106  if( fabs( (*m_verts[0])(i) -
107  (*m_verts[1])(i) +
108  (*m_verts[2])(i) -
109  (*m_verts[3])(i) )
111  {
112  Gtype = eDeformed;
113  break;
114  }
115  }
116  }
117 
119  Gtype, m_coordim, m_xmap, m_coeffs);
121  }
122  }
123 
125  const Array<OneD, const NekDouble> &coords,
126  Array<OneD, NekDouble> &Lcoords)
127  {
128  NekDouble ptdist = 1e6;
129 
130  v_FillGeom();
131 
132  // calculate local coordinate for coord
133  if(GetMetricInfo()->GetGtype() == eRegular)
134  { // Based on Spen's book, page 99
135 
136  // Point inside tetrahedron
137  PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
138 
139  // Edges
140  PointGeom er0, e10, e30, e40;
141  er0.Sub(r,*m_verts[0]);
142  e10.Sub(*m_verts[1],*m_verts[0]);
143  e30.Sub(*m_verts[3],*m_verts[0]);
144  e40.Sub(*m_verts[4],*m_verts[0]);
145 
146 
147  // Cross products (Normal times area)
148  PointGeom cp1030, cp3040, cp4010;
149  cp1030.Mult(e10,e30);
150  cp3040.Mult(e30,e40);
151  cp4010.Mult(e40,e10);
152 
153 
154  // Barycentric coordinates (relative volume)
155  NekDouble V = e40.dot(cp1030); // Pyramid Volume = {(e40)dot(e10)x(e30)}/4
156  NekDouble scaleFactor = 2.0/3.0;
157  NekDouble v1 = er0.dot(cp3040) / V; // volume1 = {(er0)dot(e30)x(e40)}/6
158  NekDouble v2 = er0.dot(cp4010) / V; // volume2 = {(er0)dot(e40)x(e10)}/6
159  NekDouble beta = v1 * scaleFactor;
160  NekDouble gamma = v2 * scaleFactor;
161  NekDouble delta = er0.dot(cp1030) / V; // volume3 = {(er0)dot(e10)x(e30)}/4
162 
163  // Make Pyramid bigger
164  Lcoords[0] = 2.0*beta - 1.0;
165  Lcoords[1] = 2.0*gamma - 1.0;
166  Lcoords[2] = 2.0*delta - 1.0;
167 
168  // Set ptdist to distance to nearest vertex
169  for(int i = 0; i < 5; ++i)
170  {
171  ptdist = min(ptdist,r.dist(*m_verts[i]));
172  }
173 
174  }
175  else
176  {
178  "inverse mapping must be set up to use this call");
179  }
180  return ptdist;
181  }
182 
184  {
185  return 5;
186  }
187 
189  {
190  return 8;
191  }
192 
194  {
195  return 5;
196  }
197 
198  int PyrGeom::v_GetDir(const int faceidx, const int facedir) const
199  {
200  if (faceidx == 0)
201  {
202  return facedir;
203  }
204  else if (faceidx == 1 || faceidx == 3)
205  {
206  return 2 * facedir;
207  }
208  else
209  {
210  return 1 + facedir;
211  }
212  }
213 
215  {
216  // find edge 0
217  int i,j;
218  unsigned int check;
219 
220  SegGeomSharedPtr edge;
221 
222  // First set up the 4 bottom edges
223  int f;
224  for (f = 1; f < 5; f++)
225  {
226  int nEdges = m_faces[f]->GetNumEdges();
227  check = 0;
228  for (i = 0; i < 4; i++)
229  {
230  for (j = 0; j < nEdges; j++)
231  {
232  if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
233  {
234  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[0])->GetEdge(i));
235  m_edges.push_back(edge);
236  check++;
237  }
238  }
239  }
240 
241  if (check < 1)
242  {
243  std::ostringstream errstrm;
244  errstrm << "Connected faces do not share an edge. Faces ";
245  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
246  ASSERTL0(false, errstrm.str());
247  }
248  else if (check > 1)
249  {
250  std::ostringstream errstrm;
251  errstrm << "Connected faces share more than one edge. Faces ";
252  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
253  ASSERTL0(false, errstrm.str());
254  }
255  }
256 
257  // Then, set up the 4 vertical edges
258  check = 0;
259  for(i = 0; i < 3; i++) //Set up the vertical edge :face(1) and face(4)
260  {
261  for(j = 0; j < 3; j++)
262  {
263  if( (m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j) )
264  {
265  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[1])->GetEdge(i));
266  m_edges.push_back(edge);
267  check++;
268  }
269  }
270  }
271  if( check < 1 )
272  {
273  std::ostringstream errstrm;
274  errstrm << "Connected faces do not share an edge. Faces ";
275  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
276  ASSERTL0(false, errstrm.str());
277  }
278  else if( check > 1)
279  {
280  std::ostringstream errstrm;
281  errstrm << "Connected faces share more than one edge. Faces ";
282  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
283  ASSERTL0(false, errstrm.str());
284  }
285 
286  // Set up vertical edges: face(1) through face(4)
287  for (f = 1; f < 4; f++)
288  {
289  check = 0;
290  for(i = 0; i < m_faces[f]->GetNumEdges(); i++)
291  {
292  for(j = 0; j < m_faces[f+1]->GetNumEdges(); j++)
293  {
294  if( (m_faces[f])->GetEid(i) == (m_faces[f+1])->GetEid(j))
295  {
296  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[f])->GetEdge(i));
297  m_edges.push_back(edge);
298  check++;
299  }
300  }
301  }
302 
303  if( check < 1 )
304  {
305  std::ostringstream errstrm;
306  errstrm << "Connected faces do not share an edge. Faces ";
307  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
308  ASSERTL0(false, errstrm.str());
309  }
310  else if( check > 1)
311  {
312  std::ostringstream errstrm;
313  errstrm << "Connected faces share more than one edge. Faces ";
314  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
315  ASSERTL0(false, errstrm.str());
316  }
317  }
318  }
319 
321  {
322  // Set up the first 2 vertices (i.e. vertex 0,1)
323  if (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0) ||
324  m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1))
325  {
326  m_verts.push_back(m_edges[0]->GetVertex(1));
327  m_verts.push_back(m_edges[0]->GetVertex(0));
328  }
329  else if (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0) ||
330  m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1))
331  {
332  m_verts.push_back(m_edges[0]->GetVertex(0));
333  m_verts.push_back(m_edges[0]->GetVertex(1));
334  }
335  else
336  {
337  std::ostringstream errstrm;
338  errstrm << "Connected edges do not share a vertex. Edges ";
339  errstrm << m_edges[0]->GetEid() << ", " << m_edges[1]->GetEid();
340  ASSERTL0(false, errstrm.str());
341  }
342 
343  // set up the other bottom vertices (i.e. vertex 2,3)
344  for(int i = 1; i < 3; i++)
345  {
346  if (m_edges[i]->GetVid(0) == m_verts[i]->GetVid())
347  {
348  m_verts.push_back(m_edges[i]->GetVertex(1));
349  }
350  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetVid())
351  {
352  m_verts.push_back(m_edges[i]->GetVertex(0));
353  }
354  else
355  {
356  std::ostringstream errstrm;
357  errstrm << "Connected edges do not share a vertex. Edges ";
358  errstrm << m_edges[i]->GetEid() << ", " << m_edges[i-1]->GetEid();
359  ASSERTL0(false, errstrm.str());
360  }
361  }
362 
363  // set up top vertex
364  if (m_edges[4]->GetVid(0) == m_verts[0]->GetVid())
365  {
366  m_verts.push_back(m_edges[4]->GetVertex(1));
367  }
368  else
369  {
370  m_verts.push_back(m_edges[4]->GetVertex(0));
371  }
372 
373  int check = 0;
374  for (int i = 5; i < 8; ++i)
375  {
376  if( (m_edges[i]->GetVid(0) == m_verts[i-4]->GetVid()
377  && m_edges[i]->GetVid(1) == m_verts[4]->GetVid())
378  ||(m_edges[i]->GetVid(1) == m_verts[i-4]->GetVid()
379  && m_edges[i]->GetVid(0) == m_verts[4]->GetVid()))
380  {
381  check++;
382  }
383  }
384  if (check != 3) {
385  std::ostringstream errstrm;
386  errstrm << "Connected edges do not share a vertex. Edges ";
387  errstrm << m_edges[3]->GetEid() << ", " << m_edges[2]->GetEid();
388  ASSERTL0(false, errstrm.str());
389  }
390  }
391 
393  {
394  // This 2D array holds the local id's of all the vertices for every
395  // edge. For every edge, they are ordered to what we define as being
396  // Forwards.
397  const unsigned int edgeVerts[kNedges][2] =
398  { {0,1}, {1,2}, {3,2}, {0,3}, {0,4}, {1,4}, {2,4}, {3,4} };
399 
400  int i;
401  for (i = 0; i < kNedges; i++)
402  {
403  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetVid())
404  {
406  }
407  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetVid())
408  {
410  }
411  else
412  {
413  ASSERTL0(false,"Could not find matching vertex for the edge");
414  }
415  }
416  }
417 
419  {
420  int f,i;
421 
422  // These arrays represent the vector of the A and B coordinate of
423  // the local elemental coordinate system where A corresponds with
424  // the coordinate direction xi_i with the lowest index i (for that
425  // particular face) Coordinate 'B' then corresponds to the other
426  // local coordinate (i.e. with the highest index)
427  Array<OneD,NekDouble> elementAaxis(m_coordim);
428  Array<OneD,NekDouble> elementBaxis(m_coordim);
429 
430  // These arrays correspond to the local coordinate
431  // system of the face itself (i.e. the Geometry2D)
432  // faceAaxis correspond to the xi_0 axis
433  // faceBaxis correspond to the xi_1 axis
434  Array<OneD,NekDouble> faceAaxis(m_coordim);
435  Array<OneD,NekDouble> faceBaxis(m_coordim);
436 
437  // This is the base vertex of the face (i.e. the Geometry2D) This
438  // corresponds to thevertex with local ID 0 of the Geometry2D
439  unsigned int baseVertex;
440 
441  // The lenght of the vectors above
442  NekDouble elementAaxis_length;
443  NekDouble elementBaxis_length;
444  NekDouble faceAaxis_length;
445  NekDouble faceBaxis_length;
446 
447  // This 2D array holds the local id's of all the vertices for every
448  // face. For every face, they are ordered in such a way that the
449  // implementation below allows a unified approach for all faces.
450  const unsigned int faceVerts[kNfaces][4] = {
451  {0,1,2,3},
452  {0,1,4,0}, // Last four elements are triangles which only
453  {1,2,4,0}, // require three vertices.
454  {3,2,4,0},
455  {0,3,4,0}
456  };
457 
458  NekDouble dotproduct1 = 0.0;
459  NekDouble dotproduct2 = 0.0;
460 
461  unsigned int orientation;
462 
463  // Loop over all the faces to set up the orientation
464  for(f = 0; f < kNqfaces + kNtfaces; f++)
465  {
466  // initialisation
467  elementAaxis_length = 0.0;
468  elementBaxis_length = 0.0;
469  faceAaxis_length = 0.0;
470  faceBaxis_length = 0.0;
471 
472  dotproduct1 = 0.0;
473  dotproduct2 = 0.0;
474 
475  baseVertex = m_faces[f]->GetVid(0);
476 
477  // We are going to construct the vectors representing the A and
478  // B axis of every face. These vectors will be constructed as a
479  // vector-representation of the edges of the face. However, for
480  // both coordinate directions, we can represent the vectors by
481  // two different edges. That's why we need to make sure that we
482  // pick the edge to which the baseVertex of the
483  // Geometry2D-representation of the face belongs...
484 
485  // Compute the length of edges on a base-face
486  if (f > 0)
487  {
488  for(i = 0; i < m_coordim; i++)
489  {
490  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
491  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
492  }
493  }
494  else
495  {
496  if( baseVertex == m_verts[ faceVerts[f][0] ]->GetVid() )
497  {
498  for(i = 0; i < m_coordim; i++)
499  {
500  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
501  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
502  }
503  }
504  else if( baseVertex == m_verts[ faceVerts[f][1] ]->GetVid() )
505  {
506  for(i = 0; i < m_coordim; i++)
507  {
508  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
509  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
510  }
511  }
512  else if( baseVertex == m_verts[ faceVerts[f][2] ]->GetVid() )
513  {
514  for(i = 0; i < m_coordim; i++)
515  {
516  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
517  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
518  }
519  }
520  else if( baseVertex == m_verts[ faceVerts[f][3] ]->GetVid() )
521  {
522  for(i = 0; i < m_coordim; i++)
523  {
524  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
525  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
526  }
527  }
528  else
529  {
530  ASSERTL0(false, "Could not find matching vertex for the face");
531  }
532  }
533 
534  // Now, construct the edge-vectors of the local coordinates of
535  // the Geometry2D-representation of the face
536  for(i = 0; i < m_coordim; i++)
537  {
538  int v = m_faces[f]->GetNumVerts()-1;
539  faceAaxis[i] = (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
540  faceBaxis[i] = (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
541 
542  elementAaxis_length += pow(elementAaxis[i],2);
543  elementBaxis_length += pow(elementBaxis[i],2);
544  faceAaxis_length += pow(faceAaxis[i],2);
545  faceBaxis_length += pow(faceBaxis[i],2);
546  }
547 
548  elementAaxis_length = sqrt(elementAaxis_length);
549  elementBaxis_length = sqrt(elementBaxis_length);
550  faceAaxis_length = sqrt(faceAaxis_length);
551  faceBaxis_length = sqrt(faceBaxis_length);
552 
553  // Calculate the inner product of both the A-axis
554  // (i.e. Elemental A axis and face A axis)
555  for(i = 0 ; i < m_coordim; i++)
556  {
557  dotproduct1 += elementAaxis[i]*faceAaxis[i];
558  }
559 
560  orientation = 0;
561  // if the innerproduct is equal to the (absolute value of the ) products of the lengths
562  // of both vectors, then, the coordinate systems will NOT be transposed
563  if( fabs(elementAaxis_length*faceAaxis_length - fabs(dotproduct1)) < NekConstants::kNekZeroTol )
564  {
565  // if the inner product is negative, both A-axis point
566  // in reverse direction
567  if(dotproduct1 < 0.0)
568  {
569  orientation += 2;
570  }
571 
572  // calculate the inner product of both B-axis
573  for(i = 0 ; i < m_coordim; i++)
574  {
575  dotproduct2 += elementBaxis[i]*faceBaxis[i];
576  }
577 
578  // if the inner product is negative, both B-axis point
579  // in reverse direction
580  if( dotproduct2 < 0.0 )
581  {
582  orientation++;
583  }
584  }
585  // The coordinate systems are transposed
586  else
587  {
588  orientation = 4;
589 
590  // Calculate the inner product between the elemental A-axis
591  // and the B-axis of the face (which are now the corresponding axis)
592  dotproduct1 = 0.0;
593  for(i = 0 ; i < m_coordim; i++)
594  {
595  dotproduct1 += elementAaxis[i]*faceBaxis[i];
596  }
597 
598  // check that both these axis are indeed parallel
599  if (fabs(elementAaxis_length*faceBaxis_length
600  - fabs(dotproduct1)) > NekConstants::kNekZeroTol)
601  {
602  cout << "Warning: Prism axes not parallel" << endl;
603  }
604 
605  // if the result is negative, both axis point in reverse
606  // directions
607  if(dotproduct1 < 0.0)
608  {
609  orientation += 2;
610  }
611 
612  // Do the same for the other two corresponding axis
613  dotproduct2 = 0.0;
614  for(i = 0 ; i < m_coordim; i++)
615  {
616  dotproduct2 += elementBaxis[i]*faceAaxis[i];
617  }
618 
619  // check that both these axis are indeed parallel
620  if (fabs(elementBaxis_length*faceAaxis_length
621  - fabs(dotproduct2)) > NekConstants::kNekZeroTol)
622  {
623  cout << "Warning: Prism axes not parallel" << endl;
624  }
625 
626  if( dotproduct2 < 0.0 )
627  {
628  orientation++;
629  }
630  }
631 
632  orientation = orientation + 5;
633 
634  // Fill the m_forient array
635  m_forient[f] = (StdRegions::Orientation) orientation;
636  }
637  }
638 
640  CurveMap &curvedEdges,
641  CurveMap &curvedFaces)
642  {
643  Geometry::v_Reset(curvedEdges, curvedFaces);
644 
645  for (int i = 0; i < 5; ++i)
646  {
647  m_faces[i]->Reset(curvedEdges, curvedFaces);
648  }
649 
650  SetUpXmap();
651  SetUpCoeffs(m_xmap->GetNcoeffs());
652  }
653 
654  /**
655  * @brief Set up the #m_xmap object by determining the order of each
656  * direction from derived faces.
657  */
659  {
660  vector<int> tmp;
661  int order0, order1;
662 
663  if (m_forient[0] < 9)
664  {
665  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
666  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
667  order0 = *max_element(tmp.begin(), tmp.end());
668  }
669  else
670  {
671  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
672  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
673  order0 = *max_element(tmp.begin(), tmp.end());
674  }
675 
676  if (m_forient[0] < 9)
677  {
678  tmp.clear();
679  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
680  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
681  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
682  order1 = *max_element(tmp.begin(), tmp.end());
683  }
684  else
685  {
686  tmp.clear();
687  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
688  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
689  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
690  order1 = *max_element(tmp.begin(), tmp.end());
691  }
692 
693  tmp.clear();
694  tmp.push_back(order0);
695  tmp.push_back(order1);
696  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(1));
697  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(2));
698  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(1));
699  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(2));
700  int order2 = *max_element(tmp.begin(), tmp.end());
701 
702 
703  const LibUtilities::BasisKey A1(
707  const LibUtilities::BasisKey A2(
711  const LibUtilities::BasisKey C(
715 
717  A1, A2, C);
718  }
719  }; //end of namespace
720 }; //end of namespace
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:92
Principle Modified Functions .
Definition: BasisType.h:51
virtual int v_GetDir(const int faceidx, const int facedir) const
Definition: PyrGeom.cpp:198
static const int kNfaces
Definition: PyrGeom.h:60
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
GeomFactorsSharedPtr m_geomFactors
Definition: Geometry.h:170
NekDouble dot(PointGeom &a)
Definition: PointGeom.cpp:192
virtual int v_GetNumEdges() const
Definition: PyrGeom.cpp:188
int GetEid(int i) const
Return the ID of edge i in this element.
Definition: Geometry3D.cpp:71
Principle Modified Functions .
Definition: BasisType.h:49
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:168
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: Geometry.cpp:307
int GetFid(int i) const
Definition: Geometry.h:329
StdRegions::StdExpansionSharedPtr GetXmap() const
Definition: Geometry.h:383
virtual int v_GetNumVerts() const
Definition: PyrGeom.cpp:183
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:93
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: PyrGeom.cpp:639
static const NekDouble kNekZeroTol
static const int kNedges
Definition: PyrGeom.h:57
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
Defines a specification for a set of points.
Definition: Points.h:58
double NekDouble
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: PyrGeom.cpp:124
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
static const int kNtfaces
Definition: PyrGeom.h:59
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:244
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
void Mult(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:177
const Geometry1DSharedPtr GetEdge(int i) const
static const int kNqfaces
Definition: PyrGeom.h:58
3D geometry information
Definition: Geometry3D.h:70
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometry is straight-sided with constant geometric factors.
int GetVid(int i) const
Definition: Geometry.h:319
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
virtual int v_GetNumFaces() const
Definition: PyrGeom.cpp:193
NekDouble dist(PointGeom &a)
Definition: PointGeom.cpp:186
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
GeomFactorsSharedPtr GetMetricInfo()
Definition: Geometry.h:299
GeomType
Indicates the type of element geometry.
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PyrGeom.cpp:658
virtual void v_GenGeomFactors()
Definition: PyrGeom.cpp:80
Geometry is curved or has non-constant factors.
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Describes the specification for a Basis.
Definition: Basis.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50