Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
72  /// Determine necessary order for standard region.
73  vector<int> tmp;
74 
75  int order0, points0, order1, points1;
76 
77  if (m_forient[0] < 9)
78  {
79  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(0));
80  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(2));
81  order0 = *max_element(tmp.begin(), tmp.end());
82 
83  tmp.clear();
84  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(0));
85  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(2));
86  points0 = *max_element(tmp.begin(), tmp.end());
87  }
88  else
89  {
90  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(1));
91  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(3));
92  order0 = *max_element(tmp.begin(), tmp.end());
93 
94  tmp.clear();
95  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(1));
96  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(3));
97  points0 = *max_element(tmp.begin(), tmp.end());
98  }
99 
100  if (m_forient[0] < 9)
101  {
102  tmp.clear();
103  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(1));
104  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(3));
105  tmp.push_back(faces[2]->GetXmap()->GetEdgeNcoeffs(2));
106  order1 = *max_element(tmp.begin(), tmp.end());
107 
108  tmp.clear();
109  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(1));
110  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(3));
111  tmp.push_back(faces[2]->GetXmap()->GetEdgeNumPoints(2));
112  points1 = *max_element(tmp.begin(), tmp.end());
113  }
114  else
115  {
116  tmp.clear();
117  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(0));
118  tmp.push_back(faces[0]->GetXmap()->GetEdgeNcoeffs(2));
119  tmp.push_back(faces[2]->GetXmap()->GetEdgeNcoeffs(2));
120  order1 = *max_element(tmp.begin(), tmp.end());
121 
122  tmp.clear();
123  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(0));
124  tmp.push_back(faces[0]->GetXmap()->GetEdgeNumPoints(2));
125  tmp.push_back(faces[2]->GetXmap()->GetEdgeNumPoints(2));
126  points1 = *max_element(tmp.begin(), tmp.end());
127  }
128 
129  tmp.clear();
130  tmp.push_back(order0);
131  tmp.push_back(order1);
132  tmp.push_back(faces[1]->GetXmap()->GetEdgeNcoeffs(1));
133  tmp.push_back(faces[1]->GetXmap()->GetEdgeNcoeffs(2));
134  tmp.push_back(faces[3]->GetXmap()->GetEdgeNcoeffs(1));
135  tmp.push_back(faces[3]->GetXmap()->GetEdgeNcoeffs(2));
136  int order2 = *max_element(tmp.begin(), tmp.end());
137 
138  tmp.clear();
139  tmp.push_back(points0);
140  tmp.push_back(points1);
141  tmp.push_back(faces[1]->GetXmap()->GetEdgeNumPoints(1));
142  tmp.push_back(faces[1]->GetXmap()->GetEdgeNumPoints(2));
143  tmp.push_back(faces[3]->GetXmap()->GetEdgeNumPoints(1));
144  tmp.push_back(faces[3]->GetXmap()->GetEdgeNumPoints(2));
145  tmp.push_back(faces[1]->GetEdge(1)->GetBasis(0)->GetNumPoints());
146  tmp.push_back(faces[1]->GetEdge(2)->GetBasis(0)->GetNumPoints());
147  tmp.push_back(faces[3]->GetEdge(1)->GetBasis(0)->GetNumPoints());
148  tmp.push_back(faces[3]->GetEdge(2)->GetBasis(0)->GetNumPoints());
149  int points2 = *max_element(tmp.begin(), tmp.end());
150 
151  const LibUtilities::BasisKey A1(
155  const LibUtilities::BasisKey A2(
159  const LibUtilities::BasisKey C(
163 
165  A1, A2, C);
166  SetUpCoeffs(m_xmap->GetNcoeffs());
167  }
168 
170  {
171 
172  }
173 
175  {
177  {
178  int i;
179  GeomType Gtype = eRegular;
180 
181  v_FillGeom();
182 
183  // check to see if expansions are linear
184  for(i = 0; i < m_coordim; ++i)
185  {
186  if (m_xmap->GetBasisNumModes(0) != 2 ||
187  m_xmap->GetBasisNumModes(1) != 2 ||
188  m_xmap->GetBasisNumModes(2) != 2 )
189  {
190  Gtype = eDeformed;
191  }
192  }
193 
194  // check to see if all quadrilateral faces are parallelograms
195  if(Gtype == eRegular)
196  {
197  // Ensure each face is a parallelogram? Check this.
198  for (i = 0; i < m_coordim; i++)
199  {
200  if( fabs( (*m_verts[0])(i) -
201  (*m_verts[1])(i) +
202  (*m_verts[2])(i) -
203  (*m_verts[3])(i) )
205  {
206  Gtype = eDeformed;
207  break;
208  }
209  }
210  }
211 
213  Gtype, m_coordim, m_xmap, m_coeffs);
215  }
216  }
217 
219  const Array<OneD, const NekDouble> &coords,
220  Array<OneD, NekDouble> &Lcoords)
221  {
222  NekDouble resid = 0.0;
223 
224  v_FillGeom();
225 
226  // calculate local coordinate for coord
227  if(GetMetricInfo()->GetGtype() == eRegular)
228  { // Based on Spen's book, page 99
229 
230  // Point inside tetrahedron
231  PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
232 
233  // Edges
234  PointGeom er0, e10, e30, e40;
235  er0.Sub(r,*m_verts[0]);
236  e10.Sub(*m_verts[1],*m_verts[0]);
237  e30.Sub(*m_verts[3],*m_verts[0]);
238  e40.Sub(*m_verts[4],*m_verts[0]);
239 
240 
241  // Cross products (Normal times area)
242  PointGeom cp1030, cp3040, cp4010;
243  cp1030.Mult(e10,e30);
244  cp3040.Mult(e30,e40);
245  cp4010.Mult(e40,e10);
246 
247 
248  // Barycentric coordinates (relative volume)
249  NekDouble V = e40.dot(cp1030); // Pyramid Volume = {(e40)dot(e10)x(e30)}/4
250  NekDouble scaleFactor = 2.0/3.0;
251  NekDouble v1 = er0.dot(cp3040) / V; // volume1 = {(er0)dot(e30)x(e40)}/6
252  NekDouble v2 = er0.dot(cp4010) / V; // volume2 = {(er0)dot(e40)x(e10)}/6
253  NekDouble beta = v1 * scaleFactor;
254  NekDouble gamma = v2 * scaleFactor;
255  NekDouble delta = er0.dot(cp1030) / V; // volume3 = {(er0)dot(e10)x(e30)}/4
256 
257 
258  // Make Pyramid bigger
259  Lcoords[0] = 2.0*beta - 1.0;
260  Lcoords[1] = 2.0*gamma - 1.0;
261  Lcoords[2] = 2.0*delta - 1.0;
262  }
263  else
264  {
266  "inverse mapping must be set up to use this call");
267  }
268  return resid;
269  }
270 
272  {
273  return 5;
274  }
275 
277  {
278  return 8;
279  }
280 
281  int PyrGeom::v_GetDir(const int faceidx, const int facedir) const
282  {
283  if (faceidx == 0)
284  {
285  return facedir;
286  }
287  else if (faceidx == 1 || faceidx == 3)
288  {
289  return 2 * facedir;
290  }
291  else
292  {
293  return 1 + facedir;
294  }
295  }
296 
298  {
299  // find edge 0
300  int i,j;
301  unsigned int check;
302 
303  SegGeomSharedPtr edge;
304 
305  // First set up the 4 bottom edges
306  int f;
307  for (f = 1; f < 5; f++)
308  {
309  int nEdges = m_faces[f]->GetNumEdges();
310  check = 0;
311  for (i = 0; i < 4; i++)
312  {
313  for (j = 0; j < nEdges; j++)
314  {
315  if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
316  {
317  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[0])->GetEdge(i));
318  m_edges.push_back(edge);
319  check++;
320  }
321  }
322  }
323 
324  if (check < 1)
325  {
326  std::ostringstream errstrm;
327  errstrm << "Connected faces do not share an edge. Faces ";
328  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
329  ASSERTL0(false, errstrm.str());
330  }
331  else if (check > 1)
332  {
333  std::ostringstream errstrm;
334  errstrm << "Connected faces share more than one edge. Faces ";
335  errstrm << (m_faces[0])->GetFid() << ", " << (m_faces[f])->GetFid();
336  ASSERTL0(false, errstrm.str());
337  }
338  }
339 
340  // Then, set up the 4 vertical edges
341  check = 0;
342  for(i = 0; i < 3; i++) //Set up the vertical edge :face(1) and face(4)
343  {
344  for(j = 0; j < 3; j++)
345  {
346  if( (m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j) )
347  {
348  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[1])->GetEdge(i));
349  m_edges.push_back(edge);
350  check++;
351  }
352  }
353  }
354  if( check < 1 )
355  {
356  std::ostringstream errstrm;
357  errstrm << "Connected faces do not share an edge. Faces ";
358  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
359  ASSERTL0(false, errstrm.str());
360  }
361  else if( check > 1)
362  {
363  std::ostringstream errstrm;
364  errstrm << "Connected faces share more than one edge. Faces ";
365  errstrm << (m_faces[1])->GetFid() << ", " << (m_faces[4])->GetFid();
366  ASSERTL0(false, errstrm.str());
367  }
368 
369  // Set up vertical edges: face(1) through face(4)
370  for (f = 1; f < 4; f++)
371  {
372  check = 0;
373  for(i = 0; i < m_faces[f]->GetNumEdges(); i++)
374  {
375  for(j = 0; j < m_faces[f+1]->GetNumEdges(); j++)
376  {
377  if( (m_faces[f])->GetEid(i) == (m_faces[f+1])->GetEid(j))
378  {
379  edge = boost::dynamic_pointer_cast<SegGeom>((m_faces[f])->GetEdge(i));
380  m_edges.push_back(edge);
381  check++;
382  }
383  }
384  }
385 
386  if( check < 1 )
387  {
388  std::ostringstream errstrm;
389  errstrm << "Connected faces do not share an edge. Faces ";
390  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
391  ASSERTL0(false, errstrm.str());
392  }
393  else if( check > 1)
394  {
395  std::ostringstream errstrm;
396  errstrm << "Connected faces share more than one edge. Faces ";
397  errstrm << (m_faces[f])->GetFid() << ", " << (m_faces[f+1])->GetFid();
398  ASSERTL0(false, errstrm.str());
399  }
400  }
401  }
402 
404  {
405  // Set up the first 2 vertices (i.e. vertex 0,1)
406  if (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0) ||
407  m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1))
408  {
409  m_verts.push_back(m_edges[0]->GetVertex(1));
410  m_verts.push_back(m_edges[0]->GetVertex(0));
411  }
412  else if (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0) ||
413  m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1))
414  {
415  m_verts.push_back(m_edges[0]->GetVertex(0));
416  m_verts.push_back(m_edges[0]->GetVertex(1));
417  }
418  else
419  {
420  std::ostringstream errstrm;
421  errstrm << "Connected edges do not share a vertex. Edges ";
422  errstrm << m_edges[0]->GetEid() << ", " << m_edges[1]->GetEid();
423  ASSERTL0(false, errstrm.str());
424  }
425 
426  // set up the other bottom vertices (i.e. vertex 2,3)
427  for(int i = 1; i < 3; i++)
428  {
429  if (m_edges[i]->GetVid(0) == m_verts[i]->GetVid())
430  {
431  m_verts.push_back(m_edges[i]->GetVertex(1));
432  }
433  else if (m_edges[i]->GetVid(1) == m_verts[i]->GetVid())
434  {
435  m_verts.push_back(m_edges[i]->GetVertex(0));
436  }
437  else
438  {
439  std::ostringstream errstrm;
440  errstrm << "Connected edges do not share a vertex. Edges ";
441  errstrm << m_edges[i]->GetEid() << ", " << m_edges[i-1]->GetEid();
442  ASSERTL0(false, errstrm.str());
443  }
444  }
445 
446  // set up top vertex
447  if (m_edges[4]->GetVid(0) == m_verts[0]->GetVid())
448  {
449  m_verts.push_back(m_edges[4]->GetVertex(1));
450  }
451  else
452  {
453  m_verts.push_back(m_edges[4]->GetVertex(0));
454  }
455 
456  int check = 0;
457  for (int i = 5; i < 8; ++i)
458  {
459  if( (m_edges[i]->GetVid(0) == m_verts[i-4]->GetVid()
460  && m_edges[i]->GetVid(1) == m_verts[4]->GetVid())
461  ||(m_edges[i]->GetVid(1) == m_verts[i-4]->GetVid()
462  && m_edges[i]->GetVid(0) == m_verts[4]->GetVid()))
463  {
464  check++;
465  }
466  }
467  if (check != 3) {
468  std::ostringstream errstrm;
469  errstrm << "Connected edges do not share a vertex. Edges ";
470  errstrm << m_edges[3]->GetEid() << ", " << m_edges[2]->GetEid();
471  ASSERTL0(false, errstrm.str());
472  }
473  }
474 
476  {
477  // This 2D array holds the local id's of all the vertices for every
478  // edge. For every edge, they are ordered to what we define as being
479  // Forwards.
480  const unsigned int edgeVerts[kNedges][2] =
481  { {0,1}, {1,2}, {3,2}, {0,3}, {0,4}, {1,4}, {2,4}, {3,4} };
482 
483  int i;
484  for (i = 0; i < kNedges; i++)
485  {
486  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetVid())
487  {
489  }
490  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetVid())
491  {
493  }
494  else
495  {
496  ASSERTL0(false,"Could not find matching vertex for the edge");
497  }
498  }
499  }
500 
502  {
503  int f,i;
504 
505  // These arrays represent the vector of the A and B coordinate of
506  // the local elemental coordinate system where A corresponds with
507  // the coordinate direction xi_i with the lowest index i (for that
508  // particular face) Coordinate 'B' then corresponds to the other
509  // local coordinate (i.e. with the highest index)
510  Array<OneD,NekDouble> elementAaxis(m_coordim);
511  Array<OneD,NekDouble> elementBaxis(m_coordim);
512 
513  // These arrays correspond to the local coordinate
514  // system of the face itself (i.e. the Geometry2D)
515  // faceAaxis correspond to the xi_0 axis
516  // faceBaxis correspond to the xi_1 axis
517  Array<OneD,NekDouble> faceAaxis(m_coordim);
518  Array<OneD,NekDouble> faceBaxis(m_coordim);
519 
520  // This is the base vertex of the face (i.e. the Geometry2D) This
521  // corresponds to thevertex with local ID 0 of the Geometry2D
522  unsigned int baseVertex;
523 
524  // The lenght of the vectors above
525  NekDouble elementAaxis_length;
526  NekDouble elementBaxis_length;
527  NekDouble faceAaxis_length;
528  NekDouble faceBaxis_length;
529 
530  // This 2D array holds the local id's of all the vertices for every
531  // face. For every face, they are ordered in such a way that the
532  // implementation below allows a unified approach for all faces.
533  const unsigned int faceVerts[kNfaces][4] = {
534  {0,1,2,3},
535  {0,1,4,0}, // Last four elements are triangles which only
536  {1,2,4,0}, // require three vertices.
537  {3,2,4,0},
538  {0,3,4,0}
539  };
540 
541  NekDouble dotproduct1 = 0.0;
542  NekDouble dotproduct2 = 0.0;
543 
544  unsigned int orientation;
545 
546  // Loop over all the faces to set up the orientation
547  for(f = 0; f < kNqfaces + kNtfaces; f++)
548  {
549  // initialisation
550  elementAaxis_length = 0.0;
551  elementBaxis_length = 0.0;
552  faceAaxis_length = 0.0;
553  faceBaxis_length = 0.0;
554 
555  dotproduct1 = 0.0;
556  dotproduct2 = 0.0;
557 
558  baseVertex = m_faces[f]->GetVid(0);
559 
560  // We are going to construct the vectors representing the A and
561  // B axis of every face. These vectors will be constructed as a
562  // vector-representation of the edges of the face. However, for
563  // both coordinate directions, we can represent the vectors by
564  // two different edges. That's why we need to make sure that we
565  // pick the edge to which the baseVertex of the
566  // Geometry2D-representation of the face belongs...
567 
568  // Compute the length of edges on a base-face
569  if (f > 0)
570  {
571  for(i = 0; i < m_coordim; i++)
572  {
573  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
574  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
575  }
576  }
577  else
578  {
579  if( baseVertex == m_verts[ faceVerts[f][0] ]->GetVid() )
580  {
581  for(i = 0; i < m_coordim; i++)
582  {
583  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
584  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
585  }
586  }
587  else if( baseVertex == m_verts[ faceVerts[f][1] ]->GetVid() )
588  {
589  for(i = 0; i < m_coordim; i++)
590  {
591  elementAaxis[i] = (*m_verts[ faceVerts[f][1] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
592  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
593  }
594  }
595  else if( baseVertex == m_verts[ faceVerts[f][2] ]->GetVid() )
596  {
597  for(i = 0; i < m_coordim; i++)
598  {
599  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
600  elementBaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][1] ])[i];
601  }
602  }
603  else if( baseVertex == m_verts[ faceVerts[f][3] ]->GetVid() )
604  {
605  for(i = 0; i < m_coordim; i++)
606  {
607  elementAaxis[i] = (*m_verts[ faceVerts[f][2] ])[i] - (*m_verts[ faceVerts[f][3] ])[i];
608  elementBaxis[i] = (*m_verts[ faceVerts[f][3] ])[i] - (*m_verts[ faceVerts[f][0] ])[i];
609  }
610  }
611  else
612  {
613  ASSERTL0(false, "Could not find matching vertex for the face");
614  }
615  }
616 
617  // Now, construct the edge-vectors of the local coordinates of
618  // the Geometry2D-representation of the face
619  for(i = 0; i < m_coordim; i++)
620  {
621  int v = m_faces[f]->GetNumVerts()-1;
622  faceAaxis[i] = (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
623  faceBaxis[i] = (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
624 
625  elementAaxis_length += pow(elementAaxis[i],2);
626  elementBaxis_length += pow(elementBaxis[i],2);
627  faceAaxis_length += pow(faceAaxis[i],2);
628  faceBaxis_length += pow(faceBaxis[i],2);
629  }
630 
631  elementAaxis_length = sqrt(elementAaxis_length);
632  elementBaxis_length = sqrt(elementBaxis_length);
633  faceAaxis_length = sqrt(faceAaxis_length);
634  faceBaxis_length = sqrt(faceBaxis_length);
635 
636  // Calculate the inner product of both the A-axis
637  // (i.e. Elemental A axis and face A axis)
638  for(i = 0 ; i < m_coordim; i++)
639  {
640  dotproduct1 += elementAaxis[i]*faceAaxis[i];
641  }
642 
643  orientation = 0;
644  // if the innerproduct is equal to the (absolute value of the ) products of the lengths
645  // of both vectors, then, the coordinate systems will NOT be transposed
646  if( fabs(elementAaxis_length*faceAaxis_length - fabs(dotproduct1)) < NekConstants::kNekZeroTol )
647  {
648  // if the inner product is negative, both A-axis point
649  // in reverse direction
650  if(dotproduct1 < 0.0)
651  {
652  orientation += 2;
653  }
654 
655  // calculate the inner product of both B-axis
656  for(i = 0 ; i < m_coordim; i++)
657  {
658  dotproduct2 += elementBaxis[i]*faceBaxis[i];
659  }
660 
661  // if the inner product is negative, both B-axis point
662  // in reverse direction
663  if( dotproduct2 < 0.0 )
664  {
665  orientation++;
666  }
667  }
668  // The coordinate systems are transposed
669  else
670  {
671  orientation = 4;
672 
673  // Calculate the inner product between the elemental A-axis
674  // and the B-axis of the face (which are now the corresponding axis)
675  dotproduct1 = 0.0;
676  for(i = 0 ; i < m_coordim; i++)
677  {
678  dotproduct1 += elementAaxis[i]*faceBaxis[i];
679  }
680 
681  // check that both these axis are indeed parallel
682  if (fabs(elementAaxis_length*faceBaxis_length
683  - fabs(dotproduct1)) > NekConstants::kNekZeroTol)
684  {
685  cout << "Warning: Prism axes not parallel" << endl;
686  }
687 
688  // if the result is negative, both axis point in reverse
689  // directions
690  if(dotproduct1 < 0.0)
691  {
692  orientation += 2;
693  }
694 
695  // Do the same for the other two corresponding axis
696  dotproduct2 = 0.0;
697  for(i = 0 ; i < m_coordim; i++)
698  {
699  dotproduct2 += elementBaxis[i]*faceAaxis[i];
700  }
701 
702  // check that both these axis are indeed parallel
703  if (fabs(elementBaxis_length*faceAaxis_length
704  - fabs(dotproduct2)) > NekConstants::kNekZeroTol)
705  {
706  cout << "Warning: Prism axes not parallel" << endl;
707  }
708 
709  if( dotproduct2 < 0.0 )
710  {
711  orientation++;
712  }
713  }
714 
715  orientation = orientation + 5;
716 
717  // Fill the m_forient array
718  m_forient[f] = (StdRegions::Orientation) orientation;
719  }
720  }
721  }; //end of namespace
722 }; //end of namespace