48 namespace SpatialDomains
56 TriGeom::TriGeom(
const int id,
59 :
Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
81 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
97 for (
int i = 0; i <
kNedges; i++)
115 return m_xmap->PhysEvaluate(Lcoord, tmp);
119 const TriGeom &face2,
bool doRot,
int dir,
123 doRot, dir, angle, tol);
130 int i, j, vmap[3] = {-1, -1, -1};
136 for (i = 0; i < 3; ++i)
138 rotPt.
Rotate((*face1[i]), dir, angle);
139 for (j = 0; j < 3; ++j)
141 if (rotPt.
dist(*face2[j]) < tol)
152 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
158 for (i = 0; i < 3; ++i)
160 cx += (*face2[i])(0) - (*face1[i])(0);
161 cy += (*face2[i])(1) - (*face1[i])(1);
162 cz += (*face2[i])(2) - (*face1[i])(2);
171 for (i = 0; i < 3; ++i)
176 for (j = 0; j < 3; ++j)
178 x1 = (*face2[j])(0) - cx;
179 y1 = (*face2[j])(1) - cy;
180 z1 = (*face2[j])(2) - cz;
181 if (
sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
182 (z1 - z) * (z1 - z)) < 1e-8)
191 if (vmap[1] == (vmap[0] + 1) % 3)
222 ASSERTL0(
false,
"Unable to determine triangle orientation");
242 if (
m_xmap->GetBasisNumModes(0) != 2 ||
243 m_xmap->GetBasisNumModes(1) != 2)
270 int nEdgeCoeffs =
m_xmap->GetTraceNcoeffs(0);
283 int N =
m_curve->m_points.size();
287 ASSERTL0(nEdgePts * (nEdgePts + 1) / 2 == N,
288 "NUMPOINTS should be a triangle number for"
289 " triangle curved face " +
294 for (i = 0; i < 3; ++i)
299 std::stringstream ss;
300 ss <<
"Curved vertex " << i <<
" of triangle " <<
m_globalID
301 <<
" is separated from expansion vertex by"
303 <<
" (dist = " << dist <<
")";
314 ASSERTL0(edgeCurve->m_points.size() == nEdgePts,
315 "Number of edge points does not correspond "
316 "to number of face points in triangle " +
319 const int offset = 3 + i * (nEdgePts - 2);
334 for (j = 0; j < nEdgePts - 2; ++j)
337 *(edgeCurve->m_points[j + 1]));
338 maxDist = dist > maxDist ? dist : maxDist;
343 for (j = 0; j < nEdgePts - 2; ++j)
346 *(edgeCurve->m_points[nEdgePts - 2 - j]));
347 maxDist = dist > maxDist ? dist : maxDist;
353 std::stringstream ss;
354 ss <<
"Curved edge " << i <<
" of triangle " <<
m_globalID
355 <<
" has a point separated from edge interior"
356 <<
" points by more than "
358 <<
" (maxdist = " << maxDist <<
")";
372 max(nEdgePts * nEdgePts,
m_xmap->GetTotPoints()));
382 for (j = 0; j < N; ++j)
384 phys[j] = (
m_curve->m_points[j]->GetPtr())[i];
387 t->BwdTrans(phys, tmp);
393 m_xmap->GetBasis(0)->GetPointsKey(),
394 m_xmap->GetBasis(1)->GetPointsKey(),
403 int npts =
m_curve->m_points.size();
412 ASSERTL0(nEdgePts * nEdgePts == npts,
413 "NUMPOINTS should be a square number for"
420 "Number of edge points does not correspond to "
421 "number of face points in triangle " +
427 for (j = 0; j < npts; ++j)
429 tmp[j] = (
m_curve->m_points[j]->GetPtr())[i];
437 m_xmap->GetBasis(0)->GetPointsKey(),
438 m_xmap->GetBasis(1)->GetPointsKey(),
448 "Only 1D/2D points distributions "
461 nEdgeCoeffs =
m_edges[i]->GetXmap()->GetNcoeffs();
465 for (k = 0; k < nEdgeCoeffs; k++)
468 signArray[k] *
m_edges[i]->GetCoeffs(j)[k];
478 boost::ignore_unused(j);
486 CurveMap::iterator it = curvedFaces.find(
m_globalID);
488 if (it != curvedFaces.end())
493 for (
int i = 0; i < 3; ++i)
495 m_edges[i]->Reset(curvedEdges, curvedFaces);
506 for (
int i = 0; i < 3; ++i)
518 int order0 =
m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes();
519 int order1 = max(order0,
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#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_eorient
bool m_setupState
Wether or not the setup routines have been run.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
GeomState m_state
Enumeration to dictate whether coefficients are filled.
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
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.
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
NekDouble dist(PointGeom &a)
return distance between this and input a
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
static const int kNedges
Get the orientation of face1.
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
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()
virtual void v_FillGeom()
PointsManagerT & PointsManager(void)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
@ eModified_B
Principle Modified Functions .
@ eOrtho_A
Principle Orthogonal Functions .
@ eOrtho_B
Principle Orthogonal Functions .
@ eModified_A
Principle Modified Functions .
static const NekDouble kVertexTheSameDouble
std::vector< PointGeomSharedPtr > PointGeomVector
std::shared_ptr< Curve > CurveSharedPtr
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< StdNodalTriExp > StdNodalTriExpSharedPtr
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1FwdDir1_Dir2FwdDir2
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir1_Dir2FwdDir2
@ eDir1FwdDir2_Dir2FwdDir1
@ eDir1FwdDir2_Dir2BwdDir1
The above copyright notice and this permission notice shall be included.
scalarT< T > sqrt(scalarT< T > in)