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