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