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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Pyramidic geometry information.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <SpatialDomains/PyrGeom.h>
37 
40 #include <StdRegions/StdPyrExp.h>
41 #include <SpatialDomains/SegGeom.h>
42 #include <SpatialDomains/TriGeom.h>
44 
45 
46 namespace Nektar
47 {
48  namespace SpatialDomains
49  {
51  {
53  }
54 
56  Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
57  {
59 
60  /// Copy the face shared pointers
61  m_faces.insert(m_faces.begin(), faces, faces+PyrGeom::kNfaces);
62 
63  /// Set up orientation vectors with correct amount of elements.
64  m_eorient.resize(kNedges);
65  m_forient.resize(kNfaces);
66 
71  SetUpXmap();
72  SetUpCoeffs(m_xmap->GetNcoeffs());
73  }
74 
76  {
77 
78  }
79 
81  {
83  {
84  int i;
85  GeomType Gtype = eRegular;
86 
87  v_FillGeom();
88 
89  // check to see if expansions are linear
90  for(i = 0; i < m_coordim; ++i)
91  {
92  if (m_xmap->GetBasisNumModes(0) != 2 ||
93  m_xmap->GetBasisNumModes(1) != 2 ||
94  m_xmap->GetBasisNumModes(2) != 2 )
95  {
96  Gtype = eDeformed;
97  }
98  }
99 
100  // check to see if all quadrilateral faces are parallelograms
101  if(Gtype == eRegular)
102  {
103  // Ensure each face is a parallelogram? Check this.
104  for (i = 0; i < m_coordim; i++)
105  {
106  if( fabs( (*m_verts[0])(i) -
107  (*m_verts[1])(i) +
108  (*m_verts[2])(i) -
109  (*m_verts[3])(i) )
111  {
112  Gtype = eDeformed;
113  break;
114  }
115  }
116  }
117 
119  Gtype, m_coordim, m_xmap, m_coeffs);
121  }
122  }
123 
125  const Array<OneD, const NekDouble> &coords,
126  Array<OneD, NekDouble> &Lcoords)
127  {
128  NekDouble ptdist = 1e6;
129 
130  v_FillGeom();
131 
132  // calculate local coordinate for coord
133  if(GetMetricInfo()->GetGtype() == eRegular)
134  { // Based on Spen's book, page 99
135 
136  // Point inside tetrahedron
137  PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
138 
139  // Edges
140  PointGeom er0, e10, e30, e40;
141  er0.Sub(r,*m_verts[0]);
142  e10.Sub(*m_verts[1],*m_verts[0]);
143  e30.Sub(*m_verts[3],*m_verts[0]);
144  e40.Sub(*m_verts[4],*m_verts[0]);
145 
146 
147  // Cross products (Normal times area)
148  PointGeom cp1030, cp3040, cp4010;
149  cp1030.Mult(e10,e30);
150  cp3040.Mult(e30,e40);
151  cp4010.Mult(e40,e10);
152 
153 
154  // Barycentric coordinates (relative volume)
155  NekDouble V = e40.dot(cp1030); // Pyramid Volume = {(e40)dot(e10)x(e30)}/4
156  NekDouble scaleFactor = 2.0/3.0;
157  NekDouble v1 = er0.dot(cp3040) / V; // volume1 = {(er0)dot(e30)x(e40)}/6
158  NekDouble v2 = er0.dot(cp4010) / V; // volume2 = {(er0)dot(e40)x(e10)}/6
159  NekDouble beta = v1 * scaleFactor;
160  NekDouble gamma = v2 * scaleFactor;
161  NekDouble delta = er0.dot(cp1030) / V; // volume3 = {(er0)dot(e10)x(e30)}/4
162 
163  // Make Pyramid bigger
164  Lcoords[0] = 2.0*beta - 1.0;
165  Lcoords[1] = 2.0*gamma - 1.0;
166  Lcoords[2] = 2.0*delta - 1.0;
167 
168  // Set ptdist to distance to nearest vertex
169  for(int i = 0; i < 5; ++i)
170  {
171  ptdist = min(ptdist,r.dist(*m_verts[i]));
172  }
173 
174  }
175  else
176  {
178  "inverse mapping must be set up to use this call");
179  }
180  return ptdist;
181  }
182 
184  {
185  return 5;
186  }
187 
189  {
190  return 8;
191  }
192 
193  int PyrGeom::v_GetDir(const int faceidx, const int facedir) const
194  {
195  if (faceidx == 0)
196  {
197  return facedir;
198  }
199  else if (faceidx == 1 || faceidx == 3)
200  {
201  return 2 * facedir;
202  }
203  else
204  {
205  return 1 + facedir;
206  }
207  }
208 
210  {
211  // find edge 0
212  int i,j;
213  unsigned int check;
214 
215  SegGeomSharedPtr edge;
216 
217  // First set up the 4 bottom edges
218  int f;
219  for (f = 1; f < 5; f++)
220  {
221  int nEdges = m_faces[f]->GetNumEdges();
222  check = 0;
223  for (i = 0; i < 4; i++)
224  {
225  for (j = 0; j < nEdges; j++)
226  {
227  if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
228  {
229  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[0])->GetEdge(i));
230  m_edges.push_back(edge);
231  check++;
232  }
233  }
234  }
235 
236  if (check < 1)
237  {
238  std::ostringstream errstrm;
239  errstrm << "Connected faces do not share an edge. Faces ";
240  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
241  ASSERTL0(false, errstrm.str());
242  }
243  else if (check > 1)
244  {
245  std::ostringstream errstrm;
246  errstrm << "Connected faces share more than one edge. Faces ";
247  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
248  ASSERTL0(false, errstrm.str());
249  }
250  }
251 
252  // Then, set up the 4 vertical edges
253  check = 0;
254  for(i = 0; i < 3; i++) //Set up the vertical edge :face(1) and face(4)
255  {
256  for(j = 0; j < 3; j++)
257  {
258  if( (m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j) )
259  {
260  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[1])->GetEdge(i));
261  m_edges.push_back(edge);
262  check++;
263  }
264  }
265  }
266  if( check < 1 )
267  {
268  std::ostringstream errstrm;
269  errstrm << "Connected faces do not share an edge. Faces ";
270  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
271  ASSERTL0(false, errstrm.str());
272  }
273  else if( check > 1)
274  {
275  std::ostringstream errstrm;
276  errstrm << "Connected faces share more than one edge. Faces ";
277  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
278  ASSERTL0(false, errstrm.str());
279  }
280 
281  // Set up vertical edges: face(1) through face(4)
282  for (f = 1; f < 4; f++)
283  {
284  check = 0;
285  for(i = 0; i < m_faces[f]->GetNumEdges(); i++)
286  {
287  for(j = 0; j < m_faces[f+1]->GetNumEdges(); j++)
288  {
289  if( (m_faces[f])->GetEid(i) == (m_faces[f+1])->GetEid(j))
290  {
291  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[f])->GetEdge(i));
292  m_edges.push_back(edge);
293  check++;
294  }
295  }
296  }
297 
298  if( check < 1 )
299  {
300  std::ostringstream errstrm;
301  errstrm << "Connected faces do not share an edge. Faces ";
302  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
303  ASSERTL0(false, errstrm.str());
304  }
305  else if( check > 1)
306  {
307  std::ostringstream errstrm;
308  errstrm << "Connected faces share more than one edge. Faces ";
309  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
310  ASSERTL0(false, errstrm.str());
311  }
312  }
313  }
314 
316  {
317  // Set up the first 2 vertices (i.e. vertex 0,1)
318  if (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0) ||
319  m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1))
320  {
321  m_verts.push_back(m_edges[0]->GetVertex(1));
322  m_verts.push_back(m_edges[0]->GetVertex(0));
323  }
324  else if (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0) ||
325  m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1))
326  {
327  m_verts.push_back(m_edges[0]->GetVertex(0));
328  m_verts.push_back(m_edges[0]->GetVertex(1));
329  }
330  else
331  {
332  std::ostringstream errstrm;
333  errstrm << "Connected edges do not share a vertex. Edges ";
334  errstrm << m_edges[0]->GetEid() << ", " << m_edges[1]->GetEid();
335  ASSERTL0(false, errstrm.str());
336  }
337 
338  // set up the other bottom vertices (i.e. vertex 2,3)
339  for(int i = 1; i < 3; i++)
340  {
341  if (m_edges[i]->GetVid(0) == m_verts[i]->GetVid())
342  {
343  m_verts.push_back(m_edges[i]->GetVertex(1));
344  }
345  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetVid())
346  {
347  m_verts.push_back(m_edges[i]->GetVertex(0));
348  }
349  else
350  {
351  std::ostringstream errstrm;
352  errstrm << "Connected edges do not share a vertex. Edges ";
353  errstrm << m_edges[i]->GetEid() << ", " << m_edges[i-1]->GetEid();
354  ASSERTL0(false, errstrm.str());
355  }
356  }
357 
358  // set up top vertex
359  if (m_edges[4]->GetVid(0) == m_verts[0]->GetVid())
360  {
361  m_verts.push_back(m_edges[4]->GetVertex(1));
362  }
363  else
364  {
365  m_verts.push_back(m_edges[4]->GetVertex(0));
366  }
367 
368  int check = 0;
369  for (int i = 5; i < 8; ++i)
370  {
371  if( (m_edges[i]->GetVid(0) == m_verts[i-4]->GetVid()
372  && m_edges[i]->GetVid(1) == m_verts[4]->GetVid())
373  ||(m_edges[i]->GetVid(1) == m_verts[i-4]->GetVid()
374  && m_edges[i]->GetVid(0) == m_verts[4]->GetVid()))
375  {
376  check++;
377  }
378  }
379  if (check != 3) {
380  std::ostringstream errstrm;
381  errstrm << "Connected edges do not share a vertex. Edges ";
382  errstrm << m_edges[3]->GetEid() << ", " << m_edges[2]->GetEid();
383  ASSERTL0(false, errstrm.str());
384  }
385  }
386 
388  {
389  // This 2D array holds the local id's of all the vertices for every
390  // edge. For every edge, they are ordered to what we define as being
391  // Forwards.
392  const unsigned int edgeVerts[kNedges][2] =
393  { {0,1}, {1,2}, {3,2}, {0,3}, {0,4}, {1,4}, {2,4}, {3,4} };
394 
395  int i;
396  for (i = 0; i < kNedges; i++)
397  {
398  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetVid())
399  {
401  }
402  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetVid())
403  {
405  }
406  else
407  {
408  ASSERTL0(false,"Could not find matching vertex for the edge");
409  }
410  }
411  }
412 
414  {
415  int f,i;
416 
417  // These arrays represent the vector of the A and B coordinate of
418  // the local elemental coordinate system where A corresponds with
419  // the coordinate direction xi_i with the lowest index i (for that
420  // particular face) Coordinate 'B' then corresponds to the other
421  // local coordinate (i.e. with the highest index)
422  Array<OneD,NekDouble> elementAaxis(m_coordim);
423  Array<OneD,NekDouble> elementBaxis(m_coordim);
424 
425  // These arrays correspond to the local coordinate
426  // system of the face itself (i.e. the Geometry2D)
427  // faceAaxis correspond to the xi_0 axis
428  // faceBaxis correspond to the xi_1 axis
429  Array<OneD,NekDouble> faceAaxis(m_coordim);
430  Array<OneD,NekDouble> faceBaxis(m_coordim);
431 
432  // This is the base vertex of the face (i.e. the Geometry2D) This
433  // corresponds to thevertex with local ID 0 of the Geometry2D
434  unsigned int baseVertex;
435 
436  // The lenght of the vectors above
437  NekDouble elementAaxis_length;
438  NekDouble elementBaxis_length;
439  NekDouble faceAaxis_length;
440  NekDouble faceBaxis_length;
441 
442  // This 2D array holds the local id's of all the vertices for every
443  // face. For every face, they are ordered in such a way that the
444  // implementation below allows a unified approach for all faces.
445  const unsigned int faceVerts[kNfaces][4] = {
446  {0,1,2,3},
447  {0,1,4,0}, // Last four elements are triangles which only
448  {1,2,4,0}, // require three vertices.
449  {3,2,4,0},
450  {0,3,4,0}
451  };
452 
453  NekDouble dotproduct1 = 0.0;
454  NekDouble dotproduct2 = 0.0;
455 
456  unsigned int orientation;
457 
458  // Loop over all the faces to set up the orientation
459  for(f = 0; f < kNqfaces + kNtfaces; f++)
460  {
461  // initialisation
462  elementAaxis_length = 0.0;
463  elementBaxis_length = 0.0;
464  faceAaxis_length = 0.0;
465  faceBaxis_length = 0.0;
466 
467  dotproduct1 = 0.0;
468  dotproduct2 = 0.0;
469 
470  baseVertex = m_faces[f]->GetVid(0);
471 
472  // We are going to construct the vectors representing the A and
473  // B axis of every face. These vectors will be constructed as a
474  // vector-representation of the edges of the face. However, for
475  // both coordinate directions, we can represent the vectors by
476  // two different edges. That's why we need to make sure that we
477  // pick the edge to which the baseVertex of the
478  // Geometry2D-representation of the face belongs...
479 
480  // Compute the length of edges on a base-face
481  if (f > 0)
482  {
483  for(i = 0; i < m_coordim; i++)
484  {
485  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
486  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
487  }
488  }
489  else
490  {
491  if( baseVertex == m_verts[ faceVerts[f][0] ]->GetVid() )
492  {
493  for(i = 0; i < m_coordim; i++)
494  {
495  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
496  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
497  }
498  }
499  else if( baseVertex == m_verts[ faceVerts[f][1] ]->GetVid() )
500  {
501  for(i = 0; i < m_coordim; i++)
502  {
503  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
504  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
505  }
506  }
507  else if( baseVertex == m_verts[ faceVerts[f][2] ]->GetVid() )
508  {
509  for(i = 0; i < m_coordim; i++)
510  {
511  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
512  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
513  }
514  }
515  else if( baseVertex == m_verts[ faceVerts[f][3] ]->GetVid() )
516  {
517  for(i = 0; i < m_coordim; i++)
518  {
519  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
520  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
521  }
522  }
523  else
524  {
525  ASSERTL0(false, "Could not find matching vertex for the face");
526  }
527  }
528 
529  // Now, construct the edge-vectors of the local coordinates of
530  // the Geometry2D-representation of the face
531  for(i = 0; i < m_coordim; i++)
532  {
533  int v = m_faces[f]->GetNumVerts()-1;
534  faceAaxis[i] = (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
535  faceBaxis[i] = (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
536 
537  elementAaxis_length += pow(elementAaxis[i],2);
538  elementBaxis_length += pow(elementBaxis[i],2);
539  faceAaxis_length += pow(faceAaxis[i],2);
540  faceBaxis_length += pow(faceBaxis[i],2);
541  }
542 
543  elementAaxis_length = sqrt(elementAaxis_length);
544  elementBaxis_length = sqrt(elementBaxis_length);
545  faceAaxis_length = sqrt(faceAaxis_length);
546  faceBaxis_length = sqrt(faceBaxis_length);
547 
548  // Calculate the inner product of both the A-axis
549  // (i.e. Elemental A axis and face A axis)
550  for(i = 0 ; i < m_coordim; i++)
551  {
552  dotproduct1 += elementAaxis[i]*faceAaxis[i];
553  }
554 
555  orientation = 0;
556  // if the innerproduct is equal to the (absolute value of the ) products of the lengths
557  // of both vectors, then, the coordinate systems will NOT be transposed
558  if( fabs(elementAaxis_length*faceAaxis_length - fabs(dotproduct1)) < NekConstants::kNekZeroTol )
559  {
560  // if the inner product is negative, both A-axis point
561  // in reverse direction
562  if(dotproduct1 < 0.0)
563  {
564  orientation += 2;
565  }
566 
567  // calculate the inner product of both B-axis
568  for(i = 0 ; i < m_coordim; i++)
569  {
570  dotproduct2 += elementBaxis[i]*faceBaxis[i];
571  }
572 
573  // if the inner product is negative, both B-axis point
574  // in reverse direction
575  if( dotproduct2 < 0.0 )
576  {
577  orientation++;
578  }
579  }
580  // The coordinate systems are transposed
581  else
582  {
583  orientation = 4;
584 
585  // Calculate the inner product between the elemental A-axis
586  // and the B-axis of the face (which are now the corresponding axis)
587  dotproduct1 = 0.0;
588  for(i = 0 ; i < m_coordim; i++)
589  {
590  dotproduct1 += elementAaxis[i]*faceBaxis[i];
591  }
592 
593  // check that both these axis are indeed parallel
594  if (fabs(elementAaxis_length*faceBaxis_length
595  - fabs(dotproduct1)) > NekConstants::kNekZeroTol)
596  {
597  cout << "Warning: Prism axes not parallel" << endl;
598  }
599 
600  // if the result is negative, both axis point in reverse
601  // directions
602  if(dotproduct1 < 0.0)
603  {
604  orientation += 2;
605  }
606 
607  // Do the same for the other two corresponding axis
608  dotproduct2 = 0.0;
609  for(i = 0 ; i < m_coordim; i++)
610  {
611  dotproduct2 += elementBaxis[i]*faceAaxis[i];
612  }
613 
614  // check that both these axis are indeed parallel
615  if (fabs(elementBaxis_length*faceAaxis_length
616  - fabs(dotproduct2)) > NekConstants::kNekZeroTol)
617  {
618  cout << "Warning: Prism axes not parallel" << endl;
619  }
620 
621  if( dotproduct2 < 0.0 )
622  {
623  orientation++;
624  }
625  }
626 
627  orientation = orientation + 5;
628 
629  // Fill the m_forient array
630  m_forient[f] = (StdRegions::Orientation) orientation;
631  }
632  }
633 
635  CurveMap &curvedEdges,
636  CurveMap &curvedFaces)
637  {
638  Geometry::v_Reset(curvedEdges, curvedFaces);
639 
640  for (int i = 0; i < 5; ++i)
641  {
642  m_faces[i]->Reset(curvedEdges, curvedFaces);
643  }
644 
645  SetUpXmap();
646  SetUpCoeffs(m_xmap->GetNcoeffs());
647  }
648 
649  /**
650  * @brief Set up the #m_xmap object by determining the order of each
651  * direction from derived faces.
652  */
654  {
655  vector<int> tmp;
656  int order0, points0, order1, points1;
657 
658  if (m_forient[0] < 9)
659  {
660  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
661  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
662  order0 = *max_element(tmp.begin(), tmp.end());
663 
664  tmp.clear();
665  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(0));
666  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(2));
667  points0 = *max_element(tmp.begin(), tmp.end());
668  }
669  else
670  {
671  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
672  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
673  order0 = *max_element(tmp.begin(), tmp.end());
674 
675  tmp.clear();
676  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(1));
677  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(3));
678  points0 = *max_element(tmp.begin(), tmp.end());
679  }
680 
681  if (m_forient[0] < 9)
682  {
683  tmp.clear();
684  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
685  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
686  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
687  order1 = *max_element(tmp.begin(), tmp.end());
688 
689  tmp.clear();
690  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(1));
691  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(3));
692  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNumPoints(2));
693  points1 = *max_element(tmp.begin(), tmp.end());
694  }
695  else
696  {
697  tmp.clear();
698  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
699  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
700  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
701  order1 = *max_element(tmp.begin(), tmp.end());
702 
703  tmp.clear();
704  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(0));
705  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNumPoints(2));
706  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNumPoints(2));
707  points1 = *max_element(tmp.begin(), tmp.end());
708  }
709 
710  tmp.clear();
711  tmp.push_back(order0);
712  tmp.push_back(order1);
713  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(1));
714  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(2));
715  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(1));
716  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(2));
717  int order2 = *max_element(tmp.begin(), tmp.end());
718 
719  tmp.clear();
720  tmp.push_back(points0);
721  tmp.push_back(points1);
722  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNumPoints(1));
723  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNumPoints(2));
724  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNumPoints(1));
725  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNumPoints(2));
726  tmp.push_back(m_faces[1]->GetEdge(1)->GetBasis(0)->GetNumPoints());
727  tmp.push_back(m_faces[1]->GetEdge(2)->GetBasis(0)->GetNumPoints());
728  tmp.push_back(m_faces[3]->GetEdge(1)->GetBasis(0)->GetNumPoints());
729  tmp.push_back(m_faces[3]->GetEdge(2)->GetBasis(0)->GetNumPoints());
730  int points2 = *max_element(tmp.begin(), tmp.end());
731 
732  const LibUtilities::BasisKey A1(
736  const LibUtilities::BasisKey A2(
740  const LibUtilities::BasisKey C(
744 
746  A1, A2, C);
747  }
748  }; //end of namespace
749 }; //end of namespace
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:92
Principle Modified Functions .
Definition: BasisType.h:51
virtual int v_GetDir(const int faceidx, const int facedir) const
Definition: PyrGeom.cpp:193
static const int kNfaces
Definition: PyrGeom.h:60
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
GeomFactorsSharedPtr m_geomFactors
Definition: Geometry.h:170
NekDouble dot(PointGeom &a)
Definition: PointGeom.cpp:192
virtual int v_GetNumEdges() const
Definition: PyrGeom.cpp:188
int GetEid(int i) const
Return the ID of edge i in this element.
Definition: Geometry3D.cpp:71
Principle Modified Functions .
Definition: BasisType.h:49
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:168
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: Geometry.cpp:307
int GetFid(int i) const
Definition: Geometry.h:329
StdRegions::StdExpansionSharedPtr GetXmap() const
Definition: Geometry.h:383
virtual int v_GetNumVerts() const
Definition: PyrGeom.cpp:183
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:93
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: PyrGeom.cpp:634
const LibUtilities::BasisSharedPtr GetBasis(const int i)
Return the j-th basis of the i-th co-ordinate dimension.
Definition: Geometry.h:475
static const NekDouble kNekZeroTol
static const int kNedges
Definition: PyrGeom.h:57
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
Defines a specification for a set of points.
Definition: Points.h:58
double NekDouble
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: PyrGeom.cpp:124
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
static const int kNtfaces
Definition: PyrGeom.h:59
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:244
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
void Mult(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:177
const Geometry1DSharedPtr GetEdge(int i) const
static const int kNqfaces
Definition: PyrGeom.h:58
3D geometry information
Definition: Geometry3D.h:70
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometry is straight-sided with constant geometric factors.
int GetVid(int i) const
Definition: Geometry.h:319
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
NekDouble dist(PointGeom &a)
Definition: PointGeom.cpp:186
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
GeomFactorsSharedPtr GetMetricInfo()
Definition: Geometry.h:299
GeomType
Indicates the type of element geometry.
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PyrGeom.cpp:653
virtual void v_GenGeomFactors()
Definition: PyrGeom.cpp:80
Geometry is curved or has non-constant factors.
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Describes the specification for a Basis.
Definition: Basis.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50