Nektar++
TetGeom.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: TetGeom.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: Tetrahedral geometry information.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <SpatialDomains/TetGeom.h>
36 
38 #include <StdRegions/StdTetExp.h>
39 #include <SpatialDomains/SegGeom.h>
40 
41 using namespace std;
42 
43 namespace Nektar
44 {
45 namespace SpatialDomains
46 {
47 const unsigned int TetGeom::VertexEdgeConnectivity[4][3] = {
48  {0, 2, 3}, {0, 1, 4}, {1, 2, 5}, {3, 4, 5}};
49 const unsigned int TetGeom::VertexFaceConnectivity[4][3] = {
50  {0, 1, 3}, {0, 1, 2}, {0, 2, 3}, {1, 2, 3}};
51 const unsigned int TetGeom::EdgeFaceConnectivity[6][2] = {
52  {0, 1}, {0, 2}, {0, 3}, {1, 3}, {1, 2}, {2, 3}};
53 
54 TetGeom::TetGeom()
55 {
56  m_shapeType = LibUtilities::eTetrahedron;
57 }
58 
59 TetGeom::TetGeom(int id, const TriGeomSharedPtr faces[])
60  : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
61 {
63  m_globalID = id;
64 
65  /// Copy the face shared pointers
66  m_faces.insert(m_faces.begin(), faces, faces + TetGeom::kNfaces);
67 
68  /// Set up orientation vectors with correct amount of elements.
69  m_eorient.resize(kNedges);
70  m_forient.resize(kNfaces);
71 
76 }
77 
79 {
80 }
81 
82 /**
83  * @brief Determines if a point specified in global coordinates is
84  * located within this tetrahedral geometry and return local caretsian
85  * coordinates
86  */
88  Array<OneD, NekDouble> &locCoord,
89  NekDouble tol,
90  NekDouble &resid)
91 {
92  //Rough check if within twice min/max point
93  if (GetMetricInfo()->GetGtype() != eRegular)
94  {
95  if (!MinMaxCheck(gloCoord))
96  {
97  return false;
98  }
99  }
100 
101  // Convert to the local (eta) coordinates.
102  resid = GetLocCoords(gloCoord, locCoord);
103 
104  // Check local coordinate is within cartesian bounds.
105  if (locCoord[0] >= -(1 + tol) && locCoord[1] >= -(1 + tol) &&
106  locCoord[2] >= -(1 + tol) &&
107  locCoord[0] + locCoord[1] + locCoord[2] <= -1 + tol)
108  {
109  return true;
110  }
111 
112  //Clamp local coords
113  ClampLocCoords(locCoord, tol);
114 
115  return false;
116 }
117 
119  Array<OneD, NekDouble> &Lcoords)
120 {
121  NekDouble ptdist = 1e6;
122 
123  // calculate local coordinates (eta) for coord
124  if (GetMetricInfo()->GetGtype() == eRegular)
125  {
126  // Point inside tetrahedron
127  PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
128 
129  // Edges
130  PointGeom er0, e10, e20, e30;
131  er0.Sub(r, *m_verts[0]);
132  e10.Sub(*m_verts[1], *m_verts[0]);
133  e20.Sub(*m_verts[2], *m_verts[0]);
134  e30.Sub(*m_verts[3], *m_verts[0]);
135 
136  // Cross products (Normal times area)
137  PointGeom cp1020, cp2030, cp3010;
138  cp1020.Mult(e10, e20);
139  cp2030.Mult(e20, e30);
140  cp3010.Mult(e30, e10);
141 
142  // Barycentric coordinates (relative volume)
143  NekDouble V = e30.dot(cp1020); // Tet Volume={(e30)dot(e10)x(e20)}/6
144  NekDouble beta = er0.dot(cp2030) / V; // volume1={(er0)dot(e20)x(e30)}/6
145  NekDouble gamma = er0.dot(cp3010) / V; // volume1={(er0)dot(e30)x(e10)}/6
146  NekDouble delta = er0.dot(cp1020) / V; // volume1={(er0)dot(e10)x(e20)}/6
147 
148  // Make tet bigger
149  Lcoords[0] = 2.0 * beta - 1.0;
150  Lcoords[1] = 2.0 * gamma - 1.0;
151  Lcoords[2] = 2.0 * delta - 1.0;
152 
153  // Set ptdist to distance to nearest vertex
154  for (int i = 0; i < 4; ++i)
155  {
156  ptdist = min(ptdist, r.dist(*m_verts[i]));
157  }
158  }
159  else
160  {
161  v_FillGeom();
162 
163  // Determine nearest point of coords to values in m_xmap
164  int npts = m_xmap->GetTotPoints();
165  Array<OneD, NekDouble> ptsx(npts), ptsy(npts), ptsz(npts);
166  Array<OneD, NekDouble> tmp1(npts), tmp2(npts);
167 
168  m_xmap->BwdTrans(m_coeffs[0], ptsx);
169  m_xmap->BwdTrans(m_coeffs[1], ptsy);
170  m_xmap->BwdTrans(m_coeffs[2], ptsz);
171 
172  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
173  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
174  const Array<OneD, const NekDouble> zc = m_xmap->GetPoints(2);
175 
176  // guess the first local coords based on nearest point
177  Vmath::Sadd(npts, -coords[0], ptsx, 1, tmp1, 1);
178  Vmath::Vmul(npts, tmp1, 1, tmp1, 1, tmp1, 1);
179  Vmath::Sadd(npts, -coords[1], ptsy, 1, tmp2, 1);
180  Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
181  Vmath::Sadd(npts, -coords[2], ptsz, 1, tmp2, 1);
182  Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
183 
184  int min_i = Vmath::Imin(npts, tmp1, 1);
185 
186  // distance from coordinate to nearest point for return value.
187  ptdist = sqrt(tmp1[min_i]);
188 
189  // Get collapsed coordinate
190  int qa = za.num_elements(), qb = zb.num_elements();
191  Lcoords[2] = zc[min_i / (qa * qb)];
192  min_i = min_i % (qa * qb);
193  Lcoords[1] = zb[min_i / qa];
194  Lcoords[0] = za[min_i % qa];
195 
196  // recover cartesian coordinate from collapsed coordinate.
197  Lcoords[1] = (1.0 + Lcoords[0]) * (1.0 - Lcoords[2]) / 2 - 1.0;
198  Lcoords[0] = (1.0 + Lcoords[0]) * (-Lcoords[1] - Lcoords[2]) / 2 - 1.0;
199 
200  // Perform newton iteration to find local coordinates
201  NekDouble resid = 0.0;
202  NewtonIterationForLocCoord(coords, ptsx, ptsy, ptsz, Lcoords, resid);
203  }
204  return ptdist;
205 }
206 
207 int TetGeom::v_GetDir(const int faceidx, const int facedir) const
208 {
209  if (faceidx == 0)
210  {
211  return facedir;
212  }
213  else if (faceidx == 1)
214  {
215  return 2 * facedir;
216  }
217  else
218  {
219  return 1 + facedir;
220  }
221 }
222 
223 int TetGeom::v_GetVertexEdgeMap(const int i, const int j) const
224 {
225  return VertexEdgeConnectivity[i][j];
226 }
227 
228 int TetGeom::v_GetVertexFaceMap(const int i, const int j) const
229 {
230  return VertexFaceConnectivity[i][j];
231 }
232 
233 int TetGeom::v_GetEdgeFaceMap(const int i, const int j) const
234 {
235  return EdgeFaceConnectivity[i][j];
236 }
237 
239 {
240 
241  // find edge 0
242  int i, j;
243  unsigned int check;
244 
245  SegGeomSharedPtr edge;
246 
247  // First set up the 3 bottom edges
248 
249  if (m_faces[0]->GetEid(0) != m_faces[1]->GetEid(0))
250  {
251  std::ostringstream errstrm;
252  errstrm << "Local edge 0 (eid=" << m_faces[0]->GetEid(0);
253  errstrm << ") on face " << m_faces[0]->GetGlobalID();
254  errstrm << " must be the same as local edge 0 (eid="
255  << m_faces[1]->GetEid(0);
256  errstrm << ") on face " << m_faces[1]->GetGlobalID();
257  ASSERTL0(false, errstrm.str());
258  }
259 
260  int faceConnected;
261  for (faceConnected = 1; faceConnected < 4; faceConnected++)
262  {
263  check = 0;
264  for (i = 0; i < 3; i++)
265  {
266  if ((m_faces[0])->GetEid(i) == (m_faces[faceConnected])->GetEid(0))
267  {
268  edge = dynamic_pointer_cast<SegGeom>(
269  (m_faces[0])->GetEdge(i));
270  m_edges.push_back(edge);
271  check++;
272  }
273  }
274 
275  if (check < 1)
276  {
277  std::ostringstream errstrm;
278  errstrm << "Face 0 does not share an edge with first edge of "
279  "adjacent face. Faces ";
280  errstrm << (m_faces[0])->GetGlobalID() << ", "
281  << (m_faces[faceConnected])->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[0])->GetGlobalID() << ", "
289  << (m_faces[faceConnected])->GetGlobalID();
290  ASSERTL0(false, errstrm.str());
291  }
292  }
293 
294  // Then, set up the 3 vertical edges
295  check = 0;
296  for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(3)
297  {
298  for (j = 0; j < 3; j++)
299  {
300  if ((m_faces[1])->GetEid(i) == (m_faces[3])->GetEid(j))
301  {
302  edge = dynamic_pointer_cast<SegGeom>(
303  (m_faces[1])->GetEdge(i));
304  m_edges.push_back(edge);
305  check++;
306  }
307  }
308  }
309  if (check < 1)
310  {
311  std::ostringstream errstrm;
312  errstrm << "Connected faces do not share an edge. Faces ";
313  errstrm << (m_faces[1])->GetGlobalID() << ", "
314  << (m_faces[3])->GetGlobalID();
315  ASSERTL0(false, errstrm.str());
316  }
317  else if (check > 1)
318  {
319  std::ostringstream errstrm;
320  errstrm << "Connected faces share more than one edge. Faces ";
321  errstrm << (m_faces[1])->GetGlobalID() << ", "
322  << (m_faces[3])->GetGlobalID();
323  ASSERTL0(false, errstrm.str());
324  }
325  // Set up vertical edges: face(1) through face(3)
326  for (faceConnected = 1; faceConnected < 3; faceConnected++)
327  {
328  check = 0;
329  for (i = 0; i < 3; i++)
330  {
331  for (j = 0; j < 3; j++)
332  {
333  if ((m_faces[faceConnected])->GetEid(i) ==
334  (m_faces[faceConnected + 1])->GetEid(j))
335  {
336  edge = dynamic_pointer_cast<SegGeom>(
337  (m_faces[faceConnected])->GetEdge(i));
338  m_edges.push_back(edge);
339  check++;
340  }
341  }
342  }
343 
344  if (check < 1)
345  {
346  std::ostringstream errstrm;
347  errstrm << "Connected faces do not share an edge. Faces ";
348  errstrm << (m_faces[faceConnected])->GetGlobalID() << ", "
349  << (m_faces[faceConnected + 1])->GetGlobalID();
350  ASSERTL0(false, errstrm.str());
351  }
352  else if (check > 1)
353  {
354  std::ostringstream errstrm;
355  errstrm << "Connected faces share more than one edge. Faces ";
356  errstrm << (m_faces[faceConnected])->GetGlobalID() << ", "
357  << (m_faces[faceConnected + 1])->GetGlobalID();
358  ASSERTL0(false, errstrm.str());
359  }
360  }
361 }
362 
364 {
365 
366  // Set up the first 2 vertices (i.e. vertex 0,1)
367  if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
368  (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
369  {
370  m_verts.push_back(m_edges[0]->GetVertex(1));
371  m_verts.push_back(m_edges[0]->GetVertex(0));
372  }
373  else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
374  (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
375  {
376  m_verts.push_back(m_edges[0]->GetVertex(0));
377  m_verts.push_back(m_edges[0]->GetVertex(1));
378  }
379  else
380  {
381  std::ostringstream errstrm;
382  errstrm << "Connected edges do not share a vertex. Edges ";
383  errstrm << m_edges[0]->GetGlobalID() << ", "
384  << m_edges[1]->GetGlobalID();
385  ASSERTL0(false, errstrm.str());
386  }
387 
388  // set up the other bottom vertices (i.e. vertex 2)
389  for (int i = 1; i < 2; i++)
390  {
391  if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
392  {
393  m_verts.push_back(m_edges[i]->GetVertex(1));
394  }
395  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
396  {
397  m_verts.push_back(m_edges[i]->GetVertex(0));
398  }
399  else
400  {
401  std::ostringstream errstrm;
402  errstrm << "Connected edges do not share a vertex. Edges ";
403  errstrm << m_edges[i]->GetGlobalID() << ", "
404  << m_edges[i - 1]->GetGlobalID();
405  ASSERTL0(false, errstrm.str());
406  }
407  }
408 
409  // set up top vertex
410  if (m_edges[3]->GetVid(0) == m_verts[0]->GetGlobalID())
411  {
412  m_verts.push_back(m_edges[3]->GetVertex(1));
413  }
414  else
415  {
416  m_verts.push_back(m_edges[3]->GetVertex(0));
417  }
418 
419  // Check the other edges match up.
420  int check = 0;
421  for (int i = 4; i < 6; ++i)
422  {
423  if ((m_edges[i]->GetVid(0) == m_verts[i - 3]->GetGlobalID() &&
424  m_edges[i]->GetVid(1) == m_verts[3]->GetGlobalID()) ||
425  (m_edges[i]->GetVid(1) == m_verts[i - 3]->GetGlobalID() &&
426  m_edges[i]->GetVid(0) == m_verts[3]->GetGlobalID()))
427  {
428  check++;
429  }
430  }
431  if (check != 2)
432  {
433  std::ostringstream errstrm;
434  errstrm << "Connected edges do not share a vertex. Edges ";
435  errstrm << m_edges[3]->GetGlobalID() << ", "
436  << m_edges[2]->GetGlobalID();
437  ASSERTL0(false, errstrm.str());
438  }
439 }
440 
442 {
443 
444  // This 2D array holds the local id's of all the vertices
445  // for every edge. For every edge, they are ordered to what we
446  // define as being Forwards
447  const unsigned int edgeVerts[kNedges][2] = {
448  {0, 1}, {1, 2}, {0, 2}, {0, 3}, {1, 3}, {2, 3}};
449 
450  int i;
451  for (i = 0; i < kNedges; i++)
452  {
453  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
454  {
456  }
457  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetGlobalID())
458  {
460  }
461  else
462  {
463  ASSERTL0(false, "Could not find matching vertex for the edge");
464  }
465  }
466 }
467 
469 {
470 
471  int f, i;
472 
473  // These arrays represent the vector of the A and B
474  // coordinate of the local elemental coordinate system
475  // where A corresponds with the coordinate direction xi_i
476  // with the lowest index i (for that particular face)
477  // Coordinate 'B' then corresponds to the other local
478  // coordinate (i.e. with the highest index)
479  Array<OneD, NekDouble> elementAaxis(m_coordim);
480  Array<OneD, NekDouble> elementBaxis(m_coordim);
481 
482  // These arrays correspond to the local coordinate
483  // system of the face itself (i.e. the Geometry2D)
484  // faceAaxis correspond to the xi_0 axis
485  // faceBaxis correspond to the xi_1 axis
488 
489  // This is the base vertex of the face (i.e. the Geometry2D)
490  // This corresponds to thevertex with local ID 0 of the
491  // Geometry2D
492  unsigned int baseVertex;
493 
494  // The lenght of the vectors above
495  NekDouble elementAaxis_length;
496  NekDouble elementBaxis_length;
497  NekDouble faceAaxis_length;
498  NekDouble faceBaxis_length;
499 
500  // This 2D array holds the local id's of all the vertices
501  // for every face. For every face, they are ordered in such
502  // a way that the implementation below allows a unified approach
503  // for all faces.
504  const unsigned int faceVerts[kNfaces][TriGeom::kNverts] = {
505  {0, 1, 2}, {0, 1, 3}, {1, 2, 3}, {0, 2, 3}};
506 
507  NekDouble dotproduct1 = 0.0;
508  NekDouble dotproduct2 = 0.0;
509 
510  unsigned int orientation;
511 
512  // Loop over all the faces to set up the orientation
513  for (f = 0; f < kNqfaces + kNtfaces; f++)
514  {
515  // initialisation
516  elementAaxis_length = 0.0;
517  elementBaxis_length = 0.0;
518  faceAaxis_length = 0.0;
519  faceBaxis_length = 0.0;
520 
521  dotproduct1 = 0.0;
522  dotproduct2 = 0.0;
523 
524  baseVertex = m_faces[f]->GetVid(0);
525 
526  // We are going to construct the vectors representing the A and B axis
527  // of every face. These vectors will be constructed as a
528  // vector-representation
529  // of the edges of the face. However, for both coordinate directions, we
530  // can
531  // represent the vectors by two different edges. That's why we need to
532  // make sure that
533  // we pick the edge to which the baseVertex of the
534  // Geometry2D-representation of the face
535  // belongs...
536 
537  // Compute the length of edges on a base-face
538  if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
539  {
540  for (i = 0; i < m_coordim; i++)
541  {
542  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
543  (*m_verts[faceVerts[f][0]])[i];
544  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
545  (*m_verts[faceVerts[f][0]])[i];
546  }
547  }
548  else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
549  {
550  for (i = 0; i < m_coordim; i++)
551  {
552  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
553  (*m_verts[faceVerts[f][0]])[i];
554  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
555  (*m_verts[faceVerts[f][1]])[i];
556  }
557  }
558  else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
559  {
560  for (i = 0; i < m_coordim; i++)
561  {
562  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
563  (*m_verts[faceVerts[f][2]])[i];
564  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
565  (*m_verts[faceVerts[f][0]])[i];
566  }
567  }
568  else
569  {
570  ASSERTL0(false, "Could not find matching vertex for the face");
571  }
572 
573  // Now, construct the edge-vectors of the local coordinates of
574  // the Geometry2D-representation of the face
575  for (i = 0; i < m_coordim; i++)
576  {
577  faceAaxis[i] =
578  (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
579  faceBaxis[i] =
580  (*m_faces[f]->GetVertex(2))[i] - (*m_faces[f]->GetVertex(0))[i];
581 
582  elementAaxis_length += pow(elementAaxis[i], 2);
583  elementBaxis_length += pow(elementBaxis[i], 2);
584  faceAaxis_length += pow(faceAaxis[i], 2);
585  faceBaxis_length += pow(faceBaxis[i], 2);
586  }
587 
588  elementAaxis_length = sqrt(elementAaxis_length);
589  elementBaxis_length = sqrt(elementBaxis_length);
590  faceAaxis_length = sqrt(faceAaxis_length);
591  faceBaxis_length = sqrt(faceBaxis_length);
592 
593  // Calculate the inner product of both the A-axis
594  // (i.e. Elemental A axis and face A axis)
595  for (i = 0; i < m_coordim; i++)
596  {
597  dotproduct1 += elementAaxis[i] * faceAaxis[i];
598  }
599 
600  NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
601  orientation = 0;
602 
603  // if the innerproduct is equal to the (absolute value of the ) products
604  // of the lengths of both vectors, then, the coordinate systems will NOT
605  // be transposed
606  if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
607  {
608  // if the inner product is negative, both A-axis point
609  // in reverse direction
610  if (dotproduct1 < 0.0)
611  {
612  orientation += 2;
613  }
614 
615  // calculate the inner product of both B-axis
616  for (i = 0; i < m_coordim; i++)
617  {
618  dotproduct2 += elementBaxis[i] * faceBaxis[i];
619  }
620 
621  norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
622 
623  // check that both these axis are indeed parallel
624  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
625  "These vectors should be parallel");
626 
627  // if the inner product is negative, both B-axis point
628  // in reverse direction
629  if (dotproduct2 < 0.0)
630  {
631  orientation++;
632  }
633  }
634  // The coordinate systems are transposed
635  else
636  {
637  orientation = 4;
638 
639  // Calculate the inner product between the elemental A-axis
640  // and the B-axis of the face (which are now the corresponding axis)
641  dotproduct1 = 0.0;
642  for (i = 0; i < m_coordim; i++)
643  {
644  dotproduct1 += elementAaxis[i] * faceBaxis[i];
645  }
646 
647  norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
648 
649  // check that both these axis are indeed parallel
650  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
651  "These vectors should be parallel");
652 
653  // if the result is negative, both axis point in reverse
654  // directions
655  if (dotproduct1 < 0.0)
656  {
657  orientation += 2;
658  }
659 
660  // Do the same for the other two corresponding axis
661  dotproduct2 = 0.0;
662  for (i = 0; i < m_coordim; i++)
663  {
664  dotproduct2 += elementBaxis[i] * faceAaxis[i];
665  }
666 
667  norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
668 
669  // check that both these axis are indeed parallel
670  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
671  "These vectors should be parallel");
672 
673  if (dotproduct2 < 0.0)
674  {
675  orientation++;
676  }
677  }
678 
679  orientation = orientation + 5;
680 
681  // Fill the m_forient array
682  m_forient[f] = (StdRegions::Orientation)orientation;
683  }
684 }
685 
686 void TetGeom::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
687 {
688  Geometry::v_Reset(curvedEdges, curvedFaces);
689 
690  for (int i = 0; i < 4; ++i)
691  {
692  m_faces[i]->Reset(curvedEdges, curvedFaces);
693  }
694 }
695 
697 {
698  if(!m_setupState)
699  {
700  for (int i = 0; i < 4; ++i)
701  {
702  m_faces[i]->Setup();
703  }
704  SetUpXmap();
705  SetUpCoeffs(m_xmap->GetNcoeffs());
706  m_setupState = true;
707  }
708 }
709 
710 /**
711  * Generate the geometry factors for this element.
712  */
714 {
715  if(!m_setupState)
716  {
718  }
719 
721  {
722  GeomType Gtype = eRegular;
723 
724  v_FillGeom();
725 
726  // check to see if expansions are linear
727  for (int i = 0; i < m_coordim; ++i)
728  {
729  if (m_xmap->GetBasisNumModes(0) != 2 ||
730  m_xmap->GetBasisNumModes(1) != 2 ||
731  m_xmap->GetBasisNumModes(2) != 2)
732  {
733  Gtype = eDeformed;
734  }
735  }
736 
738  Gtype, m_coordim, m_xmap, m_coeffs);
740  }
741 }
742 
743 /**
744  * @brief Set up the #m_xmap object by determining the order of each
745  * direction from derived faces.
746  */
748 {
749  vector<int> tmp;
750  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(0));
751  int order0 = *max_element(tmp.begin(), tmp.end());
752 
753  tmp.clear();
754  tmp.push_back(order0);
755  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(1));
756  tmp.push_back(m_faces[0]->GetXmap()->GetEdgeNcoeffs(2));
757  int order1 = *max_element(tmp.begin(), tmp.end());
758 
759  tmp.clear();
760  tmp.push_back(order0);
761  tmp.push_back(order1);
762  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(1));
763  tmp.push_back(m_faces[1]->GetXmap()->GetEdgeNcoeffs(2));
764  tmp.push_back(m_faces[3]->GetXmap()->GetEdgeNcoeffs(1));
765  int order2 = *max_element(tmp.begin(), tmp.end());
766 
769  order0,
771  const LibUtilities::BasisKey B(
773  order1,
776  const LibUtilities::BasisKey C(
778  order2,
781 
783 }
784 
785 }
786 }
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
Principle Modified Functions .
Definition: BasisType.h:50
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
static const int kNverts
Definition: TriGeom.h:74
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_GenGeomFactors()
Definition: TetGeom.cpp:713
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
static const int kNtfaces
Definition: TetGeom.h:57
static const int kNfaces
Definition: TetGeom.h:58
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
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
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: TetGeom.cpp:228
static const unsigned int EdgeFaceConnectivity[6][2]
Definition: TetGeom.h:85
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
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
static const NekDouble kNekZeroTol
static const int kNedges
Definition: TetGeom.h:55
Principle Modified Functions .
Definition: BasisType.h:49
static const unsigned int VertexFaceConnectivity[4][3]
Definition: TetGeom.h:84
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: TetGeom.cpp:747
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
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 kNqfaces
Definition: TetGeom.h:56
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:235
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: TetGeom.cpp:118
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
Definition: PointGeom.cpp:134
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: TetGeom.cpp:686
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
static const unsigned int VertexEdgeConnectivity[4][3]
Definition: TetGeom.h:83
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:298
GeomType
Indicates the type of element geometry.
#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
Describes the specification for a Basis.
Definition: Basis.h:49
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: TetGeom.cpp:223
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: TetGeom.cpp:233
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
virtual int v_GetDir(const int faceidx, const int facedir) const
Definition: TetGeom.cpp:207
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
Determines if a point specified in global coordinates is located within this tetrahedral geometry and...
Definition: TetGeom.cpp:87