Nektar++
Loading...
Searching...
No Matches
Geometry.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: Geometry.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: This file contains the base class implementation for the
32// Geometry class.
33//
34//
35////////////////////////////////////////////////////////////////////////////////
36
40
42{
43
44// static class property
46
47/**
48 * @brief Default constructor.
49 */
51 : m_coordim(0), m_geomFactorsState(eNotFilled), m_state(eNotFilled),
52 m_setupState(false), m_shapeType(LibUtilities::eNoShapeType),
53 m_globalID(-1), m_straightEdge(0)
54{
55}
56
57/**
58 * @brief Constructor when supplied a coordinate dimension.
59 */
60Geometry::Geometry(const int coordim)
61 : m_coordim(coordim), m_geomFactorsState(eNotFilled), m_state(eNotFilled),
62 m_setupState(false), m_shapeType(LibUtilities::eNoShapeType),
63 m_globalID(-1), m_straightEdge(0)
64{
65}
66
67/**
68 * @brief Check to see if a geometric factor has already been created that
69 * contains the same regular information.
70 *
71 * The principle behind this is that many regular (i.e. constant Jacobian)
72 * elements have identicial geometric factors. Memory may therefore be reduced
73 * by storing only the unique factors.
74 *
75 * @param geomFactor The GeomFactor to check.
76 *
77 * @return Either the cached GeomFactor or @p geomFactor.
78 *
79 * @todo Currently this method is disabled since the lookup is very expensive.
80 */
82 GeomFactorsSharedPtr geomFactor)
83{
84 GeomFactorsSharedPtr returnval = geomFactor;
85
86 return returnval;
87}
88
89bool SortByGlobalId(const Geometry *&lhs, const Geometry *&rhs)
90{
91 return lhs->GetGlobalID() < rhs->GetGlobalID();
92}
93
94bool GlobalIdEquality(const Geometry *&lhs, const Geometry *&rhs)
95{
96 return lhs->GetGlobalID() == rhs->GetGlobalID();
97}
98
99/**
100 * @brief Get the ID of vertex @p i of this object.
101 */
102int Geometry::v_GetVid(int i) const
103{
104 return GetVertex(i)->GetGlobalID();
105}
106
107/**
108 * @brief Get the ID of edge @p i of this object.
109 */
110int Geometry::GetEid(int i) const
111{
112 return GetEdge(i)->GetGlobalID();
113}
114
115/**
116 * @brief Get the ID of face @p i of this object.
117 */
118int Geometry::GetFid(int i) const
119{
120 return GetFace(i)->GetGlobalID();
121}
122
123/**
124 * @copydoc Geometry::GetVertex()
125 */
126PointGeom *Geometry::v_GetVertex([[maybe_unused]] int i) const
127{
129 "This function is only valid for shape type geometries");
130 return nullptr;
131}
132
133/**
134 * @copydoc Geometry::GetEdge()
135 */
136Geometry1D *Geometry::v_GetEdge([[maybe_unused]] int i) const
137{
139 "This function is only valid for shape type geometries");
140 return nullptr;
141}
142
143/**
144 * @copydoc Geometry::GetFace()
145 */
146Geometry2D *Geometry::v_GetFace([[maybe_unused]] int i) const
147{
149 "This function is only valid for shape type geometries");
150 return nullptr;
151}
152
153/**
154 * @copydoc Geometry::GetNumVerts()
155 */
157{
159 "This function is only valid for shape type geometries");
160 return 0;
161}
162
163/**
164 * @copydoc Geometry::GetEorient()
165 */
167 [[maybe_unused]] const int i) const
168{
170 "This function is not valid for this geometry.");
172}
173
174/**
175 * @copydoc Geometry::GetForient()
176 */
178 [[maybe_unused]] const int i) const
179{
181 "This function is not valid for this geometry.");
182 return StdRegions::eFwd;
183}
184
185/**
186 * @copydoc Geometry::GetNumEdges()
187 */
189{
190 return 0;
191}
192
193/**
194 * @copydoc Geometry::GetNumFaces()
195 */
197{
198 return 0;
199}
200
201/**
202 * @copydoc Geometry::GetShapeDim()
203 */
205{
207 "This function is only valid for shape type geometries");
208 return 0;
209}
210
212 [[maybe_unused]] const Array<OneD, const NekDouble> &gloCoord)
213{
214 return 0;
215}
216
218{
220 "This function is only valid for shape type geometries");
221}
222
223/**
224 * @copydoc Geometry::GetXmap()
225 */
230
231/**
232 * @copydoc Geometry::ContainsPoint(
233 * const Array<OneD, const NekDouble> &, Array<OneD, NekDouble> &,
234 * NekDouble, NekDouble&)
235 * dist is assigned value for curved elements
236 */
238 Array<OneD, NekDouble> &locCoord, NekDouble tol,
239 NekDouble &dist)
240{
241 int inside = PreliminaryCheck(gloCoord);
242 if (inside == -1)
243 {
244 dist = std::numeric_limits<double>::max();
245 return false;
246 }
247 dist = GetLocCoords(gloCoord, locCoord);
248 if (inside == 1)
249 {
250 dist = 0.;
251 return true;
252 }
253 else
254 {
256 m_xmap->LocCoordToLocCollapsed(locCoord, eta);
257 if (ClampLocCoords(eta, tol))
258 {
259 if (GetMetricInfo()->GetGtype() == eRegular)
260 {
261 dist = std::numeric_limits<double>::max();
262 }
263 return false;
264 }
265 return 3 != m_coordim ||
268 dist <= tol;
269 }
270}
271
273 [[maybe_unused]] const Array<OneD, const NekDouble> &xs,
274 [[maybe_unused]] Array<OneD, NekDouble> &xi)
275{
277 "This function has not been defined for this geometry");
278 return false;
279}
280
281/**
282 * @copydoc Geometry::GetVertexEdgeMap()
283 */
284int Geometry::v_GetVertexEdgeMap([[maybe_unused]] const int i,
285 [[maybe_unused]] const int j) const
286{
288 "This function has not been defined for this geometry");
289 return 0;
290}
291
292/**
293 * @copydoc Geometry::GetVertexFaceMap()
294 */
295int Geometry::v_GetVertexFaceMap([[maybe_unused]] const int i,
296 [[maybe_unused]] const int j) const
297{
299 "This function has not been defined for this geometry");
300 return 0;
301}
302
303/**
304 * @copydoc Geometry::GetEdgeFaceMap()
305 */
306int Geometry::v_GetEdgeFaceMap([[maybe_unused]] const int i,
307 [[maybe_unused]] const int j) const
308{
310 "This function has not been defined for this geometry");
311 return 0;
312}
313
314/**
315 * @copydoc Geometry::GetEdgeNormalToFaceVert()
316 */
317int Geometry::v_GetEdgeNormalToFaceVert([[maybe_unused]] const int i,
318 [[maybe_unused]] const int j) const
319{
321 "This function has not been defined for this geometry");
322 return 0;
323}
324
325/**
326 * @copydoc Geometry::GetDir()
327 */
328int Geometry::v_GetDir([[maybe_unused]] const int i,
329 [[maybe_unused]] const int j) const
330{
332 "This function has not been defined for this geometry");
333 return 0;
334}
335
336/**
337 * @copydoc Geometry::GetCoord()
338 */
340 [[maybe_unused]] const int i,
341 [[maybe_unused]] const Array<OneD, const NekDouble> &Lcoord)
342{
344 "This function is only valid for expansion type geometries");
345 return 0.0;
346}
347
348/**
349 * @copydoc Geometry::GetLocCoords()
350 */
352 [[maybe_unused]] const Array<OneD, const NekDouble> &coords,
353 [[maybe_unused]] Array<OneD, NekDouble> &Lcoords)
354{
356 "This function is only valid for expansion type geometries");
357 return 0.0;
358}
359
360/**
361 * @copydoc Geometry::FillGeom()
362 */
364{
366 "This function is only valid for expansion type geometries");
367}
368
369/**
370 * @copydoc Geometry::Reset()
371 */
372void Geometry::v_Reset([[maybe_unused]] CurveMap &curvedEdges,
373 [[maybe_unused]] CurveMap &curvedFaces)
374{
375
376 // Reset state
379
380 // Junk geometric factors
382}
383
385{
387 "This function is only valid for expansion type geometries");
388}
389
390/**
391 * @brief Generates the bounding box for the element.
392 *
393 * For regular elements, the vertices are sufficient to define the extent of
394 * the bounding box. For non-regular elements, the extremes of the quadrature
395 * point coordinates are used. A 10% margin is added around this computed
396 * region to account for convex hull elements where the true extent of the
397 * element may extend slightly beyond the quadrature points.
398 */
399std::array<NekDouble, 6> Geometry::GetBoundingBox()
400{
401 if (m_boundingBox.size() == 6)
402 {
403 return {{m_boundingBox[0], m_boundingBox[1], m_boundingBox[2],
405 }
406 // NekDouble minx, miny, minz, maxx, maxy, maxz;
408
409 // Always get vertexes min/max
410 PointGeom *p = GetVertex(0);
411 Array<OneD, NekDouble> x(3, 0.0);
412 p->GetCoords(x[0], x[1], x[2]);
413 for (int j = 0; j < 3; ++j)
414 {
415 min[j] = x[j];
416 max[j] = x[j];
417 }
418 for (int i = 1; i < GetNumVerts(); ++i)
419 {
420 p = GetVertex(i);
421 p->GetCoords(x[0], x[1], x[2]);
422 for (int j = 0; j < 3; ++j)
423 {
424 min[j] = (x[j] < min[j] ? x[j] : min[j]);
425 max[j] = (x[j] > max[j] ? x[j] : max[j]);
426 }
427 }
428 // If element is deformed loop over quadrature points
430 if (GetGeomFactors()->GetGtype() != eRegular)
431 {
432 marginFactor = 0.1;
433 const int nq = GetXmap()->GetTotPoints();
435 for (int j = 0; j < 3; ++j)
436 {
437 xvec[j] = Array<OneD, NekDouble>(nq, 0.0);
438 }
439 for (int j = 0; j < GetCoordim(); ++j)
440 {
441 GetXmap()->BwdTrans(m_coeffs[j], xvec[j]);
442 }
443 for (int j = 0; j < 3; ++j)
444 {
445 for (int i = 0; i < nq; ++i)
446 {
447 min[j] = (xvec[j][i] < min[j] ? xvec[j][i] : min[j]);
448 max[j] = (xvec[j][i] > max[j] ? xvec[j][i] : max[j]);
449 }
450 }
451 }
452 // Add margin to bounding box, in order to
453 // return the nearest element
454 for (int j = 0; j < 3; ++j)
455 {
456 NekDouble margin =
457 marginFactor * (max[j] - min[j]) + NekConstants::kFindDistanceMin;
458 min[j] -= margin;
459 max[j] += margin;
460 }
461
462 // save bounding box
464 for (int j = 0; j < 3; ++j)
465 {
466 m_boundingBox[j] = min[j];
467 m_boundingBox[j + 3] = max[j];
468 }
469 // Return bounding box
470 return {{min[0], min[1], min[2], max[0], max[1], max[2]}};
471}
472
474{
475 m_boundingBox = {};
476}
477
478/**
479 * @brief A fast and robust check if a given global coord is outside of a
480 * deformed element. For regular elements, this check is unnecessary.
481 *
482 * @param coords Input Cartesian global coordinates
483 *
484 * @return 1 is inside of the element.
485 * 0 maybe inside
486 * -1 outside of the element
487 */
489{
490 // bounding box check
491 if (!MinMaxCheck(gloCoord))
492 {
493 return -1;
494 }
495
496 // regular element check
497 if (GetMetricInfo()->GetGtype() == eRegular)
498 {
499 return 0;
500 }
501
502 // All left check for straight edges/plane surfaces
503 return v_AllLeftCheck(gloCoord);
504}
505
506/**
507 * @brief Check if given global coord is within the BoundingBox of the element.
508 *
509 * @param coords Input Cartesian global coordinates
510 *
511 * @return True if within distance or False otherwise.
512 */
514{
515 // Validation checks
516 ASSERTL1(gloCoord.size() >= m_coordim,
517 "Expects number of global coordinates supplied to be greater than "
518 "or equal to the mesh dimension.");
519
520 std::array<NekDouble, 6> minMax = GetBoundingBox();
521 for (int i = 0; i < m_coordim; ++i)
522 {
523 if ((gloCoord[i] < minMax[i]) || (gloCoord[i] > minMax[i + 3]))
524 {
525 return false;
526 }
527 }
528 return true;
529}
530
531/**
532 * @brief Clamp local coords to be within standard regions [-1, 1]^dim.
533 *
534 * @param Lcoords Corresponding local coordinates
535 */
537{
538 // Validation checks
539 ASSERTL1(locCoord.size() >= GetShapeDim(),
540 "Expects local coordinates to be same or "
541 "larger than shape dimension.");
542
543 // If out of range clamp locCoord to be within [-1,1]^dim
544 // since any larger value will be very oscillatory if
545 // called by 'returnNearestElmt' option in
546 // ExpList::GetExpIndex
547 bool clamp = false;
548 for (int i = 0; i < GetShapeDim(); ++i)
549 {
550 if (!std::isfinite(locCoord[i]))
551 {
552 locCoord[i] = 0.;
553 clamp = true;
554 }
555 else if (locCoord[i] < -(1. + tol))
556 {
557 locCoord[i] = -(1. + tol);
558 clamp = true;
559 }
560 else if (locCoord[i] > (1. + tol))
561 {
562 locCoord[i] = 1. + tol;
563 clamp = true;
564 }
565 }
566 return clamp;
567}
568
569} // namespace Nektar::SpatialDomains
#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....
1D geometry information
Definition Geometry1D.h:49
2D geometry information
Definition Geometry2D.h:50
Base class for shape geometry information.
Definition Geometry.h:74
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition Geometry.cpp:188
virtual StdRegions::Orientation v_GetForient(const int i) const
Returns the orientation of face i with respect to the ordering of faces in the standard element.
Definition Geometry.cpp:177
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.
Definition Geometry.cpp:339
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition Geometry.h:203
virtual PointGeom * v_GetVertex(int i) const
Returns vertex i of this object.
Definition Geometry.cpp:126
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition Geometry.h:191
virtual void v_CalculateInverseIsoParam()
Definition Geometry.cpp:217
virtual int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord)
Definition Geometry.cpp:211
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition Geometry.h:548
int PreliminaryCheck(const Array< OneD, const NekDouble > &gloCoord)
A fast and robust check if a given global coord is outside of a deformed element. For regular element...
Definition Geometry.cpp:488
virtual int v_GetEdgeNormalToFaceVert(const int i, const int j) const
Returns the standard lement edge IDs that are normal to a given face vertex.
Definition Geometry.cpp:317
int GetShapeDim() const
Get the object's shape dimension.
Definition Geometry.h:430
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition Geometry.cpp:363
Geometry()
Default constructor.
Definition Geometry.cpp:50
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition Geometry.cpp:351
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition Geometry.cpp:196
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition Geometry.cpp:226
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition Geometry.cpp:372
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition Geometry.cpp:536
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
PointGeom * GetVertex(int i) const
Returns vertex i of this object.
Definition Geometry.h:361
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition Geometry.h:279
int GetFid(int i) const
Get the ID of face i of this object.
Definition Geometry.cpp:118
virtual int v_GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
Definition Geometry.cpp:306
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition Geometry.h:197
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition Geometry.h:201
GeomState m_geomFactorsState
State of the geometric factors.
Definition Geometry.h:187
virtual Geometry2D * v_GetFace(int i) const
Returns face i of this object.
Definition Geometry.cpp:146
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition Geometry.h:306
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition Geometry.h:189
virtual StdRegions::Orientation v_GetEorient(const int i) const
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
Definition Geometry.cpp:166
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition Geometry.h:439
static GeomFactorsVector m_regGeomFactorsManager
Definition Geometry.h:180
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition Geometry.cpp:399
int GetNumVerts() const
Get the number of vertices of this object.
Definition Geometry.h:403
Geometry1D * GetEdge(int i) const
Returns edge i of this object.
Definition Geometry.h:369
Geometry2D * GetFace(int i) const
Returns face i of this object.
Definition Geometry.h:377
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition Geometry.cpp:513
virtual int v_GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
Definition Geometry.cpp:284
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition Geometry.cpp:272
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition Geometry.h:297
static GeomFactorsSharedPtr ValidateRegGeomFactor(GeomFactorsSharedPtr geomFactor)
Check to see if a geometric factor has already been created that contains the same regular informatio...
Definition Geometry.cpp:81
virtual Geometry1D * v_GetEdge(int i) const
Returns edge i of this object.
Definition Geometry.cpp:136
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition Geometry.cpp:156
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition Geometry.h:185
int m_coordim
Coordinate dimension of this geometry object.
Definition Geometry.h:183
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition Geometry.cpp:328
virtual int v_GetVid(int i) const
Get the ID of vertex i of this object.
Definition Geometry.cpp:102
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition Geometry.cpp:204
virtual int v_GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
Definition Geometry.cpp:295
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
Determine whether an element contains a particular Cartesian coordinate .
Definition Geometry.cpp:237
int GetEid(int i) const
Get the ID of edge i of this object.
Definition Geometry.cpp:110
static const NekDouble kGeomFactorsTol
static const NekDouble kFindDistanceMin
bool GlobalIdEquality(const Geometry *&lhs, const Geometry *&rhs)
Definition Geometry.cpp:94
std::vector< GeomFactorsSharedPtr > GeomFactorsVector
A vector of GeomFactor pointers.
Definition GeomFactors.h:60
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition GeomFactors.h:58
bool SortByGlobalId(const Geometry *&lhs, const Geometry *&rhs)
Less than operator to sort Geometry objects by global id when sorting STL containers.
Definition Geometry.cpp:89
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition Curve.hpp:61
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eNotFilled
Geometric information has not been generated.
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300