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