45namespace SpatialDomains
 
   48    {0, 2, 3}, {0, 1, 4}, {1, 2, 5}, {3, 4, 5}};
 
   50    {0, 1, 3}, {0, 1, 2}, {0, 2, 3}, {1, 2, 3}};
 
   52    {0, 1}, {0, 2}, {0, 3}, {1, 3}, {1, 2}, {2, 3}};
 
   54    {3, 4, 5}, {1, 2, 5}, {0, 2, 3}, {0, 1, 4}};
 
   62    : 
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
 
   90    else if (faceidx == 1)
 
  133        std::ostringstream errstrm;
 
  134        errstrm << 
"Local edge 0 (eid=" << 
m_faces[0]->GetEid(0);
 
  135        errstrm << 
") on face " << 
m_faces[0]->GetGlobalID();
 
  136        errstrm << 
" must be the same as local edge 0 (eid=" 
  138        errstrm << 
") on face " << 
m_faces[1]->GetGlobalID();
 
  143    for (faceConnected = 1; faceConnected < 4; faceConnected++)
 
  146        for (i = 0; i < 3; i++)
 
  158            std::ostringstream errstrm;
 
  159            errstrm << 
"Face 0 does not share an edge with first edge of " 
  160                       "adjacent face. Faces ";
 
  167            std::ostringstream errstrm;
 
  168            errstrm << 
"Connected faces share more than one edge. Faces ";
 
  177    for (i = 0; i < 3; i++) 
 
  179        for (j = 0; j < 3; j++)
 
  191        std::ostringstream errstrm;
 
  192        errstrm << 
"Connected faces do not share an edge. Faces ";
 
  199        std::ostringstream errstrm;
 
  200        errstrm << 
"Connected faces share more than one edge. Faces ";
 
  206    for (faceConnected = 1; faceConnected < 3; faceConnected++)
 
  209        for (i = 0; i < 3; i++)
 
  211            for (j = 0; j < 3; j++)
 
  216                    edge = dynamic_pointer_cast<SegGeom>(
 
  226            std::ostringstream errstrm;
 
  227            errstrm << 
"Connected faces do not share an edge. Faces ";
 
  234            std::ostringstream errstrm;
 
  235            errstrm << 
"Connected faces share more than one edge. Faces ";
 
  261        std::ostringstream errstrm;
 
  262        errstrm << 
"Connected edges do not share a vertex. Edges ";
 
  263        errstrm << 
m_edges[0]->GetGlobalID() << 
", " 
  269    for (
int i = 1; i < 2; i++)
 
  281            std::ostringstream errstrm;
 
  282            errstrm << 
"Connected edges do not share a vertex. Edges ";
 
  283            errstrm << 
m_edges[i]->GetGlobalID() << 
", " 
  284                    << 
m_edges[i - 1]->GetGlobalID();
 
  301    for (
int i = 4; i < 6; ++i)
 
  313        std::ostringstream errstrm;
 
  314        errstrm << 
"Connected edges do not share a vertex. Edges ";
 
  315        errstrm << 
m_edges[3]->GetGlobalID() << 
", " 
  327    const unsigned int edgeVerts[
kNedges][2] = {{0, 1}, {1, 2}, {0, 2},
 
  328                                                {0, 3}, {1, 3}, {2, 3}};
 
  344            ASSERTL0(
false, 
"Could not find matching vertex for the edge");
 
  373    unsigned int baseVertex;
 
  386        {0, 1, 2}, {0, 1, 3}, {1, 2, 3}, {0, 2, 3}};
 
  391    unsigned int orientation;
 
  397        elementAaxis_length = 0.0;
 
  398        elementBaxis_length = 0.0;
 
  399        faceAaxis_length    = 0.0;
 
  400        faceBaxis_length    = 0.0;
 
  405        baseVertex = 
