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