47 namespace SpatialDomains
56 :
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
95 if (
m_xmap->GetBasisNumModes(0) != 2 ||
96 m_xmap->GetBasisNumModes(1) != 2 ||
97 m_xmap->GetBasisNumModes(2) != 2)
131 else if (faceidx == 1 || faceidx == 3)
151 for (f = 1; f < 5; f++)
153 int nEdges =
m_faces[f]->GetNumEdges();
155 for (i = 0; i < 4; i++)
157 for (j = 0; j < nEdges; j++)
161 edge = dynamic_pointer_cast<SegGeom>(
171 std::ostringstream errstrm;
172 errstrm <<
"Connected faces do not share an edge. Faces ";
179 std::ostringstream errstrm;
180 errstrm <<
"Connected faces share more than one edge. Faces ";
189 for (i = 0; i < 3; i++)
191 for (j = 0; j < 3; j++)
195 edge = dynamic_pointer_cast<SegGeom>(
204 std::ostringstream errstrm;
205 errstrm <<
"Connected faces do not share an edge. Faces ";
212 std::ostringstream errstrm;
213 errstrm <<
"Connected faces share more than one edge. Faces ";
220 for (f = 1; f < 4; f++)
223 for (i = 0; i <
m_faces[f]->GetNumEdges(); i++)
225 for (j = 0; j <
m_faces[f + 1]->GetNumEdges(); j++)
229 edge = dynamic_pointer_cast<SegGeom>(
239 std::ostringstream errstrm;
240 errstrm <<
"Connected faces do not share an edge. Faces ";
247 std::ostringstream errstrm;
248 errstrm <<
"Connected faces share more than one edge. Faces ";
273 std::ostringstream errstrm;
274 errstrm <<
"Connected edges do not share a vertex. Edges ";
275 errstrm <<
m_edges[0]->GetGlobalID() <<
", "
281 for (
int i = 1; i < 3; i++)
293 std::ostringstream errstrm;
294 errstrm <<
"Connected edges do not share a vertex. Edges ";
295 errstrm <<
m_edges[i]->GetGlobalID() <<
", "
296 <<
m_edges[i - 1]->GetGlobalID();
312 for (
int i = 5; i < 8; ++i)
324 std::ostringstream errstrm;
325 errstrm <<
"Connected edges do not share a vertex. Edges ";
326 errstrm <<
m_edges[3]->GetGlobalID() <<
", "
337 const unsigned int edgeVerts[
kNedges][2] = {
338 {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 4}, {3, 4}};
353 ASSERTL0(
false,
"Could not find matching vertex for the edge");
379 unsigned int baseVertex;
390 const unsigned int faceVerts[
kNfaces][4] = {
400 unsigned int orientation;
406 elementAaxis_length = 0.0;
407 elementBaxis_length = 0.0;
408 faceAaxis_length = 0.0;
409 faceBaxis_length = 0.0;
414 baseVertex =
m_faces[f]->GetVid(0);
431 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
432 (*
m_verts[faceVerts[f][0]])[i];
433 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
434 (*
m_verts[faceVerts[f][0]])[i];
441 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
442 (*
m_verts[faceVerts[f][0]])[i];
443 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
444 (*
m_verts[faceVerts[f][1]])[i];
451 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
452 (*
m_verts[faceVerts[f][2]])[i];
453 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
454 (*
m_verts[faceVerts[f][0]])[i];
459 ASSERTL0(
false,
"Could not find matching vertex for the face");
468 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
469 (*
m_verts[faceVerts[f][0]])[i];
470 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
471 (*
m_verts[faceVerts[f][0]])[i];
478 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
479 (*
m_verts[faceVerts[f][0]])[i];
480 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
481 (*
m_verts[faceVerts[f][1]])[i];
488 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
489 (*
m_verts[faceVerts[f][3]])[i];
490 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
491 (*
m_verts[faceVerts[f][1]])[i];
498 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
499 (*
m_verts[faceVerts[f][3]])[i];
500 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
501 (*
m_verts[faceVerts[f][0]])[i];
506 ASSERTL0(
false,
"Could not find matching vertex for the face");
514 int v =
m_faces[f]->GetNumVerts() - 1;
520 elementAaxis_length += pow(elementAaxis[i], 2);
521 elementBaxis_length += pow(elementBaxis[i], 2);
522 faceAaxis_length += pow(faceAaxis[i], 2);
523 faceBaxis_length += pow(faceBaxis[i], 2);
526 elementAaxis_length =
sqrt(elementAaxis_length);
527 elementBaxis_length =
sqrt(elementBaxis_length);
528 faceAaxis_length =
sqrt(faceAaxis_length);
529 faceBaxis_length =
sqrt(faceBaxis_length);
535 dotproduct1 += elementAaxis[i] * faceAaxis[i];
538 NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
548 if (dotproduct1 < 0.0)
556 dotproduct2 += elementBaxis[i] * faceBaxis[i];
559 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
563 "These vectors should be parallel");
567 if (dotproduct2 < 0.0)
582 dotproduct1 += elementAaxis[i] * faceBaxis[i];
585 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
587 "These vectors should be parallel");
591 if (dotproduct1 < 0.0)
600 dotproduct2 += elementBaxis[i] * faceAaxis[i];
603 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
607 "These vectors should be parallel");
609 if (dotproduct2 < 0.0)
615 orientation = orientation + 5;
620 "Orientation of triangular face (id = " +
622 ") is inconsistent with face "+
623 boost::lexical_cast<string>(f) +
624 " of pyramid element (id = "+
626 ") since Dir2 is aligned with Dir1. Mesh setup "
627 "needs investigation");
639 for (
int i = 0; i < 5; ++i)
641 m_faces[i]->Reset(curvedEdges, curvedFaces);
652 for (
int i = 0; i < 5; ++i)
675 order0 = *max_element(tmp.begin(), tmp.end());
681 order0 = *max_element(tmp.begin(), tmp.end());
690 order1 = *max_element(tmp.begin(), tmp.end());
698 order1 = *max_element(tmp.begin(), tmp.end());
702 tmp.push_back(order0);
703 tmp.push_back(order1);
708 int order2 = *max_element(tmp.begin(), tmp.end());
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
std::vector< StdRegions::Orientation > m_forient
std::vector< StdRegions::Orientation > m_eorient
bool m_setupState
Wether or not the setup routines have been run.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
int GetVid(int i) const
Get the ID of vertex i of this object.
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
int GetGlobalID(void) const
Get the ID of this object.
LibUtilities::ShapeType m_shapeType
Type of shape.
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
GeomState m_geomFactorsState
State of the geometric factors.
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
int m_coordim
Coordinate dimension of this geometry object.
int GetEid(int i) const
Get the ID of edge i of this object.
void SetUpFaceOrientation()
void SetUpEdgeOrientation()
void SetUpLocalVertices()
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
static const int kNtfaces
static const int kNqfaces
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
virtual void v_GenGeomFactors()
@ eGaussRadauMAlpha2Beta0
Gauss Radau pinned at x=-1, .
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModifiedPyr_C
Principle Modified Functions .
@ eModified_A
Principle Modified Functions .
static const NekDouble kNekZeroTol
std::unordered_map< int, CurveSharedPtr > CurveMap
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< SegGeom > SegGeomSharedPtr
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
@ ePtsFilled
Geometric information has been generated.
@ eDir1FwdDir2_Dir2FwdDir1
The above copyright notice and this permission notice shall be included.
scalarT< T > sqrt(scalarT< T > in)