Nektar++
Loading...
Searching...
No Matches
PrismGeom.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: PrismGeom.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: Prismatic geometry definition.
32//
33////////////////////////////////////////////////////////////////////////////////
34
42
44{
45
46const unsigned int PrismGeom::VertexEdgeConnectivity[6][3] = {
47 {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7}, {4, 5, 8}, {6, 7, 8}};
48const unsigned int PrismGeom::VertexFaceConnectivity[6][3] = {
49 {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4}, {1, 2, 4}, {2, 3, 4}};
50const unsigned int PrismGeom::EdgeFaceConnectivity[9][2] = {
51 {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2}, {2, 3}, {3, 4}, {2, 4}};
52const unsigned int PrismGeom::EdgeNormalToFaceVert[5][4] = {
53 {4, 5, 6, 7}, {1, 3, 8, -1}, {0, 2, 4, 7}, {1, 3, 8, -1}, {0, 2, 5, 6}};
54
60
65
66PrismGeom::PrismGeom(int id, std::array<Geometry2D *, kNfaces> faces)
67 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
68{
70 m_globalID = id;
71 m_faces = faces;
72
73 /// Set up local objects.
78}
79
80int PrismGeom::v_GetDir(const int faceidx, const int facedir) const
81{
82 if (faceidx == 0)
83 {
84 return facedir;
85 }
86 else if (faceidx == 1 || faceidx == 3)
87 {
88 return 2 * facedir;
89 }
90 else
91 {
92 return 1 + facedir;
93 }
94}
95
97{
98 if (!m_setupState)
99 {
101 }
102
104 {
105 GeomType Gtype = eRegular;
106
107 v_FillGeom();
108
109 // check to see if expansions are linear
110 m_straightEdge = 1;
111 if (m_xmap->GetBasisNumModes(0) != 2 ||
112 m_xmap->GetBasisNumModes(1) != 2 ||
113 m_xmap->GetBasisNumModes(2) != 2)
114 {
115 Gtype = eDeformed;
116 m_straightEdge = 0;
117 }
118
119 // check to see if all quadrilateral faces are parallelograms
120 if (Gtype == eRegular)
121 {
123 for (int i = 0; i < 3; ++i)
124 {
126 NekDouble A = (*m_verts[0])(i);
127 NekDouble B = (*m_verts[1])(i);
128 NekDouble C = (*m_verts[2])(i);
129 NekDouble D = (*m_verts[3])(i);
130 NekDouble E = (*m_verts[4])(i);
131 NekDouble F = (*m_verts[5])(i);
132 m_isoParameter[i][0] = 0.25 * (B + C + E + F);
133
134 m_isoParameter[i][1] = 0.25 * (-A + B + C - D); // xi1
135 m_isoParameter[i][2] = 0.25 * (-B + C - E + F); // xi2
136 m_isoParameter[i][3] = 0.25 * (-A - D + E + F); // xi3
137
138 m_isoParameter[i][4] = 0.25 * (A - B + C - D); // xi1*xi2
139 m_isoParameter[i][5] = 0.25 * (A - D - E + F); // xi2*xi3
140 NekDouble tmp = fabs(m_isoParameter[i][1]) +
141 fabs(m_isoParameter[i][2]) +
142 fabs(m_isoParameter[i][3]);
144 for (int d = 4; d < 6; ++d)
145 {
146 if (fabs(m_isoParameter[i][d]) > tmp)
147 {
148 Gtype = eDeformed;
149 }
150 }
151 }
152 }
153
154 if (Gtype == eRegular)
155 {
157 }
158
160 Gtype, m_coordim, m_xmap, m_coeffs);
162 }
163}
164
165int PrismGeom::v_GetVertexEdgeMap(const int i, const int j) const
166{
167 return VertexEdgeConnectivity[i][j];
168}
169
170int PrismGeom::v_GetVertexFaceMap(const int i, const int j) const
171{
172 return VertexFaceConnectivity[i][j];
173}
174
175int PrismGeom::v_GetEdgeFaceMap(const int i, const int j) const
176{
177 return EdgeFaceConnectivity[i][j];
178}
179
180int PrismGeom::v_GetEdgeNormalToFaceVert(const int i, const int j) const
181{
182 return EdgeNormalToFaceVert[i][j];
183}
184
186{
187 // find edge 0
188 int i, j;
189 unsigned int check;
190
191 // First set up the 4 bottom edges
192 int f; // Connected face index
193 for (f = 1; f < 5; f++)
194 {
195 int nEdges = m_faces[f]->GetNumEdges();
196 check = 0;
197 for (i = 0; i < 4; i++)
198 {
199 for (j = 0; j < nEdges; j++)
200 {
201 if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
202 {
203 m_edges[f - 1] =
204 static_cast<SegGeom *>((m_faces[0])->GetEdge(i));
205 check++;
206 }
207 }
208 }
209
210 if (check < 1)
211 {
212 std::ostringstream errstrm;
213 errstrm << "Connected faces do not share an edge. Faces ";
214 errstrm << (m_faces[0])->GetGlobalID() << ", "
215 << (m_faces[f])->GetGlobalID();
216 ASSERTL0(false, errstrm.str());
217 }
218 else if (check > 1)
219 {
220 std::ostringstream errstrm;
221 errstrm << "Connected faces share more than one edge. Faces ";
222 errstrm << (m_faces[0])->GetGlobalID() << ", "
223 << (m_faces[f])->GetGlobalID();
224 ASSERTL0(false, errstrm.str());
225 }
226 }
227
228 // Then, set up the 4 vertical edges
229 check = 0;
230 for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(4)
231 {
232 for (j = 0; j < 4; j++)
233 {
234 if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
235 {
236 m_edges[4] = static_cast<SegGeom *>((m_faces[1])->GetEdge(i));
237 check++;
238 }
239 }
240 }
241 if (check < 1)
242 {
243 std::ostringstream errstrm;
244 errstrm << "Connected faces do not share an edge. Faces ";
245 errstrm << (m_faces[1])->GetGlobalID() << ", "
246 << (m_faces[4])->GetGlobalID();
247 ASSERTL0(false, errstrm.str());
248 }
249 else if (check > 1)
250 {
251 std::ostringstream errstrm;
252 errstrm << "Connected faces share more than one edge. Faces ";
253 errstrm << (m_faces[1])->GetGlobalID() << ", "
254 << (m_faces[4])->GetGlobalID();
255 ASSERTL0(false, errstrm.str());
256 }
257 // Set up vertical edges: face(1) through face(4)
258 for (f = 1; f < 4; f++)
259 {
260 check = 0;
261 for (i = 0; i < m_faces[f]->GetNumEdges(); i++)
262 {
263 for (j = 0; j < m_faces[f + 1]->GetNumEdges(); j++)
264 {
265 if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
266 {
267 m_edges[f + 4] =
268 static_cast<SegGeom *>((m_faces[f])->GetEdge(i));
269 check++;
270 }
271 }
272 }
273
274 if (check < 1)
275 {
276 std::ostringstream errstrm;
277 errstrm << "Connected faces do not share an edge. Faces ";
278 errstrm << (m_faces[f])->GetGlobalID() << ", "
279 << (m_faces[f + 1])->GetGlobalID();
280 ASSERTL0(false, errstrm.str());
281 }
282 else if (check > 1)
283 {
284 std::ostringstream errstrm;
285 errstrm << "Connected faces share more than one edge. Faces ";
286 errstrm << (m_faces[f])->GetGlobalID() << ", "
287 << (m_faces[f + 1])->GetGlobalID();
288 ASSERTL0(false, errstrm.str());
289 }
290 }
291
292 // Finally, set up the 1 top edge
293 check = 0;
294 for (i = 0; i < 4; i++)
295 {
296 for (j = 0; j < 4; j++)
297 {
298 if ((m_faces[2])->GetEid(i) == (m_faces[4])->GetEid(j))
299 {
300 m_edges[8] = static_cast<SegGeom *>((m_faces[2])->GetEdge(i));
301 check++;
302 }
303 }
304 }
305
306 if (check < 1)
307 {
308 std::ostringstream errstrm;
309 errstrm << "Connected faces do not share an edge. Faces ";
310 errstrm << (m_faces[1])->GetGlobalID() << ", "
311 << (m_faces[3])->GetGlobalID();
312 ASSERTL0(false, errstrm.str());
313 }
314 else if (check > 1)
315 {
316 std::ostringstream errstrm;
317 errstrm << "Connected faces share more than one edge. Faces ";
318 errstrm << (m_faces[1])->GetGlobalID() << ", "
319 << (m_faces[3])->GetGlobalID();
320 ASSERTL0(false, errstrm.str());
321 }
322}
323
325{
326
327 // Set up the first 2 vertices (i.e. vertex 0,1)
328 if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
329 (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
330 {
331 m_verts[0] = m_edges[0]->GetVertex(1);
332 m_verts[1] = m_edges[0]->GetVertex(0);
333 }
334 else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
335 (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
336 {
337 m_verts[0] = m_edges[0]->GetVertex(0);
338 m_verts[1] = m_edges[0]->GetVertex(1);
339 }
340 else
341 {
342 std::ostringstream errstrm;
343 errstrm << "Connected edges do not share a vertex. Edges ";
344 errstrm << m_edges[0]->GetGlobalID() << ", "
345 << m_edges[1]->GetGlobalID();
346 ASSERTL0(false, errstrm.str());
347 }
348
349 // set up the other bottom vertices (i.e. vertex 2,3)
350 for (int i = 1; i < 3; i++)
351 {
352 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
353 {
354 m_verts[i + 1] = m_edges[i]->GetVertex(1);
355 }
356 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
357 {
358 m_verts[i + 1] = m_edges[i]->GetVertex(0);
359 }
360 else
361 {
362 std::ostringstream errstrm;
363 errstrm << "Connected edges do not share a vertex. Edges ";
364 errstrm << m_edges[i]->GetGlobalID() << ", "
365 << m_edges[i - 1]->GetGlobalID();
366 ASSERTL0(false, errstrm.str());
367 }
368 }
369
370 // set up top vertices
371 // First, set up vertices 4,5
372 if ((m_edges[8]->GetVid(0) == m_edges[4]->GetVid(0)) ||
373 (m_edges[8]->GetVid(0) == m_edges[4]->GetVid(1)))
374 {
375 m_verts[4] = m_edges[8]->GetVertex(0);
376 m_verts[5] = m_edges[8]->GetVertex(1);
377 }
378 else if ((m_edges[8]->GetVid(1) == m_edges[4]->GetVid(0)) ||
379 (m_edges[8]->GetVid(1) == m_edges[4]->GetVid(1)))
380 {
381 m_verts[4] = m_edges[8]->GetVertex(1);
382 m_verts[5] = m_edges[8]->GetVertex(0);
383 }
384 else
385 {
386 std::ostringstream errstrm;
387 errstrm << "Connected edges do not share a vertex. Edges ";
388 errstrm << m_edges[8]->GetGlobalID();
389 ASSERTL0(false, errstrm.str());
390 }
391}
392
394{
395
396 // This 2D array holds the local id's of all the vertices
397 // for every edge. For every edge, they are ordered to what we
398 // define as being Forwards
399 const unsigned int edgeVerts[kNedges][2] = {
400 {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 5}, {3, 5}, {4, 5}};
401
402 int i;
403 for (i = 0; i < kNedges; i++)
404 {
405 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
406 {
408 }
409 else if (m_edges[i]->GetVid(0) ==
410 m_verts[edgeVerts[i][1]]->GetGlobalID())
411 {
413 }
414 else
415 {
416 ASSERTL0(false, "Could not find matching vertex for the edge");
417 }
418 }
419}
420
422{
423 int f, i;
424
425 // These arrays represent the vector of the A and B
426 // coordinate of the local elemental coordinate system
427 // where A corresponds with the coordinate direction xi_i
428 // with the lowest index i (for that particular face)
429 // Coordinate 'B' then corresponds to the other local
430 // coordinate (i.e. with the highest index)
431 Array<OneD, NekDouble> elementAaxis(m_coordim);
432 Array<OneD, NekDouble> elementBaxis(m_coordim);
433
434 // These arrays correspond to the local coordinate
435 // system of the face itself (i.e. the Geometry2D)
436 // faceAaxis correspond to the xi_0 axis
437 // faceBaxis correspond to the xi_1 axis
440
441 // This is the base vertex of the face (i.e. the Geometry2D)
442 // This corresponds to thevertex with local ID 0 of the
443 // Geometry2D
444 unsigned int baseVertex;
445
446 // The lenght of the vectors above
447 NekDouble elementAaxis_length;
448 NekDouble elementBaxis_length;
449 NekDouble faceAaxis_length;
450 NekDouble faceBaxis_length;
451
452 // This 2D array holds the local id's of all the vertices
453 // for every face. For every face, they are ordered in such
454 // a way that the implementation below allows a unified approach
455 // for all faces.
456 const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {
457 {0, 1, 2, 3},
458 {0, 1, 4, 0}, // This is triangle requires only three vertices
459 {1, 2, 5, 4},
460 {3, 2, 5, 0}, // This is triangle requires only three vertices
461 {0, 3, 5, 4},
462 };
463
464 NekDouble dotproduct1 = 0.0;
465 NekDouble dotproduct2 = 0.0;
466
467 unsigned int orientation;
468
469 // Loop over all the faces to set up the orientation
470 for (f = 0; f < kNqfaces + kNtfaces; f++)
471 {
472 // initialisation
473 elementAaxis_length = 0.0;
474 elementBaxis_length = 0.0;
475 faceAaxis_length = 0.0;
476 faceBaxis_length = 0.0;
477
478 dotproduct1 = 0.0;
479 dotproduct2 = 0.0;
480
481 baseVertex = m_faces[f]->GetVid(0);
482
483 // We are going to construct the vectors representing the A and B axis
484 // of every face. These vectors will be constructed as a
485 // vector-representation
486 // of the edges of the face. However, for both coordinate directions, we
487 // can
488 // represent the vectors by two different edges. That's why we need to
489 // make sure that
490 // we pick the edge to which the baseVertex of the
491 // Geometry2D-representation of the face
492 // belongs...
493
494 // Compute the length of edges on a base-face
495 if (f == 1 || f == 3)
496 { // Face is a Triangle
497 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
498 {
499 for (i = 0; i < m_coordim; i++)
500 {
501 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
502 (*m_verts[faceVerts[f][0]])[i];
503 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
504 (*m_verts[faceVerts[f][0]])[i];
505 }
506 }
507 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
508 {
509 for (i = 0; i < m_coordim; i++)
510 {
511 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
512 (*m_verts[faceVerts[f][0]])[i];
513 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
514 (*m_verts[faceVerts[f][1]])[i];
515 }
516 }
517 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
518 {
519 for (i = 0; i < m_coordim; i++)
520 {
521 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
522 (*m_verts[faceVerts[f][2]])[i];
523 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
524 (*m_verts[faceVerts[f][0]])[i];
525 }
526 }
527 else
528 {
529 ASSERTL0(false, "Could not find matching vertex for the face");
530 }
531 }
532 else
533 { // Face is a Quad
534 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
535 {
536 for (i = 0; i < m_coordim; i++)
537 {
538 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
539 (*m_verts[faceVerts[f][0]])[i];
540 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
541 (*m_verts[faceVerts[f][0]])[i];
542 }
543 }
544 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
545 {
546 for (i = 0; i < m_coordim; i++)
547 {
548 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
549 (*m_verts[faceVerts[f][0]])[i];
550 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
551 (*m_verts[faceVerts[f][1]])[i];
552 }
553 }
554 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
555 {
556 for (i = 0; i < m_coordim; i++)
557 {
558 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
559 (*m_verts[faceVerts[f][3]])[i];
560 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
561 (*m_verts[faceVerts[f][1]])[i];
562 }
563 }
564 else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
565 {
566 for (i = 0; i < m_coordim; i++)
567 {
568 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
569 (*m_verts[faceVerts[f][3]])[i];
570 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
571 (*m_verts[faceVerts[f][0]])[i];
572 }
573 }
574 else
575 {
576 ASSERTL0(false, "Could not find matching vertex for the face");
577 }
578 }
579 // Now, construct the edge-vectors of the local coordinates of
580 // the Geometry2D-representation of the face
581 for (i = 0; i < m_coordim; i++)
582 {
583 int v = m_faces[f]->GetNumVerts() - 1;
584 faceAaxis[i] =
585 (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
586 faceBaxis[i] =
587 (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
588
589 elementAaxis_length += pow(elementAaxis[i], 2);
590 elementBaxis_length += pow(elementBaxis[i], 2);
591 faceAaxis_length += pow(faceAaxis[i], 2);
592 faceBaxis_length += pow(faceBaxis[i], 2);
593 }
594
595 elementAaxis_length = sqrt(elementAaxis_length);
596 elementBaxis_length = sqrt(elementBaxis_length);
597 faceAaxis_length = sqrt(faceAaxis_length);
598 faceBaxis_length = sqrt(faceBaxis_length);
599
600 // Calculate the inner product of both the A-axis
601 // (i.e. Elemental A axis and face A axis)
602 for (i = 0; i < m_coordim; i++)
603 {
604 dotproduct1 += elementAaxis[i] * faceAaxis[i];
605 }
606
607 orientation = 0;
608
609 // if the innerproduct is equal to the (absolute value of the ) products
610 // of the lengths of both vectors, then, the coordinate systems will NOT
611 // be transposed
612 if (fabs(elementAaxis_length * faceAaxis_length - fabs(dotproduct1)) <
614 {
615 // if the inner product is negative, both A-axis point
616 // in reverse direction
617 if (dotproduct1 < 0.0)
618 {
619 orientation += 2;
620 }
621
622 // calculate the inner product of both B-axis
623 for (i = 0; i < m_coordim; i++)
624 {
625 dotproduct2 += elementBaxis[i] * faceBaxis[i];
626 }
627
628 ASSERTL1(fabs(fabs(dotproduct2 / elementBaxis_length /
629 faceBaxis_length) -
631 "These vectors should be parallel");
632
633 // if the inner product is negative, both B-axis point
634 // in reverse direction
635 if (dotproduct2 < 0.0)
636 {
637 orientation++;
638 }
639 }
640 // The coordinate systems are transposed
641 else
642 {
643 orientation = 4;
644
645 // Calculate the inner product between the elemental A-axis
646 // and the B-axis of the face (which are now the corresponding axis)
647 dotproduct1 = 0.0;
648 for (i = 0; i < m_coordim; i++)
649 {
650 dotproduct1 += elementAaxis[i] * faceBaxis[i];
651 }
652
653 // check that both these axis are indeed parallel
654 ASSERTL1(fabs(fabs(dotproduct1) / elementAaxis_length /
655 faceBaxis_length -
657 "These vectors should be parallel");
658
659 // if the result is negative, both axis point in reverse
660 // directions
661 if (dotproduct1 < 0.0)
662 {
663 orientation += 2;
664 }
665
666 // Do the same for the other two corresponding axis
667 dotproduct2 = 0.0;
668 for (i = 0; i < m_coordim; i++)
669 {
670 dotproduct2 += elementBaxis[i] * faceAaxis[i];
671 }
672
673 ASSERTL1(fabs(fabs(dotproduct2) / elementBaxis_length /
674 faceAaxis_length -
676 "These vectors should be parallel");
677
678 if (dotproduct2 < 0.0)
679 {
680 orientation++;
681 }
682 }
683
684 orientation = orientation + 5;
685
686 if ((f == 1) || (f == 3)) // check triange orientation
687 {
688 ASSERTL0(
690 "Orientation of triangular face (id = " +
691 std::to_string(m_faces[f]->GetGlobalID()) +
692 ") is inconsistent with face " + std::to_string(f) +
693 " of prism element (id = " + std::to_string(m_globalID) +
694 ") since Dir2 is aligned with Dir1. Mesh setup "
695 "needs investigation");
696 }
697
698 // Fill the m_forient array
699 m_forient[f] = (StdRegions::Orientation)orientation;
700 }
701}
702
703void PrismGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
704{
705 Geometry::v_Reset(curvedEdges, curvedFaces);
706
707 for (int i = 0; i < 5; ++i)
708 {
709 m_faces[i]->Reset(curvedEdges, curvedFaces);
710 }
711}
712
714{
715 if (!m_setupState)
716 {
717 for (int i = 0; i < 5; ++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 std::vector<int> tmp;
734 int order0, order1;
735
736 if (m_forient[0] < 9)
737 {
738 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
739 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
740 order0 = *std::max_element(tmp.begin(), tmp.end());
741 }
742 else
743 {
744 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
745 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
746 order0 = *std::max_element(tmp.begin(), tmp.end());
747 }
748
749 if (m_forient[0] < 9)
750 {
751 tmp.clear();
752 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
753 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
754 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
755 order1 = *std::max_element(tmp.begin(), tmp.end());
756 }
757 else
758 {
759 tmp.clear();
760 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
761 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
762 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
763 order1 = *std::max_element(tmp.begin(), tmp.end());
764 }
765
766 tmp.clear();
767 tmp.push_back(order0);
768 tmp.push_back(order1);
769 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
770 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
771 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
772 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
773 int order2 = *std::max_element(tmp.begin(), tmp.end());
774
775 std::array<LibUtilities::BasisKey, 3> basis = {
778 LibUtilities::PointsKey(order0 + 1,
782 LibUtilities::PointsKey(order1 + 1,
787 LibUtilities::eGaussRadauMAlpha1Beta0))};
788
789 m_xmap = GetStdPrismFactory().CreateInstance(basis);
790}
791
792/**
793 * @brief Put all quadrature information into face/edge structure and
794 * backward transform.
795 *
796 * Note verts, edges, and faces are listed according to anticlockwise
797 * convention but points in _coeffs have to be in array format from left
798 * to right.
799 */
801{
802 if (m_state == ePtsFilled)
803 {
804 return;
805 }
806
807 int i, j, k;
808
809 for (i = 0; i < kNfaces; i++)
810 {
811 m_faces[i]->FillGeom();
812
813 int nFaceCoeffs = m_faces[i]->GetXmap()->GetNcoeffs();
814
815 Array<OneD, unsigned int> mapArray(nFaceCoeffs);
816 Array<OneD, int> signArray(nFaceCoeffs);
817
818 if (m_forient[i] < 9)
819 {
820 m_xmap->GetTraceToElementMap(
821 i, mapArray, signArray, m_forient[i],
822 m_faces[i]->GetXmap()->GetTraceNcoeffs(0),
823 m_faces[i]->GetXmap()->GetTraceNcoeffs(1));
824 }
825 else
826 {
827 m_xmap->GetTraceToElementMap(
828 i, mapArray, signArray, m_forient[i],
829 m_faces[i]->GetXmap()->GetTraceNcoeffs(1),
830 m_faces[i]->GetXmap()->GetTraceNcoeffs(0));
831 }
832
833 for (j = 0; j < m_coordim; j++)
834 {
835 const Array<OneD, const NekDouble> &coeffs =
836 m_faces[i]->GetCoeffs(j);
837
838 for (k = 0; k < nFaceCoeffs; k++)
839 {
840 NekDouble v = signArray[k] * coeffs[k];
841 m_coeffs[j][mapArray[k]] = v;
842 }
843 }
844 }
845
847}
848
849} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
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:50
bool m_setupState
Wether or not the setup routines have been run.
Definition Geometry.h:193
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition Geometry.h:191
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition Geometry.h:704
int GetVid(int i) const
Returns global id of vertex i of this object.
Definition Geometry.h:353
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:372
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
PointGeom * GetVertex(int i) const
Returns vertex i of this object.
Definition Geometry.h:361
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:439
Geometry1D * GetEdge(int i) const
Returns edge i of this object.
Definition Geometry.h:369
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:110
void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces) override
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
static const unsigned int EdgeFaceConnectivity[9][2]
Definition PrismGeom.h:127
int v_GetVertexFaceMap(const int i, const int j) const override
Returns the standard element face IDs that are connected to a given vertex.
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
int v_GetDir(const int faceidx, const int facedir) const override
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition PrismGeom.cpp:80
void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
std::array< PointGeom *, kNverts > m_verts
Definition PrismGeom.h:112
static const unsigned int EdgeNormalToFaceVert[5][4]
Definition PrismGeom.h:128
int v_GetEdgeFaceMap(const int i, const int j) const override
Returns the standard element edge IDs that are connected to a given face.
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.
static const unsigned int VertexFaceConnectivity[6][3]
Definition PrismGeom.h:126
std::array< Geometry2D *, kNfaces > m_faces
Definition PrismGeom.h:114
std::array< StdRegions::Orientation, kNedges > m_eorient
Definition PrismGeom.h:115
std::array< StdRegions::Orientation, kNfaces > m_forient
Definition PrismGeom.h:116
static const unsigned int VertexEdgeConnectivity[6][3]
Definition PrismGeom.h:125
std::array< SegGeom *, kNedges > m_edges
Definition PrismGeom.h:113
int v_GetVertexEdgeMap(const int i, const int j) const override
Returns the standard element edge IDs that are connected to a given vertex.
A simple factory for Xmap objects that is based on the element type, the basis and quadrature selecti...
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48
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.
XmapFactory< StdRegions::StdPrismExp, 3 > & GetStdPrismFactory()
Definition PrismGeom.cpp:55
@ ePtsFilled
Geometric information has been generated.
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290