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