46 namespace SpatialDomains
49 const unsigned int HexGeom::VertexEdgeConnectivity[8][3] = {
50 {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7},
51 {4, 8, 11}, {5, 8, 9}, {6, 9, 10}, {7, 10, 11}};
52 const unsigned int HexGeom::VertexFaceConnectivity[8][3] = {
53 {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4},
54 {1, 4, 5}, {1, 2, 5}, {2, 3, 5}, {3, 4, 5}};
55 const unsigned int HexGeom::EdgeFaceConnectivity[12][2] = {
56 {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2},
57 {2, 3}, {3, 4}, {1, 5}, {2, 5}, {3, 5}, {4, 5}};
65 :
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
104 if (
m_xmap->GetBasisNumModes(0) != 2 ||
105 m_xmap->GetBasisNumModes(1) != 2 ||
106 m_xmap->GetBasisNumModes(2) != 2)
128 if (fabs((*
m_verts[faceVerts[f][0]])(i) -
129 (*
m_verts[faceVerts[f][1]])(i) +
130 (*
m_verts[faceVerts[f][2]])(i) -
131 (*
m_verts[faceVerts[f][3]])(i)) >
169 if (faceidx == 0 || faceidx == 5)
173 else if (faceidx == 1 || faceidx == 3)
193 for (f = 1; f < 5; f++)
196 for (i = 0; i < 4; i++)
198 for (j = 0; j < 4; j++)
202 edge = std::dynamic_pointer_cast<SegGeom>(
212 std::ostringstream errstrm;
213 errstrm <<
"Connected faces do not share an edge. Faces ";
220 std::ostringstream errstrm;
221 errstrm <<
"Connected faces share more than one edge. Faces ";
230 for (i = 0; i < 4; i++)
232 for (j = 0; j < 4; j++)
236 edge = std::dynamic_pointer_cast<SegGeom>(
245 std::ostringstream errstrm;
246 errstrm <<
"Connected faces do not share an edge. Faces ";
253 std::ostringstream errstrm;
254 errstrm <<
"Connected faces share more than one edge. Faces ";
259 for (f = 1; f < 4; f++)
262 for (i = 0; i < 4; i++)
264 for (j = 0; j < 4; j++)
268 edge = std::dynamic_pointer_cast<SegGeom>(
278 std::ostringstream errstrm;
279 errstrm <<
"Connected faces do not share an edge. Faces ";
286 std::ostringstream errstrm;
287 errstrm <<
"Connected faces share more than one edge. Faces ";
295 for (f = 1; f < 5; f++)
298 for (i = 0; i < 4; i++)
300 for (j = 0; j < 4; j++)
304 edge = std::dynamic_pointer_cast<SegGeom>(
314 std::ostringstream errstrm;
315 errstrm <<
"Connected faces do not share an edge. Faces ";
322 std::ostringstream errstrm;
323 errstrm <<
"Connected faces share more than one edge. Faces ";
348 std::ostringstream errstrm;
349 errstrm <<
"Connected edges do not share a vertex. Edges ";
350 errstrm <<
m_edges[0]->GetGlobalID() <<
", "
357 for (i = 1; i < 3; i++)
369 std::ostringstream errstrm;
370 errstrm <<
"Connected edges do not share a vertex. Edges ";
371 errstrm <<
m_edges[i]->GetGlobalID() <<
", "
372 <<
m_edges[i - 1]->GetGlobalID();
393 std::ostringstream errstrm;
394 errstrm <<
"Connected edges do not share a vertex. Edges ";
395 errstrm <<
m_edges[8]->GetGlobalID() <<
", "
401 for (i = 9; i < 11; i++)
413 std::ostringstream errstrm;
414 errstrm <<
"Connected edges do not share a vertex. Edges ";
415 errstrm <<
m_edges[i]->GetGlobalID() <<
", "
416 <<
m_edges[i - 1]->GetGlobalID();
445 unsigned int baseVertex;
467 unsigned int orientation;
473 elementAaxis_length = 0.0;
474 elementBaxis_length = 0.0;
475 faceAaxis_length = 0.0;
476 faceBaxis_length = 0.0;
481 baseVertex =
m_faces[f]->GetVid(0);
497 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
498 (*
m_verts[faceVerts[f][0]])[i];
499 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
500 (*
m_verts[faceVerts[f][0]])[i];
507 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
508 (*
m_verts[faceVerts[f][0]])[i];
509 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
510 (*
m_verts[faceVerts[f][1]])[i];
517 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
518 (*
m_verts[faceVerts[f][3]])[i];
519 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
520 (*
m_verts[faceVerts[f][1]])[i];
527 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
528 (*
m_verts[faceVerts[f][3]])[i];
529 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
530 (*
m_verts[faceVerts[f][0]])[i];
535 ASSERTL0(
false,
"Could not find matching vertex for the face");
547 elementAaxis_length += pow(elementAaxis[i], 2);
548 elementBaxis_length += pow(elementBaxis[i], 2);
549 faceAaxis_length += pow(faceAaxis[i], 2);
550 faceBaxis_length += pow(faceBaxis[i], 2);
553 elementAaxis_length =
sqrt(elementAaxis_length);
554 elementBaxis_length =
sqrt(elementBaxis_length);
555 faceAaxis_length =
sqrt(faceAaxis_length);
556 faceBaxis_length =
sqrt(faceBaxis_length);
562 dotproduct1 += elementAaxis[i] * faceAaxis[i];
565 NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
575 if (dotproduct1 < 0.0)
583 dotproduct2 += elementBaxis[i] * faceBaxis[i];
586 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
590 "These vectors should be parallel");
594 if (dotproduct2 < 0.0)
609 dotproduct1 += elementAaxis[i] * faceBaxis[i];
612 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
616 "These vectors should be parallel");
620 if (dotproduct1 < 0.0)
629 dotproduct2 += elementBaxis[i] * faceAaxis[i];
632 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
636 "These vectors should be parallel");
638 if (dotproduct2 < 0.0)
644 orientation = orientation + 5;
656 const unsigned int edgeVerts[
kNedges][2] = {{0, 1},
682 ASSERTL0(
false,
"Could not find matching vertex for the edge");
691 for (
int i = 0; i < 6; ++i)
693 m_faces[i]->Reset(curvedEdges, curvedFaces);
704 for (
int i = 0; i < 6; ++i)
746 int order0 = *max_element(tmp1.begin(), tmp1.end());
772 int order1 = *max_element(tmp1.begin(), tmp1.end());
798 int order2 = *max_element(tmp1.begin(), tmp1.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.
virtual int v_GetVertexFaceMap(const int i, const int j) const
Returns the standard element face IDs that are connected to a given vertex.
static const unsigned int EdgeFaceConnectivity[12][2]
void SetUpFaceOrientation()
static const unsigned int VertexEdgeConnectivity[8][3]
virtual int v_GetEdgeFaceMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given face.
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 unsigned int VertexFaceConnectivity[8][3]
virtual int v_GetVertexEdgeMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given vertex.
void SetUpEdgeOrientation()
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
void SetUpLocalVertices()
virtual void v_GenGeomFactors()
static const int kNqfaces
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
static const int kNtfaces
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModified_A
Principle Modified Functions .
static const NekDouble kNekZeroTol
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
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
@ ePtsFilled
Geometric information has been generated.
The above copyright notice and this permission notice shall be included.
scalarT< T > sqrt(scalarT< T > in)