46 namespace SpatialDomains
49 const unsigned int PrismGeom::VertexEdgeConnectivity[6][3] = {
50 {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7}, {4, 5, 8}, {6, 7, 8}};
51 const unsigned int PrismGeom::VertexFaceConnectivity[6][3] = {
52 {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4}, {1, 2, 4}, {2, 3, 4}};
53 const unsigned int PrismGeom::EdgeFaceConnectivity[9][2] = {
54 {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2}, {2, 3}, {3, 4}, {2, 4}};
55 const unsigned int PrismGeom::EdgeNormalToFaceVert[5][4] = {
56 {4, 5, 6, 7}, {1, 3, 8, -1}, {0, 2, 4, 7}, {1, 3, 8, -1}, {0, 2, 5, 6}};
58 PrismGeom::PrismGeom()
64 :
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
93 else if (faceidx == 1 || faceidx == 3)
120 if (
m_xmap->GetBasisNumModes(0) != 2 ||
121 m_xmap->GetBasisNumModes(1) != 2 ||
122 m_xmap->GetBasisNumModes(2) != 2)
132 const unsigned int faceVerts[3][4] = {
133 {0, 1, 2, 3}, {1, 2, 5, 4}, {0, 3, 5, 4}};
135 for (f = 0; f < 3; f++)
140 if (fabs((*
m_verts[faceVerts[f][0]])(i) -
141 (*
m_verts[faceVerts[f][1]])(i) +
142 (*
m_verts[faceVerts[f][2]])(i) -
143 (*
m_verts[faceVerts[f][3]])(i)) >
194 for (f = 1; f < 5; f++)
196 int nEdges =
m_faces[f]->GetNumEdges();
198 for (i = 0; i < 4; i++)
200 for (j = 0; j < nEdges; j++)
214 std::ostringstream errstrm;
215 errstrm <<
"Connected faces do not share an edge. Faces ";
222 std::ostringstream errstrm;
223 errstrm <<
"Connected faces share more than one edge. Faces ";
232 for (i = 0; i < 3; i++)
234 for (j = 0; j < 4; j++)
246 std::ostringstream errstrm;
247 errstrm <<
"Connected faces do not share an edge. Faces ";
254 std::ostringstream errstrm;
255 errstrm <<
"Connected faces share more than one edge. Faces ";
261 for (f = 1; f < 4; f++)
264 for (i = 0; i <
m_faces[f]->GetNumEdges(); i++)
266 for (j = 0; j <
m_faces[f + 1]->GetNumEdges(); j++)
280 std::ostringstream errstrm;
281 errstrm <<
"Connected faces do not share an edge. Faces ";
288 std::ostringstream errstrm;
289 errstrm <<
"Connected faces share more than one edge. Faces ";
298 for (i = 0; i < 4; i++)
300 for (j = 0; j < 4; j++)
313 std::ostringstream errstrm;
314 errstrm <<
"Connected faces do not share an edge. Faces ";
321 std::ostringstream errstrm;
322 errstrm <<
"Connected faces share more than one edge. Faces ";
347 std::ostringstream errstrm;
348 errstrm <<
"Connected edges do not share a vertex. Edges ";
349 errstrm <<
m_edges[0]->GetGlobalID() <<
", "
355 for (
int i = 1; i < 3; i++)
367 std::ostringstream errstrm;
368 errstrm <<
"Connected edges do not share a vertex. Edges ";
369 errstrm <<
m_edges[i]->GetGlobalID() <<
", "
370 <<
m_edges[i - 1]->GetGlobalID();
391 std::ostringstream errstrm;
392 errstrm <<
"Connected edges do not share a vertex. Edges ";
393 errstrm <<
m_edges[8]->GetGlobalID();
404 const unsigned int edgeVerts[
kNedges][2] = {
405 {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 5}, {3, 5}, {4, 5}};
421 ASSERTL0(
false,
"Could not find matching vertex for the edge");
449 unsigned int baseVertex;
472 unsigned int orientation;
478 elementAaxis_length = 0.0;
479 elementBaxis_length = 0.0;
480 faceAaxis_length = 0.0;
481 faceBaxis_length = 0.0;
486 baseVertex =
m_faces[f]->GetVid(0);
500 if (f == 1 || f == 3)
506 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
507 (*
m_verts[faceVerts[f][0]])[i];
508 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
509 (*
m_verts[faceVerts[f][0]])[i];
516 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
517 (*
m_verts[faceVerts[f][0]])[i];
518 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
519 (*
m_verts[faceVerts[f][1]])[i];
526 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
527 (*
m_verts[faceVerts[f][2]])[i];
528 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
529 (*
m_verts[faceVerts[f][0]])[i];
534 ASSERTL0(
false,
"Could not find matching vertex for the face");
543 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
544 (*
m_verts[faceVerts[f][0]])[i];
545 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
546 (*
m_verts[faceVerts[f][0]])[i];
553 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
554 (*
m_verts[faceVerts[f][0]])[i];
555 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
556 (*
m_verts[faceVerts[f][1]])[i];
563 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
564 (*
m_verts[faceVerts[f][3]])[i];
565 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
566 (*
m_verts[faceVerts[f][1]])[i];
573 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
574 (*
m_verts[faceVerts[f][3]])[i];
575 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
576 (*
m_verts[faceVerts[f][0]])[i];
581 ASSERTL0(
false,
"Could not find matching vertex for the face");
588 int v =
m_faces[f]->GetNumVerts() - 1;
594 elementAaxis_length += pow(elementAaxis[i], 2);
595 elementBaxis_length += pow(elementBaxis[i], 2);
596 faceAaxis_length += pow(faceAaxis[i], 2);
597 faceBaxis_length += pow(faceBaxis[i], 2);
600 elementAaxis_length =
sqrt(elementAaxis_length);
601 elementBaxis_length =
sqrt(elementBaxis_length);
602 faceAaxis_length =
sqrt(faceAaxis_length);
603 faceBaxis_length =
sqrt(faceBaxis_length);
609 dotproduct1 += elementAaxis[i] * faceAaxis[i];
617 if (fabs(elementAaxis_length * faceAaxis_length - fabs(dotproduct1)) <
622 if (dotproduct1 < 0.0)
630 dotproduct2 += elementBaxis[i] * faceBaxis[i];
633 ASSERTL1(fabs(fabs(dotproduct2 / elementBaxis_length /
636 "These vectors should be parallel");
640 if (dotproduct2 < 0.0)
655 dotproduct1 += elementAaxis[i] * faceBaxis[i];
659 ASSERTL1(fabs(fabs(dotproduct1) / elementAaxis_length /
662 "These vectors should be parallel");
666 if (dotproduct1 < 0.0)
675 dotproduct2 += elementBaxis[i] * faceAaxis[i];
678 ASSERTL1(fabs(fabs(dotproduct2) / elementBaxis_length /
681 "These vectors should be parallel");
683 if (dotproduct2 < 0.0)
689 orientation = orientation + 5;
691 if ((f == 1) || (f == 3))
695 "Orientation of triangular face (id = " +
697 ") is inconsistent with face " +
698 boost::lexical_cast<string>(f) +
699 " of prism element (id = " +
701 ") since Dir2 is aligned with Dir1. Mesh setup "
702 "needs investigation");
714 for (
int i = 0; i < 5; ++i)
716 m_faces[i]->Reset(curvedEdges, curvedFaces);
724 for (
int i = 0; i < 5; ++i)
747 order0 = *max_element(tmp.begin(), tmp.end());
753 order0 = *max_element(tmp.begin(), tmp.end());
762 order1 = *max_element(tmp.begin(), tmp.end());
770 order1 = *max_element(tmp.begin(), tmp.end());
774 tmp.push_back(order0);
775 tmp.push_back(order1);
780 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.
std::vector< StdRegions::Orientation > m_forient
virtual void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
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 void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces) override
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
static const unsigned int EdgeFaceConnectivity[9][2]
virtual int v_GetVertexFaceMap(const int i, const int j) const override
Returns the standard element face IDs that are connected to a given vertex.
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 override
Returns the element coordinate direction corresponding to a given face coordinate direction.
void SetUpEdgeOrientation()
void SetUpFaceOrientation()
virtual void v_GenGeomFactors() override
static const unsigned int EdgeNormalToFaceVert[5][4]
virtual int v_GetEdgeFaceMap(const int i, const int j) const override
Returns the standard element edge IDs that are connected to a given face.
virtual int v_GetEdgeNormalToFaceVert(const int i, const int j) const override
Returns the standard lement edge IDs that are normal to a given face vertex.
static const unsigned int VertexFaceConnectivity[6][3]
virtual void v_Setup() override
static const int kNqfaces
void SetUpLocalVertices()
static const unsigned int VertexEdgeConnectivity[6][3]
virtual int v_GetVertexEdgeMap(const int i, const int j) const override
Returns the standard element edge IDs that are connected to a given vertex.
static const int kNtfaces
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModified_B
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)