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