m_faces[f]->GetVid(0);
 
  423                elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
 
  424                                  (*
m_verts[faceVerts[f][0]])[i];
 
  425                elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
 
  426                                  (*
m_verts[faceVerts[f][0]])[i];
 
  433                elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
 
  434                                  (*
m_verts[faceVerts[f][0]])[i];
 
  435                elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
 
  436                                  (*
m_verts[faceVerts[f][1]])[i];
 
  443                elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
 
  444                                  (*
m_verts[faceVerts[f][2]])[i];
 
  445                elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
 
  446                                  (*
m_verts[faceVerts[f][0]])[i];
 
  451            ASSERTL0(
false, 
"Could not find matching vertex for the face");
 
  463            elementAaxis_length += pow(elementAaxis[i], 2);
 
  464            elementBaxis_length += pow(elementBaxis[i], 2);
 
  465            faceAaxis_length += pow(faceAaxis[i], 2);
 
  466            faceBaxis_length += pow(faceBaxis[i], 2);
 
  469        elementAaxis_length = 
sqrt(elementAaxis_length);
 
  470        elementBaxis_length = 
sqrt(elementBaxis_length);
 
  471        faceAaxis_length    = 
sqrt(faceAaxis_length);
 
  472        faceBaxis_length    = 
sqrt(faceBaxis_length);
 
  478            dotproduct1 += elementAaxis[i] * faceAaxis[i];
 
  482            fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
 
  492            if (dotproduct1 < 0.0)
 
  500                dotproduct2 += elementBaxis[i] * faceBaxis[i];
 
  503            norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
 
  507                     "These vectors should be parallel");
 
  511            if (dotproduct2 < 0.0)
 
  526                dotproduct1 += elementAaxis[i] * faceBaxis[i];
 
  529            norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
 
  533                     "These vectors should be parallel");
 
  537            if (dotproduct1 < 0.0)
 
  546                dotproduct2 += elementBaxis[i] * faceAaxis[i];
 
  549            norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
 
  553                     "These vectors should be parallel");
 
  555            if (dotproduct2 < 0.0)
 
  561        orientation = orientation + 5;
 
  564                 "Orientation of triangular face (id = " +
 
  566                     ") is inconsistent with face " +
 
  567                     boost::lexical_cast<string>(f) + 
" of tet element (id = " +
 
  569                     ") since Dir2 is aligned with Dir1. Mesh setup " 
  570                     "needs investigation");
 
  581    for (
int i = 0; i < 4; ++i)
 
  583        m_faces[i]->Reset(curvedEdges, curvedFaces);
 
  591        for (
int i = 0; i < 4; ++i)
 
  619        if (
m_xmap->GetBasisNumModes(0) != 2 ||
 
  620            m_xmap->GetBasisNumModes(1) != 2 ||
 
  621            m_xmap->GetBasisNumModes(2) != 2)
 
  630            for (
int i = 0; i < 3; ++i)
 
  660    int order0 = *max_element(tmp.begin(), tmp.end());
 
  663    tmp.push_back(order0);
 
  666    int order1 = *max_element(tmp.begin(), tmp.end());
 
  669    tmp.push_back(order0);
 
  670    tmp.push_back(order1);
 
  674    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_CalculateInverseIsoParam() override
 
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.
 
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
 
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_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_GetVertexFaceMap(const int i, const int j) const override
Returns the standard element face IDs that are connected to a given vertex.
 
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 unsigned int VertexEdgeConnectivity[4][3]
 
static const unsigned int EdgeNormalToFaceVert[4][3]
 
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.
 
void SetUpFaceOrientation()
 
void SetUpEdgeOrientation()
 
void SetUpLocalVertices()
 
virtual void v_GenGeomFactors() override
 
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 int kNqfaces
 
static const int kNtfaces
 
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.
 
static const unsigned int VertexFaceConnectivity[4][3]
 
virtual void v_Setup() override
 
static const unsigned int EdgeFaceConnectivity[6][2]
 
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
 
@ eModified_B
Principle Modified Functions .
 
@ eModified_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
 
@ ePtsFilled
Geometric information has been generated.
 
std::shared_ptr< TriGeom > TriGeomSharedPtr
 
@ eDir1FwdDir2_Dir2FwdDir1
 
The above copyright notice and this permission notice shall be included.
 
scalarT< T > sqrt(scalarT< T > in)