48 namespace SpatialDomains
58 :
Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
81 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
96 for (
int i = 0; i <
kNedges; i++)
114 return m_xmap->PhysEvaluate(Lcoord, tmp);
131 int i, j, vmap[3] = {-1, -1, -1};
137 for (i = 0; i < 3; ++i)
139 rotPt.
Rotate((*face1[i]), dir, angle);
140 for (j = 0; j < 3; ++j)
142 if (rotPt.
dist(*face2[j]) < tol)
153 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
159 for (i = 0; i < 3; ++i)
161 cx += (*face2[i])(0) - (*face1[i])(0);
162 cy += (*face2[i])(1) - (*face1[i])(1);
163 cz += (*face2[i])(2) - (*face1[i])(2);
172 for (i = 0; i < 3; ++i)
177 for (j = 0; j < 3; ++j)
179 x1 = (*face2[j])(0) - cx;
180 y1 = (*face2[j])(1) - cy;
181 z1 = (*face2[j])(2) - cz;
182 if (
sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
183 (z1 - z) * (z1 - z)) < 1e-8)
192 if (vmap[1] == (vmap[0] + 1) % 3)
223 ASSERTL0(
false,
"Unable to determine triangle orientation");
243 if (
m_xmap->GetBasisNumModes(0) != 2 ||
244 m_xmap->GetBasisNumModes(1) != 2)
271 int nEdgeCoeffs =
m_xmap->GetTraceNcoeffs(0);
284 int N =
m_curve->m_points.size();
288 ASSERTL0(nEdgePts * (nEdgePts + 1) / 2 == N,
289 "NUMPOINTS should be a triangle number for"
290 " triangle curved face " +
295 for (i = 0; i < 3; ++i)
300 std::stringstream ss;
301 ss <<
"Curved vertex " << i <<
" of triangle " <<
m_globalID
302 <<
" is separated from expansion vertex by"
304 <<
" (dist = " << dist <<
")";
315 ASSERTL0(edgeCurve->m_points.size() == nEdgePts,
316 "Number of edge points does not correspond "
317 "to number of face points in triangle " +
320 const int offset = 3 + i * (nEdgePts - 2);
336 for (j = 0; j < nEdgePts - 2; ++j)
339 *(edgeCurve->m_points[j + 1]));
340 maxDist = dist > maxDist ? dist : maxDist;
345 for (j = 0; j < nEdgePts - 2; ++j)
348 *(edgeCurve->m_points[nEdgePts - 2 - j]));
349 maxDist = dist > maxDist ? dist : maxDist;
355 std::stringstream ss;
356 ss <<
"Curved edge " << i <<
" of triangle " <<
m_globalID
357 <<
" has a point separated from edge interior"
358 <<
" points by more than "
360 <<
" (maxdist = " << maxDist <<
")";
368 nEdgePts, LibUtilities::eGaussRadauMAlpha1Beta0);
374 max(nEdgePts * nEdgePts,
m_xmap->GetTotPoints()));
384 for (j = 0; j < N; ++j)
386 phys[j] = (
m_curve->m_points[j]->GetPtr())[i];
389 t->BwdTrans(phys, tmp);
393 P0, P1, tmp,
m_xmap->GetBasis(0)->GetPointsKey(),
394 m_xmap->GetBasis(1)->GetPointsKey(), phys);
402 int npts =
m_curve->m_points.size();
411 ASSERTL0(nEdgePts * nEdgePts == npts,
412 "NUMPOINTS should be a square number for"
419 "Number of edge points does not correspond to "
420 "number of face points in triangle " +
426 for (j = 0; j < npts; ++j)
428 tmp[j] = (
m_curve->m_points[j]->GetPtr())[i];
434 m_xmap->GetBasis(0)->GetPointsKey(),
435 m_xmap->GetBasis(1)->GetPointsKey(),
444 ASSERTL0(
false,
"Only 1D/2D points distributions "
457 nEdgeCoeffs =
m_edges[i]->GetXmap()->GetNcoeffs();
461 for (k = 0; k < nEdgeCoeffs; k++)
464 signArray[k] *
m_edges[i]->GetCoeffs(j)[k];
474 boost::ignore_unused(j);
476 return i == 0 ? 0 : 1;
482 CurveMap::iterator it = curvedFaces.find(
m_globalID);
484 if (it != curvedFaces.end())
489 for (
int i = 0; i < 3; ++i)
491 m_edges[i]->Reset(curvedEdges, curvedFaces);
502 for (
int i = 0; i < 3; ++i)
514 int order0 =
m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes();
516 max(order0, max(
m_edges[1]->
GetXmap()->GetBasis(0)->GetNumModes(),
#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.
virtual void v_GenGeomFactors() override
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 StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
virtual void v_FillGeom() override
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord) override
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
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...
virtual void v_Setup() override
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
@ 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)