Nektar++
HexGeom.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: HexGeom.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: Hexahedral geometry definition.
32//
33////////////////////////////////////////////////////////////////////////////////
34
41
42using namespace std;
43
45{
46
47const unsigned int HexGeom::VertexEdgeConnectivity[8][3] = {
48 {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7},
49 {4, 8, 11}, {5, 8, 9}, {6, 9, 10}, {7, 10, 11}};
50const unsigned int HexGeom::VertexFaceConnectivity[8][3] = {
51 {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4},
52 {1, 4, 5}, {1, 2, 5}, {2, 3, 5}, {3, 4, 5}};
53const unsigned int HexGeom::EdgeFaceConnectivity[12][2] = {
54 {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2},
55 {2, 3}, {3, 4}, {1, 5}, {2, 5}, {3, 5}, {4, 5}};
56const unsigned int HexGeom::EdgeNormalToFaceVert[6][4] = {
57 {4, 5, 6, 7}, {1, 3, 9, 11}, {0, 2, 8, 10},
58 {1, 3, 9, 11}, {0, 2, 8, 10}, {4, 5, 6, 7}};
59
61{
63}
64
65HexGeom::HexGeom(int id, const QuadGeomSharedPtr faces[])
66 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
67{
69 m_globalID = id;
70
71 /// Copy the face shared pointers
72 m_faces.insert(m_faces.begin(), faces, faces + HexGeom::kNfaces);
73
74 /// Set up orientation vectors with correct amount of elements.
75 m_eorient.resize(kNedges);
76 m_forient.resize(kNfaces);
77
82}
83
85{
86}
87
89{
90 if (!m_setupState)
91 {
93 }
94
96 {
97 GeomType Gtype = eRegular;
98
99 v_FillGeom();
100
101 // check to see if expansions are linear
102 m_straightEdge = 1;
103 if (m_xmap->GetBasisNumModes(0) != 2 ||
104 m_xmap->GetBasisNumModes(1) != 2 ||
105 m_xmap->GetBasisNumModes(2) != 2)
106 {
107 Gtype = eDeformed;
108 m_straightEdge = 0;
109 }
110
111 // check to see if all faces are parallelograms
112 if (Gtype == eRegular)
113 {
115 for (int i = 0; i < 3; ++i)
116 {
118 NekDouble A = (*m_verts[0])(i);
119 NekDouble B = (*m_verts[1])(i);
120 NekDouble C = (*m_verts[2])(i);
121 NekDouble D = (*m_verts[3])(i);
122 NekDouble E = (*m_verts[4])(i);
123 NekDouble F = (*m_verts[5])(i);
124 NekDouble G = (*m_verts[6])(i);
125 NekDouble H = (*m_verts[7])(i);
126 m_isoParameter[i][0] =
127 0.125 * (A + B + C + D + E + F + G + H); // 1
128
129 m_isoParameter[i][1] =
130 0.125 * (-A + B + C - D - E + F + G - H); // xi1
131 m_isoParameter[i][2] =
132 0.125 * (-A - B + C + D - E - F + G + H); // xi2
133 m_isoParameter[i][3] =
134 0.125 * (-A - B - C - D + E + F + G + H); // xi3
135
136 m_isoParameter[i][4] =
137 0.125 * (A - B + C - D + E - F + G - H); // xi1*xi2
138 m_isoParameter[i][5] =
139 0.125 * (A + B - C - D - E - F + G + H); // xi2*xi3
140 m_isoParameter[i][6] =
141 0.125 * (A - B - C + D - E + F + G - H); // xi1*xi3
142
143 m_isoParameter[i][7] =
144 0.125 * (-A + B - C + D + E - F + G - H); // xi1*xi2*xi3
145 NekDouble tmp = fabs(m_isoParameter[i][1]) +
146 fabs(m_isoParameter[i][2]) +
147 fabs(m_isoParameter[i][3]);
149 for (int d = 4; d < 8; ++d)
150 {
151 if (fabs(m_isoParameter[i][d]) > tmp)
152 {
153 Gtype = eDeformed;
154 }
155 }
156 }
157 }
158
159 if (Gtype == eRegular)
160 {
162 }
163
165 Gtype, m_coordim, m_xmap, m_coeffs);
167 }
168}
169
170int HexGeom::v_GetVertexEdgeMap(const int i, const int j) const
171{
172 return VertexEdgeConnectivity[i][j];
173}
174
175int HexGeom::v_GetVertexFaceMap(const int i, const int j) const
176{
177 return VertexFaceConnectivity[i][j];
178}
179
180int HexGeom::v_GetEdgeFaceMap(const int i, const int j) const
181{
182 return EdgeFaceConnectivity[i][j];
183}
184
185int HexGeom::v_GetEdgeNormalToFaceVert(const int i, const int j) const
186{
187 return EdgeNormalToFaceVert[i][j];
188}
189
190int HexGeom::v_GetDir(const int faceidx, const int facedir) const
191{
192 if (faceidx == 0 || faceidx == 5)
193 {
194 return facedir;
195 }
196 else if (faceidx == 1 || faceidx == 3)
197 {
198 return 2 * facedir;
199 }
200 else
201 {
202 return 1 + facedir;
203 }
204}
205
207{
208 // find edge 0
209 int i, j;
210 unsigned int check;
211
212 SegGeomSharedPtr edge;
213
214 // First set up the 4 bottom edges
215 int f;
216 for (f = 1; f < 5; f++)
217 {
218 check = 0;
219 for (i = 0; i < 4; i++)
220 {
221 for (j = 0; j < 4; j++)
222 {
223 if ((m_faces[0])->GetEid(i) == (m_faces[f])->GetEid(j))
224 {
225 edge = std::dynamic_pointer_cast<SegGeom>(
226 (m_faces[0])->GetEdge(i));
227 m_edges.push_back(edge);
228 check++;
229 }
230 }
231 }
232
233 if (check < 1)
234 {
235 std::ostringstream errstrm;
236 errstrm << "Connected faces do not share an edge. Faces ";
237 errstrm << (m_faces[0])->GetGlobalID() << ", "
238 << (m_faces[f])->GetGlobalID();
239 ASSERTL0(false, errstrm.str());
240 }
241 else if (check > 1)
242 {
243 std::ostringstream errstrm;
244 errstrm << "Connected faces share more than one edge. Faces ";
245 errstrm << (m_faces[0])->GetGlobalID() << ", "
246 << (m_faces[f])->GetGlobalID();
247 ASSERTL0(false, errstrm.str());
248 }
249 }
250
251 // Then, set up the 4 vertical edges
252 check = 0;
253 for (i = 0; i < 4; i++)
254 {
255 for (j = 0; j < 4; j++)
256 {
257 if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
258 {
259 edge = std::dynamic_pointer_cast<SegGeom>(
260 (m_faces[1])->GetEdge(i));
261 m_edges.push_back(edge);
262 check++;
263 }
264 }
265 }
266 if (check < 1)
267 {
268 std::ostringstream errstrm;
269 errstrm << "Connected faces do not share an edge. Faces ";
270 errstrm << (m_faces[1])->GetGlobalID() << ", "
271 << (m_faces[4])->GetGlobalID();
272 ASSERTL0(false, errstrm.str());
273 }
274 else if (check > 1)
275 {
276 std::ostringstream errstrm;
277 errstrm << "Connected faces share more than one edge. Faces ";
278 errstrm << (m_faces[1])->GetGlobalID() << ", "
279 << (m_faces[4])->GetGlobalID();
280 ASSERTL0(false, errstrm.str());
281 }
282 for (f = 1; f < 4; f++)
283 {
284 check = 0;
285 for (i = 0; i < 4; i++)
286 {
287 for (j = 0; j < 4; j++)
288 {
289 if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
290 {
291 edge = std::dynamic_pointer_cast<SegGeom>(
292 (m_faces[f])->GetEdge(i));
293 m_edges.push_back(edge);
294 check++;
295 }
296 }
297 }
298
299 if (check < 1)
300 {
301 std::ostringstream errstrm;
302 errstrm << "Connected faces do not share an edge. Faces ";
303 errstrm << (m_faces[f])->GetGlobalID() << ", "
304 << (m_faces[f + 1])->GetGlobalID();
305 ASSERTL0(false, errstrm.str());
306 }
307 else if (check > 1)
308 {
309 std::ostringstream errstrm;
310 errstrm << "Connected faces share more than one edge. Faces ";
311 errstrm << (m_faces[f])->GetGlobalID() << ", "
312 << (m_faces[f + 1])->GetGlobalID();
313 ASSERTL0(false, errstrm.str());
314 }
315 }
316
317 // Finally, set up the 4 top vertices
318 for (f = 1; f < 5; f++)
319 {
320 check = 0;
321 for (i = 0; i < 4; i++)
322 {
323 for (j = 0; j < 4; j++)
324 {
325 if ((m_faces[5])->GetEid(i) == (m_faces[f])->GetEid(j))
326 {
327 edge = std::dynamic_pointer_cast<SegGeom>(
328 (m_faces[5])->GetEdge(i));
329 m_edges.push_back(edge);
330 check++;
331 }
332 }
333 }
334
335 if (check < 1)
336 {
337 std::ostringstream errstrm;
338 errstrm << "Connected faces do not share an edge. Faces ";
339 errstrm << (m_faces[5])->GetGlobalID() << ", "
340 << (m_faces[f])->GetGlobalID();
341 ASSERTL0(false, errstrm.str());
342 }
343 else if (check > 1)
344 {
345 std::ostringstream errstrm;
346 errstrm << "Connected faces share more than one edge. Faces ";
347 errstrm << (m_faces[5])->GetGlobalID() << ", "
348 << (m_faces[f])->GetGlobalID();
349 ASSERTL0(false, errstrm.str());
350 }
351 }
352}
353
355{
356 // Set up the first 2 vertices (i.e. vertex 0,1)
357 if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
358 (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
359 {
360 m_verts.push_back(m_edges[0]->GetVertex(1));
361 m_verts.push_back(m_edges[0]->GetVertex(0));
362 }
363 else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
364 (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
365 {
366 m_verts.push_back(m_edges[0]->GetVertex(0));
367 m_verts.push_back(m_edges[0]->GetVertex(1));
368 }
369 else
370 {
371 std::ostringstream errstrm;
372 errstrm << "Connected edges do not share a vertex. Edges ";
373 errstrm << m_edges[0]->GetGlobalID() << ", "
374 << m_edges[1]->GetGlobalID();
375 ASSERTL0(false, errstrm.str());
376 }
377
378 // set up the other bottom vertices (i.e. vertex 2,3)
379 int i;
380 for (i = 1; i < 3; i++)
381 {
382 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
383 {
384 m_verts.push_back(m_edges[i]->GetVertex(1));
385 }
386 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
387 {
388 m_verts.push_back(m_edges[i]->GetVertex(0));
389 }
390 else
391 {
392 std::ostringstream errstrm;
393 errstrm << "Connected edges do not share a vertex. Edges ";
394 errstrm << m_edges[i]->GetGlobalID() << ", "
395 << m_edges[i - 1]->GetGlobalID();
396 ASSERTL0(false, errstrm.str());
397 }
398 }
399
400 // set up top vertices
401 // First, set up vertices 4,5
402 if ((m_edges[8]->GetVid(0) == m_edges[9]->GetVid(0)) ||
403 (m_edges[8]->GetVid(0) == m_edges[9]->GetVid(1)))
404 {
405 m_verts.push_back(m_edges[8]->GetVertex(1));
406 m_verts.push_back(m_edges[8]->GetVertex(0));
407 }
408 else if ((m_edges[8]->GetVid(1) == m_edges[9]->GetVid(0)) ||
409 (m_edges[8]->GetVid(1) == m_edges[9]->GetVid(1)))
410 {
411 m_verts.push_back(m_edges[8]->GetVertex(0));
412 m_verts.push_back(m_edges[8]->GetVertex(1));
413 }
414 else
415 {
416 std::ostringstream errstrm;
417 errstrm << "Connected edges do not share a vertex. Edges ";
418 errstrm << m_edges[8]->GetGlobalID() << ", "
419 << m_edges[9]->GetGlobalID();
420 ASSERTL0(false, errstrm.str());
421 }
422
423 // set up the other top vertices (i.e. vertex 6,7)
424 for (i = 9; i < 11; i++)
425 {
426 if (m_edges[i]->GetVid(0) == m_verts[i - 4]->GetGlobalID())
427 {
428 m_verts.push_back(m_edges[i]->GetVertex(1));
429 }
430 else if (m_edges[i]->GetVid(1) == m_verts[i - 4]->GetGlobalID())
431 {
432 m_verts.push_back(m_edges[i]->GetVertex(0));
433 }
434 else
435 {
436 std::ostringstream errstrm;
437 errstrm << "Connected edges do not share a vertex. Edges ";
438 errstrm << m_edges[i]->GetGlobalID() << ", "
439 << m_edges[i - 1]->GetGlobalID();
440 ASSERTL0(false, errstrm.str());
441 }
442 }
443}
444
446{
447 int f, i;
448
449 // These arrays represent the vector of the A and B
450 // coordinate of the local elemental coordinate system
451 // where A corresponds with the coordinate direction xi_i
452 // with the lowest index i (for that particular face)
453 // Coordinate 'B' then corresponds to the other local
454 // coordinate (i.e. with the highest index)
455 Array<OneD, NekDouble> elementAaxis(m_coordim);
456 Array<OneD, NekDouble> elementBaxis(m_coordim);
457
458 // These arrays correspond to the local coordinate
459 // system of the face itself (i.e. the Geometry2D)
460 // faceAaxis correspond to the xi_0 axis
461 // faceBaxis correspond to the xi_1 axis
464
465 // This is the base vertex of the face (i.e. the Geometry2D)
466 // This corresponds to thevertex with local ID 0 of the
467 // Geometry2D
468 unsigned int baseVertex;
469
470 // The lenght of the vectors above
471 NekDouble elementAaxis_length;
472 NekDouble elementBaxis_length;
473 NekDouble faceAaxis_length;
474 NekDouble faceBaxis_length;
475
476 // This 2D array holds the local id's of all the vertices
477 // for every face. For every face, they are ordered in such
478 // a way that the implementation below allows a unified approach
479 // for all faces.
480 const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {
481 {0, 1, 2, 3}, {0, 1, 5, 4}, {1, 2, 6, 5},
482 {3, 2, 6, 7}, {0, 3, 7, 4}, {4, 5, 6, 7}};
483
484 NekDouble dotproduct1 = 0.0;
485 NekDouble dotproduct2 = 0.0;
486
487 unsigned int orientation;
488
489 // Loop over all the faces to set up the orientation
490 for (f = 0; f < kNqfaces + kNtfaces; f++)
491 {
492 // initialisation
493 elementAaxis_length = 0.0;
494 elementBaxis_length = 0.0;
495 faceAaxis_length = 0.0;
496 faceBaxis_length = 0.0;
497
498 dotproduct1 = 0.0;
499 dotproduct2 = 0.0;
500
501 baseVertex = m_faces[f]->GetVid(0);
502
503 // We are going to construct the vectors representing the A and B axis
504 // of every face. These vectors will be constructed as a
505 // vector-representation
506 // of the edges of the face. However, for both coordinate directions, we
507 // can
508 // represent the vectors by two different edges. That's why we need to
509 // make sure that
510 // we pick the edge to which the baseVertex of the
511 // Geometry2D-representation of the face
512 // belongs...
513 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
514 {
515 for (i = 0; i < m_coordim; i++)
516 {
517 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
518 (*m_verts[faceVerts[f][0]])[i];
519 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
520 (*m_verts[faceVerts[f][0]])[i];
521 }
522 }
523 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
524 {
525 for (i = 0; i < m_coordim; i++)
526 {
527 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
528 (*m_verts[faceVerts[f][0]])[i];
529 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
530 (*m_verts[faceVerts[f][1]])[i];
531 }
532 }
533 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
534 {
535 for (i = 0; i < m_coordim; i++)
536 {
537 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
538 (*m_verts[faceVerts[f][3]])[i];
539 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
540 (*m_verts[faceVerts[f][1]])[i];
541 }
542 }
543 else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
544 {
545 for (i = 0; i < m_coordim; i++)
546 {
547 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
548 (*m_verts[faceVerts[f][3]])[i];
549 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
550 (*m_verts[faceVerts[f][0]])[i];
551 }
552 }
553 else
554 {
555 ASSERTL0(false, "Could not find matching vertex for the face");
556 }
557
558 // Now, construct the edge-vectors of the local coordinates of
559 // the Geometry2D-representation of the face
560 for (i = 0; i < m_coordim; i++)
561 {
562 faceAaxis[i] =
563 (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
564 faceBaxis[i] =
565 (*m_faces[f]->GetVertex(3))[i] - (*m_faces[f]->GetVertex(0))[i];
566
567 elementAaxis_length += pow(elementAaxis[i], 2);
568 elementBaxis_length += pow(elementBaxis[i], 2);
569 faceAaxis_length += pow(faceAaxis[i], 2);
570 faceBaxis_length += pow(faceBaxis[i], 2);
571 }
572
573 elementAaxis_length = sqrt(elementAaxis_length);
574 elementBaxis_length = sqrt(elementBaxis_length);
575 faceAaxis_length = sqrt(faceAaxis_length);
576 faceBaxis_length = sqrt(faceBaxis_length);
577
578 // Calculate the inner product of both the A-axis
579 // (i.e. Elemental A axis and face A axis)
580 for (i = 0; i < m_coordim; i++)
581 {
582 dotproduct1 += elementAaxis[i] * faceAaxis[i];
583 }
584
585 NekDouble norm =
586 fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
587 orientation = 0;
588
589 // if the innerproduct is equal to the (absolute value of the ) products
590 // of the lengths of both vectors, then, the coordinate systems will NOT
591 // be transposed
592 if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
593 {
594 // if the inner product is negative, both A-axis point
595 // in reverse direction
596 if (dotproduct1 < 0.0)
597 {
598 orientation += 2;
599 }
600
601 // calculate the inner product of both B-axis
602 for (i = 0; i < m_coordim; i++)
603 {
604 dotproduct2 += elementBaxis[i] * faceBaxis[i];
605 }
606
607 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
608
609 // check that both these axis are indeed parallel
610 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
611 "These vectors should be parallel");
612
613 // if the inner product is negative, both B-axis point
614 // in reverse direction
615 if (dotproduct2 < 0.0)
616 {
617 orientation++;
618 }
619 }
620 // The coordinate systems are transposed
621 else
622 {
623 orientation = 4;
624
625 // Calculate the inner product between the elemental A-axis
626 // and the B-axis of the face (which are now the corresponding axis)
627 dotproduct1 = 0.0;
628 for (i = 0; i < m_coordim; i++)
629 {
630 dotproduct1 += elementAaxis[i] * faceBaxis[i];
631 }
632
633 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
634
635 // check that both these axis are indeed parallel
636 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
637 "These vectors should be parallel");
638
639 // if the result is negative, both axis point in reverse
640 // directions
641 if (dotproduct1 < 0.0)
642 {
643 orientation += 2;
644 }
645
646 // Do the same for the other two corresponding axis
647 dotproduct2 = 0.0;
648 for (i = 0; i < m_coordim; i++)
649 {
650 dotproduct2 += elementBaxis[i] * faceAaxis[i];
651 }
652
653 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
654
655 // check that both these axis are indeed parallel
656 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
657 "These vectors should be parallel");
658
659 if (dotproduct2 < 0.0)
660 {
661 orientation++;
662 }
663 }
664
665 orientation = orientation + 5;
666 // Fill the m_forient array
667 m_forient[f] = (StdRegions::Orientation)orientation;
668 }
669}
670
672{
673
674 // This 2D array holds the local id's of all the vertices
675 // for every edge. For every edge, they are ordered to what we
676 // define as being Forwards
677 const unsigned int edgeVerts[kNedges][2] = {{0, 1}, {1, 2}, {2, 3}, {3, 0},
678 {0, 4}, {1, 5}, {2, 6}, {3, 7},
679 {4, 5}, {5, 6}, {6, 7}, {7, 4}};
680
681 int i;
682 for (i = 0; i < kNedges; i++)
683 {
684 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
685 {
687 }
688 else if (m_edges[i]->GetVid(0) ==
689 m_verts[edgeVerts[i][1]]->GetGlobalID())
690 {
692 }
693 else
694 {
695 ASSERTL0(false, "Could not find matching vertex for the edge");
696 }
697 }
698}
699
700void HexGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
701{
702 Geometry::v_Reset(curvedEdges, curvedFaces);
703
704 for (int i = 0; i < 6; ++i)
705 {
706 m_faces[i]->Reset(curvedEdges, curvedFaces);
707 }
708
709 SetUpXmap();
710 SetUpCoeffs(m_xmap->GetNcoeffs());
711}
712
714{
715 if (!m_setupState)
716 {
717 for (int i = 0; i < 6; ++i)
718 {
719 m_faces[i]->Setup();
720 }
721 SetUpXmap();
722 SetUpCoeffs(m_xmap->GetNcoeffs());
723 m_setupState = true;
724 }
725}
726
727/**
728 * @brief Set up the #m_xmap object by determining the order of each
729 * direction from derived faces.
730 */
732{
733 // Determine necessary order for standard region. This can almost certainly
734 // be simplified but works for now!
735 vector<int> tmp1;
736
737 if (m_forient[0] < 9)
738 {
739 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
740 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
741 }
742 else
743 {
744 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
745 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
746 }
747
748 if (m_forient[5] < 9)
749 {
750 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
751 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
752 }
753 else
754 {
755 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
756 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
757 }
758
759 int order0 = *max_element(tmp1.begin(), tmp1.end());
760
761 tmp1.clear();
762
763 if (m_forient[0] < 9)
764 {
765 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
766 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
767 }
768 else
769 {
770 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
771 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
772 }
773
774 if (m_forient[5] < 9)
775 {
776 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
777 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
778 }
779 else
780 {
781 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
782 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
783 }
784
785 int order1 = *max_element(tmp1.begin(), tmp1.end());
786
787 tmp1.clear();
788
789 if (m_forient[1] < 9)
790 {
791 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
792 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(3));
793 }
794 else
795 {
796 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(0));
797 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
798 }
799
800 if (m_forient[3] < 9)
801 {
802 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
803 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(3));
804 }
805 else
806 {
807 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(0));
808 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
809 }
810
811 int order2 = *max_element(tmp1.begin(), tmp1.end());
812
815 LibUtilities::PointsKey(order0 + 1,
819 LibUtilities::PointsKey(order1 + 1,
823 LibUtilities::PointsKey(order2 + 1,
825
827}
828
829} // 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_GetDir(const int faceidx, const int facedir) const override
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: HexGeom.cpp:190
static const unsigned int EdgeFaceConnectivity[12][2]
Definition: HexGeom.h:80
static const int kNfaces
Definition: HexGeom.h:58
static const unsigned int VertexEdgeConnectivity[8][3]
Definition: HexGeom.h:78
static const unsigned int VertexFaceConnectivity[8][3]
Definition: HexGeom.h:79
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: HexGeom.cpp:700
void v_GenGeomFactors() override
Definition: HexGeom.cpp:88
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: HexGeom.cpp:185
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: HexGeom.cpp:180
static const int kNqfaces
Definition: HexGeom.h:56
static const int kNedges
Definition: HexGeom.h:55
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: HexGeom.cpp:170
static const unsigned int EdgeNormalToFaceVert[6][4]
Definition: HexGeom.h:81
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: HexGeom.cpp:175
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: HexGeom.cpp:731
static const int kNtfaces
Definition: HexGeom.h:57
static const int kNverts
Definition: QuadGeom.h:73
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
static const NekDouble kNekZeroTol
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: HexGeom.h:45
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::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294