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