Nektar++
PyrGeom.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PyrGeom.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Pyramidic geometry information.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <SpatialDomains/PyrGeom.h>
36 
39 #include <StdRegions/StdPyrExp.h>
40 #include <SpatialDomains/SegGeom.h>
41 #include <SpatialDomains/TriGeom.h>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47 namespace SpatialDomains
48 {
49 
50 PyrGeom::PyrGeom()
51 {
52  m_shapeType = LibUtilities::ePyramid;
53 }
54 
55 PyrGeom::PyrGeom(int id, const Geometry2DSharedPtr faces[])
56  : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
57 {
59  m_globalID = id;
60 
61  /// Copy the face shared pointers
62  m_faces.insert(m_faces.begin(), faces, faces + PyrGeom::kNfaces);
63 
64  /// Set up orientation vectors with correct amount of elements.
65  m_eorient.resize(kNedges);
66  m_forient.resize(kNfaces);
67 
72 }
73 
75 {
76 }
77 
78 
80  Array<OneD, NekDouble> &locCoord,
81  NekDouble tol,
82  NekDouble &resid)
83 {
84  //Rough check if within twice min/max point
85  if (GetMetricInfo()->GetGtype() != eRegular)
86  {
87  if (!MinMaxCheck(gloCoord))
88  {
89  return false;
90  }
91  }
92 
93  // Convert to the local Cartesian coordinates.
94  resid = GetLocCoords(gloCoord, locCoord);
95 
96  // Check local coordinate is within std region bounds.
97  if (locCoord[0] >= -(1 + tol) && locCoord[1] >= -(1 + tol) &&
98  locCoord[2] >= -(1 + tol) && locCoord[0] + locCoord[2] <= tol &&
99  locCoord[1] + locCoord[2] <= tol)
100  {
101  return true;
102  }
103 
104  //Clamp local coords
105  ClampLocCoords(locCoord, tol);
106 
107  return false;
108 }
109 
111 {
112  if(!m_setupState)
113  {
115  }
116 
118  {
119  int i;
120  GeomType Gtype = eRegular;
121 
122  v_FillGeom();
123 
124  // check to see if expansions are linear
125  for (i = 0; i < m_coordim; ++i)
126  {
127  if (m_xmap->GetBasisNumModes(0) != 2 ||
128  m_xmap->GetBasisNumModes(1) != 2 ||
129  m_xmap->GetBasisNumModes(2) != 2)
130  {
131  Gtype = eDeformed;
132  }
133  }
134 
135  // check to see if all quadrilateral faces are parallelograms
136  if (Gtype == eRegular)
137  {
138  // Ensure each face is a parallelogram? Check this.
139  for (i = 0; i < m_coordim; i++)
140  {
141  if (fabs((*m_verts[0])(i) - (*m_verts[1])(i) +
142  (*m_verts[2])(i) - (*m_verts[3])(i)) >
144  {
145  Gtype = eDeformed;
146  break;
147  }
148  }
149  }
150 
152  Gtype, m_coordim, m_xmap, m_coeffs);
154  }
155 }
156 
158  Array<OneD, NekDouble> &Lcoords)
159 {
160  NekDouble ptdist = 1e6;
161 
162  v_FillGeom();
163 
164  // calculate local coordinate for coord
165  if (GetMetricInfo()->GetGtype() == eRegular&&0) // This method does not currently work and so is disabled
166  { // Based on Spen's book, page 99
167 
168  // Point inside tetrahedron
169  PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
170 
171  // Edges
172  PointGeom er0, e10, e30, e40;
173  er0.Sub(r, *m_verts[0]);
174  e10.Sub(*m_verts[1], *m_verts[0]);
175  e30.Sub(*m_verts[3], *m_verts[0]);
176  e40.Sub(*m_verts[4], *m_verts[0]);
177 
178  // Cross products (Normal times area)
179  PointGeom cp1030, cp3040, cp4010;
180  cp1030.Mult(e10, e30);
181  cp3040.Mult(e30, e40);
182  cp4010.Mult(e40, e10);
183 
184  // Barycentric coordinates (relative volume)
185  NekDouble V =
186  e40.dot(cp1030); // Pyramid Volume = {(e40)dot(e10)x(e30)}/4
187  NekDouble scaleFactor = 2.0 / 3.0;
188  NekDouble v1 = er0.dot(cp3040) / V; // volume1 = {(er0)dot(e30)x(e40)}/6
189  NekDouble v2 = er0.dot(cp4010) / V; // volume2 = {(er0)dot(e40)x(e10)}/6
190  NekDouble beta = v1 * scaleFactor;
191  NekDouble gamma = v2 * scaleFactor;
192  NekDouble delta =
193  er0.dot(cp1030) / V; // volume3 = {(er0)dot(e10)x(e30)}/4
194 
195  // Make Pyramid bigger
196  Lcoords[0] = 2.0 * beta - 1.0;
197  Lcoords[1] = 2.0 * gamma - 1.0;
198  Lcoords[2] = 2.0 * delta - 1.0;
199 
200  // Set ptdist to distance to nearest vertex
201  for (int i = 0; i < 5; ++i)
202  {
203  ptdist = min(ptdist, r.dist(*m_verts[i]));
204  }
205  }
206  else
207  {
208 
209  v_FillGeom();
210 
211  // Determine nearest point of coords to values in m_xmap
212  int npts = m_xmap->GetTotPoints();
213  Array<OneD, NekDouble> ptsx(npts), ptsy(npts), ptsz(npts);
214  Array<OneD, NekDouble> tmp1(npts), tmp2(npts);
215 
216  m_xmap->BwdTrans(m_coeffs[0], ptsx);
217  m_xmap->BwdTrans(m_coeffs[1], ptsy);
218  m_xmap->BwdTrans(m_coeffs[2], ptsz);
219 
220  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
221  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
222  const Array<OneD, const NekDouble> zc = m_xmap->GetPoints(2);
223 
224  // guess the first local coords based on nearest point
225  Vmath::Sadd(npts, -coords[0], ptsx, 1, tmp1, 1);
226  Vmath::Vmul(npts, tmp1, 1, tmp1, 1, tmp1, 1);
227  Vmath::Sadd(npts, -coords[1], ptsy, 1, tmp2, 1);
228  Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
229  Vmath::Sadd(npts, -coords[2], ptsz, 1, tmp2, 1);
230  Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
231 
232  int min_i = Vmath::Imin(npts, tmp1, 1);
233 
234  // distance from coordinate to nearest point for return value.
235  ptdist = sqrt(tmp1[min_i]);
236 
237  // Get collapsed coordinate
238  int qa = za.num_elements(), qb = zb.num_elements();
239  Lcoords[2] = zc[min_i / (qa * qb)];
240  min_i = min_i % (qa * qb);
241  Lcoords[1] = zb[min_i / qa];
242  Lcoords[0] = za[min_i % qa];
243 
244  // recover cartesian coordinate from collapsed coordinate.
245  Lcoords[0] = (1.0 + Lcoords[0]) * (1.0 - Lcoords[2]) / 2 - 1.0;
246  Lcoords[1] = (1.0 + Lcoords[1]) * (1.0 - Lcoords[2]) / 2 - 1.0;
247 
248  // Perform newton iteration to find local coordinates
249  NekDouble resid = 0.0;
250  NewtonIterationForLocCoord(coords, ptsx, ptsy, ptsz, Lcoords, resid);
251  }
252 
253  return ptdist;
254 }
255 
256 int PyrGeom::v_GetDir(const int faceidx, const int facedir) const
257 {
258  if (faceidx == 0)
259  {
260  return facedir;
261  }
262  else if (faceidx == 1 || faceidx == 3)
263  {
264  return 2 * facedir;
265  }
266  else
267  {
268  return 1 + facedir;
269  }
270 }
271 
273 {
274  // find edge 0
275  int i, j;
276  unsigned int check;
277 
278  SegGeomSharedPtr edge;
279 
280  // First set up the 4 bottom edges
281  int f;
282  for (f = 1; f < 5; f++)
283  {
284  int nEdges = m_faces[f]->GetNumEdges();
285  check = 0;
286  for (i = 0; i < 4; i++)
287  {
288  for (j = 0; j < nEdges; j++)
289  {
290  if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
291  {
292  edge = dynamic_pointer_cast<SegGeom>(
293  (m_faces[0])->GetEdge(i));
294  m_edges.push_back(edge);
295  check++;
296  }
297  }
298  }
299 
300  if (check < 1)
301  {
302  std::ostringstream errstrm;
303  errstrm << "Connected faces do not share an edge. Faces ";
304  errstrm << (m_faces[0])->GetGlobalID() << ", "
305  << (m_faces[f])->GetGlobalID();
306  ASSERTL0(false, errstrm.str());
307  }
308  else if (check > 1)
309  {
310  std::ostringstream errstrm;
311  errstrm << "Connected faces share more than one edge. Faces ";
312  errstrm << (m_faces[0])->GetGlobalID() << ", "
313  << (m_faces[f])->GetGlobalID();
314  ASSERTL0(false, errstrm.str());
315  }
316  }
317 
318  // Then, set up the 4 vertical edges
319  check = 0;
320  for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(4)
321  {
322  for (j = 0; j < 3; j++)
323  {
324  if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
325  {
326  edge = dynamic_pointer_cast<SegGeom>(
327  (m_faces[1])->GetEdge(i));
328  m_edges.push_back(edge);
329  check++;
330  }
331  }
332  }
333  if (check < 1)
334  {
335  std::ostringstream errstrm;
336  errstrm << "Connected faces do not share an edge. Faces ";
337  errstrm << (m_faces[1])->GetGlobalID() << ", "
338  << (m_faces[4])->GetGlobalID();
339  ASSERTL0(false, errstrm.str());
340  }
341  else if (check > 1)
342  {
343  std::ostringstream errstrm;
344  errstrm << "Connected faces share more than one edge. Faces ";
345  errstrm << (m_faces[1])->GetGlobalID() << ", "
346  << (m_faces[4])->GetGlobalID();
347  ASSERTL0(false, errstrm.str());
348  }
349 
350  // Set up vertical edges: face(1) through face(4)
351  for (f = 1; f < 4; f++)
352  {
353  check = 0;
354  for (i = 0; i < m_faces[f]->GetNumEdges(); i++)
355  {
356  for (j = 0; j < m_faces[f + 1]->GetNumEdges(); j++)
357  {
358  if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
359  {
360  edge = dynamic_pointer_cast<SegGeom>(
361  (m_faces[f])->GetEdge(i));
362  m_edges.push_back(edge);
363  check++;
364  }
365  }
366  }
367 
368  if (check < 1)
369  {
370  std::ostringstream errstrm;
371  errstrm << "Connected faces do not share an edge. Faces ";
372  errstrm << (m_faces[f])->GetGlobalID() << ", "
373  << (m_faces[f + 1])->GetGlobalID();
374  ASSERTL0(false, errstrm.str());
375  }
376  else if (check > 1)
377  {
378  std::ostringstream errstrm;
379  errstrm << "Connected faces share more than one edge. Faces ";
380  errstrm << (m_faces[f])->GetGlobalID() << ", "
381  << (m_faces[f + 1])->GetGlobalID();
382  ASSERTL0(false, errstrm.str());
383  }
384  }
385 }
386 
388 {
389  // Set up the first 2 vertices (i.e. vertex 0,1)
390  if (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0) ||
391  m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1))
392  {
393  m_verts.push_back(m_edges[0]->GetVertex(1));
394  m_verts.push_back(m_edges[0]->GetVertex(0));
395  }
396  else if (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0) ||
397  m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1))
398  {
399  m_verts.push_back(m_edges[0]->GetVertex(0));
400  m_verts.push_back(m_edges[0]->GetVertex(1));
401  }
402  else
403  {
404  std::ostringstream errstrm;
405  errstrm << "Connected edges do not share a vertex. Edges ";
406  errstrm << m_edges[0]->GetGlobalID() << ", "
407  << m_edges[1]->GetGlobalID();
408  ASSERTL0(false, errstrm.str());
409  }
410 
411  // set up the other bottom vertices (i.e. vertex 2,3)
412  for (int i = 1; i < 3; i++)
413  {
414  if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
415  {
416  m_verts.push_back(m_edges[i]->GetVertex(1));
417  }
418  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
419  {
420  m_verts.push_back(m_edges[i]->GetVertex(0));
421  }
422  else
423  {
424  std::ostringstream errstrm;
425  errstrm << "Connected edges do not share a vertex. Edges ";
426  errstrm << m_edges[i]->GetGlobalID() << ", "
427  << m_edges[i - 1]->GetGlobalID();
428  ASSERTL0(false, errstrm.str());
429  }
430  }
431 
432  // set up top vertex
433  if (m_edges[4]->GetVid(0) == m_verts[0]->GetGlobalID())
434  {
435  m_verts.push_back(m_edges[4]->GetVertex(1));
436  }
437  else
438  {
439  m_verts.push_back(m_edges[4]->GetVertex(0));
440  }
441 
442  int check = 0;
443  for (int i = 5; i < 8; ++i)
444  {
445  if ((m_edges[i]->GetVid(0) == m_verts[i - 4]->GetGlobalID() &&
446  m_edges[i]->GetVid(1) == m_verts[4]->GetGlobalID()) ||
447  (m_edges[i]->GetVid(1) == m_verts[i - 4]->GetGlobalID() &&
448  m_edges[i]->GetVid(0) == m_verts[4]->GetGlobalID()))
449  {
450  check++;
451  }
452  }
453  if (check != 3)
454  {
455  std::ostringstream errstrm;
456  errstrm << "Connected edges do not share a vertex. Edges ";
457  errstrm << m_edges[3]->GetGlobalID() << ", "
458  << m_edges[2]->GetGlobalID();
459  ASSERTL0(false, errstrm.str());
460  }
461 }
462 
464 {
465  // This 2D array holds the local id's of all the vertices for every
466  // edge. For every edge, they are ordered to what we define as being
467  // Forwards.
468  const unsigned int edgeVerts[kNedges][2] = {
469  {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 4}, {3, 4}};
470 
471  int i;
472  for (i = 0; i < kNedges; i++)
473  {
474  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
475  {
477  }
478  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetGlobalID())
479  {
481  }
482  else
483  {
484  ASSERTL0(false, "Could not find matching vertex for the edge");
485  }
486  }
487 }
488 
490 {
491  int f, i;
492 
493  // These arrays represent the vector of the A and B coordinate of
494  // the local elemental coordinate system where A corresponds with
495  // the coordinate direction xi_i with the lowest index i (for that
496  // particular face) Coordinate 'B' then corresponds to the other
497  // local coordinate (i.e. with the highest index)
498  Array<OneD, NekDouble> elementAaxis(m_coordim);
499  Array<OneD, NekDouble> elementBaxis(m_coordim);
500 
501  // These arrays correspond to the local coordinate
502  // system of the face itself (i.e. the Geometry2D)
503  // faceAaxis correspond to the xi_0 axis
504  // faceBaxis correspond to the xi_1 axis
507 
508  // This is the base vertex of the face (i.e. the Geometry2D) This
509  // corresponds to thevertex with local ID 0 of the Geometry2D
510  unsigned int baseVertex;
511 
512  // The lenght of the vectors above
513  NekDouble elementAaxis_length;
514  NekDouble elementBaxis_length;
515  NekDouble faceAaxis_length;
516  NekDouble faceBaxis_length;
517 
518  // This 2D array holds the local id's of all the vertices for every
519  // face. For every face, they are ordered in such a way that the
520  // implementation below allows a unified approach for all faces.
521  const unsigned int faceVerts[kNfaces][4] = {
522  {0, 1, 2, 3},
523  {0, 1, 4, 0}, // Last four elements are triangles which only
524  {1, 2, 4, 0}, // require three vertices.
525  {3, 2, 4, 0},
526  {0, 3, 4, 0}};
527 
528  NekDouble dotproduct1 = 0.0;
529  NekDouble dotproduct2 = 0.0;
530 
531  unsigned int orientation;
532 
533  // Loop over all the faces to set up the orientation
534  for (f = 0; f < kNqfaces + kNtfaces; f++)
535  {
536  // initialisation
537  elementAaxis_length = 0.0;
538  elementBaxis_length = 0.0;
539  faceAaxis_length = 0.0;
540  faceBaxis_length = 0.0;
541 
542  dotproduct1 = 0.0;
543  dotproduct2 = 0.0;
544 
545  baseVertex = m_faces[f]->GetVid(0);
546 
547  // We are going to construct the vectors representing the A and
548  // B axis of every face. These vectors will be constructed as a
549  // vector-representation of the edges of the face. However, for
550  // both coordinate directions, we can represent the vectors by
551  // two different edges. That's why we need to make sure that we
552  // pick the edge to which the baseVertex of the
553  // Geometry2D-representation of the face belongs...
554 
555  // Compute the length of edges on a base-face
556  if (f > 0)
557  {
558  if (baseVertex == m_verts[faceVerts[f][0]]->GetVid())
559  {
560  for (i = 0; i < m_coordim; i++)
561  {
562  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
563  (*m_verts[faceVerts[f][0]])[i];
564  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
565  (*m_verts[faceVerts[f][0]])[i];
566  }
567  }
568  else if (baseVertex == m_verts[faceVerts[f][1]]->GetVid())
569  {
570  for (i = 0; i < m_coordim; i++)
571  {
572  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
573  (*m_verts[faceVerts[f][0]])[i];
574  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
575  (*m_verts[faceVerts[f][1]])[i];
576  }
577  }
578  else if (baseVertex == m_verts[faceVerts[f][2]]->GetVid())
579  {
580  for (i = 0; i < m_coordim; i++)
581  {
582  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
583  (*m_verts[faceVerts[f][2]])[i];
584  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
585  (*m_verts[faceVerts[f][0]])[i];
586  }
587  }
588  else
589  {
590  ASSERTL0(false, "Could not find matching vertex for the face");
591  }
592  }
593  else
594  {
595  if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
596  {
597  for (i = 0; i < m_coordim; i++)
598  {
599  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
600  (*m_verts[faceVerts[f][0]])[i];
601  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
602  (*m_verts[faceVerts[f][0]])[i];
603  }
604  }
605  else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
606  {
607  for (i = 0; i < m_coordim; i++)
608  {
609  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
610  (*m_verts[faceVerts[f][0]])[i];
611  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
612  (*m_verts[faceVerts[f][1]])[i];
613  }
614  }
615  else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
616  {
617  for (i = 0; i < m_coordim; i++)
618  {
619  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
620  (*m_verts[faceVerts[f][3]])[i];
621  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
622  (*m_verts[faceVerts[f][1]])[i];
623  }
624  }
625  else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
626  {
627  for (i = 0; i < m_coordim; i++)
628  {
629  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
630  (*m_verts[faceVerts[f][3]])[i];
631  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
632  (*m_verts[faceVerts[f][0]])[i];
633  }
634  }
635  else
636  {
637  ASSERTL0(false, "Could not find matching vertex for the face");
638  }
639  }
640 
641  // Now, construct the edge-vectors of the local coordinates of
642  // the Geometry2D-representation of the face
643  for (i = 0; i < m_coordim; i++)
644  {
645  int v = m_faces[f]->GetNumVerts() - 1;
646  faceAaxis[i] =
647  (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
648  faceBaxis[i] =
649  (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
650 
651  elementAaxis_length += pow(elementAaxis[i], 2);
652  elementBaxis_length += pow(elementBaxis[i], 2);
653  faceAaxis_length += pow(faceAaxis[i], 2);
654  faceBaxis_length += pow(faceBaxis[i], 2);
655  }
656 
657  elementAaxis_length = sqrt(elementAaxis_length);
658  elementBaxis_length = sqrt(elementBaxis_length);
659  faceAaxis_length = sqrt(faceAaxis_length);
660  faceBaxis_length = sqrt(faceBaxis_length);
661 
662  // Calculate the inner product of both the A-axis
663  // (i.e. Elemental A axis and face A axis)
664  for (i = 0; i < m_coordim; i++)
665  {
666  dotproduct1 += elementAaxis[i] * faceAaxis[i];
667  }
668 
669  NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
670  orientation = 0;
671 
672  // if the innerproduct is equal to the (absolute value of the ) products
673  // of the lengths of both vectors, then, the coordinate systems will NOT
674  // be transposed
675  if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
676  {
677  // if the inner product is negative, both A-axis point
678  // in reverse direction
679  if (dotproduct1 < 0.0)
680  {
681  orientation += 2;
682  }
683 
684  // calculate the inner product of both B-axis
685  for (i = 0; i < m_coordim; i++)
686  {
687  dotproduct2 += elementBaxis[i] * faceBaxis[i];
688  }
689 
690  norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
691 
692  // check that both these axis are indeed parallel
693  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
694  "These vectors should be parallel");
695 
696  // if the inner product is negative, both B-axis point
697  // in reverse direction
698  if (dotproduct2 < 0.0)
699  {
700  orientation++;
701  }
702  }
703  // The coordinate systems are transposed
704  else
705  {
706  orientation = 4;
707 
708  // Calculate the inner product between the elemental A-axis
709  // and the B-axis of the face (which are now the corresponding axis)
710  dotproduct1 = 0.0;
711  for (i = 0; i < m_coordim; i++)
712  {
713  dotproduct1 += elementAaxis[i] * faceBaxis[i];
714  }
715 
716  norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
717  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
718  "These vectors should be parallel");
719 
720  // if the result is negative, both axis point in reverse
721  // directions
722  if (dotproduct1 < 0.0)
723  {
724  orientation += 2;
725  }
726 
727  // Do the same for the other two corresponding axis
728  dotproduct2 = 0.0;
729  for (i = 0; i < m_coordim; i++)
730  {
731  dotproduct2 += elementBaxis[i] * faceAaxis[i];
732  }
733 
734  norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
735 
736  // check that both these axis are indeed parallel
737  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
738  "These vectors should be parallel");
739 
740  if (dotproduct2 < 0.0)
741  {
742  orientation++;
743  }
744  }
745 
746  orientation = orientation + 5;
747 
748  // Fill the m_forient array
749  m_forient[f] = (StdRegions::Orientation)orientation;
750  }
751 }
752 
753 void PyrGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
754 {
755  Geometry::v_Reset(curvedEdges, curvedFaces);
756 
757  for (int i = 0; i < 5; ++i)
758  {
759  m_faces[i]->Reset(curvedEdges, curvedFaces);
760  }
761 
762  SetUpXmap();
763  SetUpCoeffs(m_xmap->GetNcoeffs());
764 }
765 
767 {
768  if(!m_setupState)
769  {
770  for (int i = 0; i < 5; ++i)
771  {
772  m_faces[i]->Setup();
773  }
774  SetUpXmap();
775  SetUpCoeffs(m_xmap->GetNcoeffs());
776  m_setupState = true;
777  }
778 }
779 
780 /**
781  * @brief Set up the #m_xmap object by determining the order of each
782  * direction from derived faces.
783  */
785 {
786  vector<int> tmp;
787  int order0, order1;
788 
789  if (m_forient[0] < 9)
790  {
791  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
792  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
793  order0 = *max_element(tmp.begin(), tmp.end());
794  }
795  else
796  {
797  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
798  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
799  order0 = *max_element(tmp.begin(), tmp.end());
800  }
801 
802  if (m_forient[0] < 9)
803  {
804  tmp.clear();
805  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
806  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(3));
807  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
808  order1 = *max_element(tmp.begin(), tmp.end());
809  }
810  else
811  {
812  tmp.clear();
813  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
814  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
815  tmp.push_back(m_faces[2]->GetXmap()->GetEdgeNcoeffs(2));
816  order1 = *max_element(tmp.begin(), tmp.end());
817  }
818 
819  tmp.clear();
820  tmp.push_back(order0);
821  tmp.push_back(order1);
822  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(1));
823  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(2));
824  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(1));
825  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(2));
826  int order2 = *max_element(tmp.begin(), tmp.end());
827 
828 
829  const LibUtilities::BasisKey A1(
833  const LibUtilities::BasisKey A2(
837  const LibUtilities::BasisKey C(
841 
843  A1, A2, C);
844 }
845 
846 }
847 }
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:86
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:145
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
static const int kNfaces
Definition: PyrGeom.h:57
Principle Modified Functions .
Definition: BasisType.h:52
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:421
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:62
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within twice the min/max distance of the element.
Definition: Geometry.cpp:435
NekDouble dot(PointGeom &a)
retun the dot product between this and input a
Definition: PointGeom.cpp:188
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:850
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
Principle Modified Functions .
Definition: BasisType.h:48
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:125
STL namespace.
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:335
void ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:478
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:187
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:314
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:87
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: PyrGeom.cpp:753
void NewtonIterationForLocCoord(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, const Array< OneD, const NekDouble > &ptsz, Array< OneD, NekDouble > &Lcoords, NekDouble &resid)
Definition: Geometry3D.cpp:81
static const NekDouble kNekZeroTol
static const int kNedges
Definition: PyrGeom.h:54
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition: PyrGeom.cpp:157
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:193
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:137
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:318
static const int kNtfaces
Definition: PyrGeom.h:56
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:235
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
Definition: PointGeom.cpp:134
static const int kNqfaces
Definition: PyrGeom.h:55
3D geometry information
Definition: Geometry3D.h:67
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
Geometry is straight-sided with constant geometric factors.
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition: Geometry.h:534
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:351
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:197
NekDouble dist(PointGeom &a)
return distance between this and input a
Definition: PointGeom.cpp:181
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:643
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:298
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:784
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
Definition: PyrGeom.cpp:79
virtual void v_GenGeomFactors()
Definition: PyrGeom.cpp:110
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62
Geometry is curved or has non-constant factors.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:343
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:183
virtual int v_GetDir(const int faceidx, const int facedir) const
Definition: PyrGeom.cpp:256
Describes the specification for a Basis.
Definition: Basis.h:49
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186