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 
38 #include <StdRegions/StdPrismExp.h>
39 #include <SpatialDomains/SegGeom.h>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace SpatialDomains
47 {
48 
49 const unsigned int PrismGeom::VertexEdgeConnectivity[6][3] = {
50  {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7}, {4, 5, 8}, {6, 7, 8}};
51 const unsigned int PrismGeom::VertexFaceConnectivity[6][3] = {
52  {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4}, {1, 2, 4}, {2, 3, 4}};
53 const unsigned int PrismGeom::EdgeFaceConnectivity[9][2] = {
54  {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2}, {2, 3}, {3, 4}, {2, 4}};
55 
56 PrismGeom::PrismGeom()
57 {
58  m_shapeType = LibUtilities::ePrism;
59 }
60 
61 PrismGeom::PrismGeom(int id, const Geometry2DSharedPtr faces[])
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 
85 int 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  int i, f;
111  GeomType Gtype = eRegular;
112 
113  v_FillGeom();
114 
115  // check to see if expansions are linear
116  for (i = 0; i < m_coordim; ++i)
117  {
118  if (m_xmap->GetBasisNumModes(0) != 2 ||
119  m_xmap->GetBasisNumModes(1) != 2 ||
120  m_xmap->GetBasisNumModes(2) != 2)
121  {
122  Gtype = eDeformed;
123  }
124  }
125 
126  // check to see if all quadrilateral faces are parallelograms
127  if (Gtype == eRegular)
128  {
129  // Vertex ids of quad faces
130  const unsigned int faceVerts[3][4] = {
131  {0, 1, 2, 3}, {1, 2, 5, 4}, {0, 3, 5, 4}};
132 
133  for (f = 0; f < 3; f++)
134  {
135  // Ensure each face is a parallelogram? Check this.
136  for (i = 0; i < m_coordim; i++)
137  {
138  if (fabs((*m_verts[faceVerts[f][0]])(i) -
139  (*m_verts[faceVerts[f][1]])(i) +
140  (*m_verts[faceVerts[f][2]])(i) -
141  (*m_verts[faceVerts[f][3]])(i)) >
143  {
144  Gtype = eDeformed;
145  break;
146  }
147  }
148 
149  if (Gtype == eDeformed)
150  {
151  break;
152  }
153  }
154  }
155 
157  Gtype, m_coordim, m_xmap, m_coeffs);
159  }
160 }
161 
162 int PrismGeom::v_GetVertexEdgeMap(const int i, const int j) const
163 {
164  return VertexEdgeConnectivity[i][j];
165 }
166 
167 int PrismGeom::v_GetVertexFaceMap(const int i, const int j) const
168 {
169  return VertexFaceConnectivity[i][j];
170 }
171 
172 int PrismGeom::v_GetEdgeFaceMap(const int i, const int j) const
173 {
174  return EdgeFaceConnectivity[i][j];
175 }
176 
178 {
179  // find edge 0
180  int i, j;
181  unsigned int check;
182 
183  SegGeomSharedPtr edge;
184 
185  // First set up the 4 bottom edges
186  int f; // Connected face index
187  for (f = 1; f < 5; f++)
188  {
189  int nEdges = m_faces[f]->GetNumEdges();
190  check = 0;
191  for (i = 0; i < 4; i++)
192  {
193  for (j = 0; j < nEdges; j++)
194  {
195  if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
196  {
197  edge = dynamic_pointer_cast<SegGeom>(
198  (m_faces[0])->GetEdge(i));
199  m_edges.push_back(edge);
200  check++;
201  }
202  }
203  }
204 
205  if (check < 1)
206  {
207  std::ostringstream errstrm;
208  errstrm << "Connected faces do not share an edge. Faces ";
209  errstrm << (m_faces[0])->GetGlobalID() << ", "
210  << (m_faces[f])->GetGlobalID();
211  ASSERTL0(false, errstrm.str());
212  }
213  else if (check > 1)
214  {
215  std::ostringstream errstrm;
216  errstrm << "Connected faces share more than one edge. Faces ";
217  errstrm << (m_faces[0])->GetGlobalID() << ", "
218  << (m_faces[f])->GetGlobalID();
219  ASSERTL0(false, errstrm.str());
220  }
221  }
222 
223  // Then, set up the 4 vertical edges
224  check = 0;
225  for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(4)
226  {
227  for (j = 0; j < 4; j++)
228  {
229  if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
230  {
231  edge = dynamic_pointer_cast<SegGeom>(
232  (m_faces[1])->GetEdge(i));
233  m_edges.push_back(edge);
234  check++;
235  }
236  }
237  }
238  if (check < 1)
239  {
240  std::ostringstream errstrm;
241  errstrm << "Connected faces do not share an edge. Faces ";
242  errstrm << (m_faces[1])->GetGlobalID() << ", "
243  << (m_faces[4])->GetGlobalID();
244  ASSERTL0(false, errstrm.str());
245  }
246  else if (check > 1)
247  {
248  std::ostringstream errstrm;
249  errstrm << "Connected faces share more than one edge. Faces ";
250  errstrm << (m_faces[1])->GetGlobalID() << ", "
251  << (m_faces[4])->GetGlobalID();
252  ASSERTL0(false, errstrm.str());
253  }
254  // Set up vertical edges: face(1) through face(4)
255  for (f = 1; f < 4; f++)
256  {
257  check = 0;
258  for (i = 0; i < m_faces[f]->GetNumEdges(); i++)
259  {
260  for (j = 0; j < m_faces[f + 1]->GetNumEdges(); j++)
261  {
262  if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
263  {
264  edge = dynamic_pointer_cast<SegGeom>(
265  (m_faces[f])->GetEdge(i));
266  m_edges.push_back(edge);
267  check++;
268  }
269  }
270  }
271 
272  if (check < 1)
273  {
274  std::ostringstream errstrm;
275  errstrm << "Connected faces do not share an edge. Faces ";
276  errstrm << (m_faces[f])->GetGlobalID() << ", "
277  << (m_faces[f + 1])->GetGlobalID();
278  ASSERTL0(false, errstrm.str());
279  }
280  else if (check > 1)
281  {
282  std::ostringstream errstrm;
283  errstrm << "Connected faces share more than one edge. Faces ";
284  errstrm << (m_faces[f])->GetGlobalID() << ", "
285  << (m_faces[f + 1])->GetGlobalID();
286  ASSERTL0(false, errstrm.str());
287  }
288  }
289 
290  // Finally, set up the 1 top edge
291  check = 0;
292  for (i = 0; i < 4; i++)
293  {
294  for (j = 0; j < 4; j++)
295  {
296  if ((m_faces[2])->GetEid(i) == (m_faces[4])->GetEid(j))
297  {
298  edge = dynamic_pointer_cast<SegGeom>(
299  (m_faces[2])->GetEdge(i));
300  m_edges.push_back(edge);
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.push_back(m_edges[0]->GetVertex(1));
332  m_verts.push_back(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.push_back(m_edges[0]->GetVertex(0));
338  m_verts.push_back(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.push_back(m_edges[i]->GetVertex(1));
355  }
356  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
357  {
358  m_verts.push_back(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.push_back(m_edges[8]->GetVertex(0));
376  m_verts.push_back(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.push_back(m_edges[8]->GetVertex(1));
382  m_verts.push_back(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) == m_verts[edgeVerts[i][1]]->GetGlobalID())
410  {
412  }
413  else
414  {
415  ASSERTL0(false, "Could not find matching vertex for the edge");
416  }
417  }
418 }
419 
421 {
422  int f, i;
423 
424  // These arrays represent the vector of the A and B
425  // coordinate of the local elemental coordinate system
426  // where A corresponds with the coordinate direction xi_i
427  // with the lowest index i (for that particular face)
428  // Coordinate 'B' then corresponds to the other local
429  // coordinate (i.e. with the highest index)
430  Array<OneD, NekDouble> elementAaxis(m_coordim);
431  Array<OneD, NekDouble> elementBaxis(m_coordim);
432 
433  // These arrays correspond to the local coordinate
434  // system of the face itself (i.e. the Geometry2D)
435  // faceAaxis correspond to the xi_0 axis
436  // faceBaxis correspond to the xi_1 axis
439 
440  // This is the base vertex of the face (i.e. the Geometry2D)
441  // This corresponds to thevertex with local ID 0 of the
442  // Geometry2D
443  unsigned int baseVertex;
444 
445  // The lenght of the vectors above
446  NekDouble elementAaxis_length;
447  NekDouble elementBaxis_length;
448  NekDouble faceAaxis_length;
449  NekDouble faceBaxis_length;
450 
451  // This 2D array holds the local id's of all the vertices
452  // for every face. For every face, they are ordered in such
453  // a way that the implementation below allows a unified approach
454  // for all faces.
455  const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {
456  {0, 1, 2, 3},
457  {0, 1, 4, 0}, // This is triangle requires only three vertices
458  {1, 2, 5, 4},
459  {3, 2, 5, 0}, // This is triangle requires only three vertices
460  {0, 3, 5, 4},
461  };
462 
463  NekDouble dotproduct1 = 0.0;
464  NekDouble dotproduct2 = 0.0;
465 
466  unsigned int orientation;
467 
468  // Loop over all the faces to set up the orientation
469  for (f = 0; f < kNqfaces + kNtfaces; f++)
470  {
471  // initialisation
472  elementAaxis_length = 0.0;
473  elementBaxis_length = 0.0;
474  faceAaxis_length = 0.0;
475  faceBaxis_length = 0.0;
476 
477  dotproduct1 = 0.0;
478  dotproduct2 = 0.0;
479 
480  baseVertex = m_faces[f]->GetVid(0);
481 
482  // We are going to construct the vectors representing the A and B axis
483  // of every face. These vectors will be constructed as a
484  // vector-representation
485  // of the edges of the face. However, for both coordinate directions, we
486  // can
487  // represent the vectors by two different edges. That's why we need to
488  // make sure that
489  // we pick the edge to which the baseVertex of the
490  // Geometry2D-representation of the face
491  // belongs...
492 
493  // Compute the length of edges on a base-face
494  if (f == 1 || f == 3)
495  { // Face is a Triangle
496  if (baseVertex == m_verts[faceVerts[f][0]]->GetVid())
497  {
498  for (i = 0; i < m_coordim; i++)
499  {
500  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
501  (*m_verts[faceVerts[f][0]])[i];
502  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
503  (*m_verts[faceVerts[f][0]])[i];
504  }
505  }
506  else if (baseVertex == m_verts[faceVerts[f][1]]->GetVid())
507  {
508  for (i = 0; i < m_coordim; i++)
509  {
510  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
511  (*m_verts[faceVerts[f][0]])[i];
512  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
513  (*m_verts[faceVerts[f][1]])[i];
514  }
515  }
516  else if (baseVertex == m_verts[faceVerts[f][2]]->GetVid())
517  {
518  for (i = 0; i < m_coordim; i++)
519  {
520  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
521  (*m_verts[faceVerts[f][2]])[i];
522  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
523  (*m_verts[faceVerts[f][0]])[i];
524  }
525  }
526  else
527  {
528  ASSERTL0(false, "Could not find matching vertex for the face");
529  }
530  }
531  else
532  { // Face is a Quad
533  if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
534  {
535  for (i = 0; i < m_coordim; i++)
536  {
537  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
538  (*m_verts[faceVerts[f][0]])[i];
539  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
540  (*m_verts[faceVerts[f][0]])[i];
541  }
542  }
543  else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
544  {
545  for (i = 0; i < m_coordim; i++)
546  {
547  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
548  (*m_verts[faceVerts[f][0]])[i];
549  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
550  (*m_verts[faceVerts[f][1]])[i];
551  }
552  }
553  else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
554  {
555  for (i = 0; i < m_coordim; i++)
556  {
557  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
558  (*m_verts[faceVerts[f][3]])[i];
559  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
560  (*m_verts[faceVerts[f][1]])[i];
561  }
562  }
563  else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
564  {
565  for (i = 0; i < m_coordim; i++)
566  {
567  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
568  (*m_verts[faceVerts[f][3]])[i];
569  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
570  (*m_verts[faceVerts[f][0]])[i];
571  }
572  }
573  else
574  {
575  ASSERTL0(false, "Could not find matching vertex for the face");
576  }
577  }
578  // Now, construct the edge-vectors of the local coordinates of
579  // the Geometry2D-representation of the face
580  for (i = 0; i < m_coordim; i++)
581  {
582  int v = m_faces[f]->GetNumVerts() - 1;
583  faceAaxis[i] =
584  (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
585  faceBaxis[i] =
586  (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
587 
588  elementAaxis_length += pow(elementAaxis[i], 2);
589  elementBaxis_length += pow(elementBaxis[i], 2);
590  faceAaxis_length += pow(faceAaxis[i], 2);
591  faceBaxis_length += pow(faceBaxis[i], 2);
592  }
593 
594  elementAaxis_length = sqrt(elementAaxis_length);
595  elementBaxis_length = sqrt(elementBaxis_length);
596  faceAaxis_length = sqrt(faceAaxis_length);
597  faceBaxis_length = sqrt(faceBaxis_length);
598 
599  // Calculate the inner product of both the A-axis
600  // (i.e. Elemental A axis and face A axis)
601  for (i = 0; i < m_coordim; i++)
602  {
603  dotproduct1 += elementAaxis[i] * faceAaxis[i];
604  }
605 
606  orientation = 0;
607 
608  // if the innerproduct is equal to the (absolute value of the ) products
609  // of the lengths of both vectors, then, the coordinate systems will NOT
610  // be transposed
611  if (fabs(elementAaxis_length * faceAaxis_length - fabs(dotproduct1)) <
613  {
614  // if the inner product is negative, both A-axis point
615  // in reverse direction
616  if (dotproduct1 < 0.0)
617  {
618  orientation += 2;
619  }
620 
621  // calculate the inner product of both B-axis
622  for (i = 0; i < m_coordim; i++)
623  {
624  dotproduct2 += elementBaxis[i] * faceBaxis[i];
625  }
626 
627  ASSERTL1(
628  fabs(fabs(dotproduct2 / elementBaxis_length / faceBaxis_length)
629  - 1.0) < NekConstants::kNekZeroTol,
630  "These vectors should be parallel");
631 
632  // if the inner product is negative, both B-axis point
633  // in reverse direction
634  if (dotproduct2 < 0.0)
635  {
636  orientation++;
637  }
638  }
639  // The coordinate systems are transposed
640  else
641  {
642  orientation = 4;
643 
644  // Calculate the inner product between the elemental A-axis
645  // and the B-axis of the face (which are now the corresponding axis)
646  dotproduct1 = 0.0;
647  for (i = 0; i < m_coordim; i++)
648  {
649  dotproduct1 += elementAaxis[i] * faceBaxis[i];
650  }
651 
652  // check that both these axis are indeed parallel
653  ASSERTL1(
654  fabs(fabs(dotproduct1) / elementAaxis_length / faceBaxis_length
655  - 1.0) < NekConstants::kNekZeroTol,
656  "These vectors should be parallel");
657 
658  // if the result is negative, both axis point in reverse
659  // directions
660  if (dotproduct1 < 0.0)
661  {
662  orientation += 2;
663  }
664 
665  // Do the same for the other two corresponding axis
666  dotproduct2 = 0.0;
667  for (i = 0; i < m_coordim; i++)
668  {
669  dotproduct2 += elementBaxis[i] * faceAaxis[i];
670  }
671 
672  ASSERTL1(
673  fabs(fabs(dotproduct2) / elementBaxis_length / faceAaxis_length
674  - 1.0) < NekConstants::kNekZeroTol,
675  "These vectors should be parallel");
676 
677  if (dotproduct2 < 0.0)
678  {
679  orientation++;
680  }
681  }
682 
683  orientation = orientation + 5;
684 
685  if((f == 1)||(f == 3)) // check triange orientation
686  {
688  "Orientation of triangular face (id = " +
689  boost::lexical_cast<string>(m_faces[f]->GetGlobalID()) +
690  ") is inconsistent with face "+
691  boost::lexical_cast<string>(f) +
692  " of prism element (id = "+
693  boost::lexical_cast<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 
703 void 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  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 = *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 = *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 = *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 = *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 = *max_element(tmp.begin(), tmp.end());
774 
777  order0,
779  const LibUtilities::BasisKey B(
781  order1,
783  const LibUtilities::BasisKey C(
785  order2,
788 
790 }
791 
792 }
793 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:60
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
3D geometry information
Definition: Geometry3D.h:68
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:365
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:84
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:195
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:348
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:658
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:356
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:137
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:355
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:319
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:199
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:203
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:189
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:191
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:426
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:187
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:185
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:145
virtual int v_GetVertexFaceMap(const int i, const int j) const
Returns the standard element face IDs that are connected to a given vertex.
Definition: PrismGeom.cpp:167
static const unsigned int EdgeFaceConnectivity[9][2]
Definition: PrismGeom.h:79
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PrismGeom.cpp:731
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: PrismGeom.cpp:85
static const unsigned int VertexFaceConnectivity[6][3]
Definition: PrismGeom.h:78
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: PrismGeom.cpp:703
virtual int v_GetEdgeFaceMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given face.
Definition: PrismGeom.cpp:172
static const unsigned int VertexEdgeConnectivity[6][3]
Definition: PrismGeom.h:77
virtual int v_GetVertexEdgeMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given vertex.
Definition: PrismGeom.cpp:162
static const int kNverts
Definition: QuadGeom.h:77
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
@ 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:62
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:1
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267