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