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