Nektar++
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 Default destructor.
69 */
71{
72}
73
74/**
75 * @brief Check to see if a geometric factor has already been created that
76 * contains the same regular information.
77 *
78 * The principle behind this is that many regular (i.e. constant Jacobian)
79 * elements have identicial geometric factors. Memory may therefore be reduced
80 * by storing only the unique factors.
81 *
82 * @param geomFactor The GeomFactor to check.
83 *
84 * @return Either the cached GeomFactor or @p geomFactor.
85 *
86 * @todo Currently this method is disabled since the lookup is very expensive.
87 */
89 GeomFactorsSharedPtr geomFactor)
90{
91 GeomFactorsSharedPtr returnval = geomFactor;
92
93/// \todo should this '#if 0' statement be removed?
94#if 0
95 bool found = false;
96 if (geomFactor->GetGtype() == eRegular)
97 {
98 for (GeomFactorsVectorIter iter = m_regGeomFactorsManager.begin();
99 iter != m_regGeomFactorsManager.end();
100 ++iter)
101 {
102 if (**iter == *geomFactor)
103 {
104 returnval = *iter;
105 found = true;
106 break;
107 }
108 }
109
110 if (!found)
111 {
112 m_regGeomFactorsManager.push_back(geomFactor);
113 returnval = geomFactor;
114 }
115 }
116#endif
117 return returnval;
118}
119
120bool SortByGlobalId(const std::shared_ptr<Geometry> &lhs,
121 const std::shared_ptr<Geometry> &rhs)
122{
123 return lhs->GetGlobalID() < rhs->GetGlobalID();
124}
125
126bool GlobalIdEquality(const std::shared_ptr<Geometry> &lhs,
127 const std::shared_ptr<Geometry> &rhs)
128{
129 return lhs->GetGlobalID() == rhs->GetGlobalID();
130}
131
132/**
133 * @brief Get the ID of vertex @p i of this object.
134 */
135int Geometry::GetVid(int i) const
136{
137 return GetVertex(i)->GetGlobalID();
138}
139
140/**
141 * @brief Get the ID of edge @p i of this object.
142 */
143int Geometry::GetEid(int i) const
144{
145 return GetEdge(i)->GetGlobalID();
146}
147
148/**
149 * @brief Get the ID of face @p i of this object.
150 */
151int Geometry::GetFid(int i) const
152{
153 return GetFace(i)->GetGlobalID();
154}
155
156/**
157 * @copydoc Geometry::GetEdge()
158 */
159Geometry1DSharedPtr Geometry::v_GetEdge([[maybe_unused]] int i) const
160{
162 "This function is only valid for shape type geometries");
163 return Geometry1DSharedPtr();
164}
165
166/**
167 * @copydoc Geometry::GetFace()
168 */
169Geometry2DSharedPtr Geometry::v_GetFace([[maybe_unused]] int i) const
170{
172 "This function is only valid for shape type geometries");
173 return Geometry2DSharedPtr();
174}
175
176/**
177 * @copydoc Geometry::GetNumVerts()
178 */
180{
182 "This function is only valid for shape type geometries");
183 return 0;
184}
185
186/**
187 * @copydoc Geometry::GetEorient()
188 */
190 [[maybe_unused]] const int i) const
191{
193 "This function is not valid for this geometry.");
195}
196
197/**
198 * @copydoc Geometry::GetForient()
199 */
201 [[maybe_unused]] const int i) const
202{
204 "This function is not valid for this geometry.");
205 return StdRegions::eFwd;
206}
207
208/**
209 * @copydoc Geometry::GetNumEdges()
210 */
212{
213 return 0;
214}
215
216/**
217 * @copydoc Geometry::GetNumFaces()
218 */
220{
221 return 0;
222}
223
224/**
225 * @copydoc Geometry::GetShapeDim()
226 */
228{
230 "This function is only valid for shape type geometries");
231 return 0;
232}
233
235 [[maybe_unused]] const Array<OneD, const NekDouble> &gloCoord)
236{
237 return 0;
238}
239
241{
243 "This function is only valid for shape type geometries");
244}
245
246/**
247 * @copydoc Geometry::GetXmap()
248 */
250{
251 return m_xmap;
252}
253
254/**
255 * @copydoc Geometry::ContainsPoint(
256 * const Array<OneD, const NekDouble> &, Array<OneD, NekDouble> &,
257 * NekDouble, NekDouble&)
258 * dist is assigned value for curved elements
259 */
261 Array<OneD, NekDouble> &locCoord, NekDouble tol,
262 NekDouble &dist)
263{
264 int inside = PreliminaryCheck(gloCoord);
265 if (inside == -1)
266 {
267 dist = std::numeric_limits<double>::max();
268 return false;
269 }
270 dist = GetLocCoords(gloCoord, locCoord);
271 if (inside == 1)
272 {
273 dist = 0.;
274 return true;
275 }
276 else
277 {
279 m_xmap->LocCoordToLocCollapsed(locCoord, eta);
280 return !ClampLocCoords(eta, tol);
281 }
282}
283
285 [[maybe_unused]] const Array<OneD, const NekDouble> &xs,
286 [[maybe_unused]] Array<OneD, NekDouble> &xi)
287{
289 "This function has not been defined for this geometry");
290 return false;
291}
292
293/**
294 * @copydoc Geometry::GetVertexEdgeMap()
295 */
296int Geometry::v_GetVertexEdgeMap([[maybe_unused]] const int i,
297 [[maybe_unused]] const int j) const
298{
300 "This function has not been defined for this geometry");
301 return 0;
302}
303
304/**
305 * @copydoc Geometry::GetVertexFaceMap()
306 */
307int Geometry::v_GetVertexFaceMap([[maybe_unused]] const int i,
308 [[maybe_unused]] const int j) const
309{
311 "This function has not been defined for this geometry");
312 return 0;
313}
314
315/**
316 * @copydoc Geometry::GetEdgeFaceMap()
317 */
318int Geometry::v_GetEdgeFaceMap([[maybe_unused]] const int i,
319 [[maybe_unused]] const int j) const
320{
322 "This function has not been defined for this geometry");
323 return 0;
324}
325
326/**
327 * @copydoc Geometry::GetEdgeNormalToFaceVert()
328 */
329int Geometry::v_GetEdgeNormalToFaceVert([[maybe_unused]] const int i,
330 [[maybe_unused]] const int j) const
331{
333 "This function has not been defined for this geometry");
334 return 0;
335}
336
337/**
338 * @copydoc Geometry::GetDir()
339 */
340int Geometry::v_GetDir([[maybe_unused]] const int i,
341 [[maybe_unused]] const int j) const
342{
344 "This function has not been defined for this geometry");
345 return 0;
346}
347
348/**
349 * @copydoc Geometry::GetCoord()
350 */
352 [[maybe_unused]] const int i,
353 [[maybe_unused]] const Array<OneD, const NekDouble> &Lcoord)
354{
356 "This function is only valid for expansion type geometries");
357 return 0.0;
358}
359
360/**
361 * @copydoc Geometry::GetLocCoords()
362 */
364 [[maybe_unused]] const Array<OneD, const NekDouble> &coords,
365 [[maybe_unused]] Array<OneD, NekDouble> &Lcoords)
366{
368 "This function is only valid for expansion type geometries");
369 return 0.0;
370}
371
372/**
373 * @copydoc Geometry::FillGeom()
374 */
376{
378 "This function is only valid for expansion type geometries");
379}
380
381/**
382 * @copydoc Geometry::Reset()
383 */
384void Geometry::v_Reset([[maybe_unused]] CurveMap &curvedEdges,
385 [[maybe_unused]] CurveMap &curvedFaces)
386{
387
388 // Reset state
391
392 // Junk geometric factors
394}
395
397{
399 "This function is only valid for expansion type geometries");
400}
401
402/**
403 * @brief Generates the bounding box for the element.
404 *
405 * For regular elements, the vertices are sufficient to define the extent of
406 * the bounding box. For non-regular elements, the extremes of the quadrature
407 * point coordinates are used. A 10% margin is added around this computed
408 * region to account for convex hull elements where the true extent of the
409 * element may extend slightly beyond the quadrature points.
410 */
411std::array<NekDouble, 6> Geometry::GetBoundingBox()
412{
413 if (m_boundingBox.size() == 6)
414 {
415 return {{m_boundingBox[0], m_boundingBox[1], m_boundingBox[2],
417 }
418 // NekDouble minx, miny, minz, maxx, maxy, maxz;
419 Array<OneD, NekDouble> min(3), max(3);
420
421 // Always get vertexes min/max
423 Array<OneD, NekDouble> x(3, 0.0);
424 p->GetCoords(x[0], x[1], x[2]);
425 for (int j = 0; j < 3; ++j)
426 {
427 min[j] = x[j];
428 max[j] = x[j];
429 }
430 for (int i = 1; i < GetNumVerts(); ++i)
431 {
432 p = GetVertex(i);
433 p->GetCoords(x[0], x[1], x[2]);
434 for (int j = 0; j < 3; ++j)
435 {
436 min[j] = (x[j] < min[j] ? x[j] : min[j]);
437 max[j] = (x[j] > max[j] ? x[j] : max[j]);
438 }
439 }
440 // If element is deformed loop over quadrature points
442 if (GetGeomFactors()->GetGtype() != eRegular)
443 {
444 marginFactor = 0.1;
445 const int nq = GetXmap()->GetTotPoints();
447 for (int j = 0; j < 3; ++j)
448 {
449 x[j] = Array<OneD, NekDouble>(nq, 0.0);
450 }
451 for (int j = 0; j < GetCoordim(); ++j)
452 {
453 GetXmap()->BwdTrans(m_coeffs[j], x[j]);
454 }
455 for (int j = 0; j < 3; ++j)
456 {
457 for (int i = 0; i < nq; ++i)
458 {
459 min[j] = (x[j][i] < min[j] ? x[j][i] : min[j]);
460 max[j] = (x[j][i] > max[j] ? x[j][i] : max[j]);
461 }
462 }
463 }
464 // Add margin to bounding box, in order to
465 // return the nearest element
466 for (int j = 0; j < 3; ++j)
467 {
468 NekDouble margin =
469 marginFactor * (max[j] - min[j]) + NekConstants::kFindDistanceMin;
470 min[j] -= margin;
471 max[j] += margin;
472 }
473
474 // save bounding box
476 for (int j = 0; j < 3; ++j)
477 {
478 m_boundingBox[j] = min[j];
479 m_boundingBox[j + 3] = max[j];
480 }
481 // Return bounding box
482 return {{min[0], min[1], min[2], max[0], max[1], max[2]}};
483}
484
486{
487 m_boundingBox = {};
488}
489
490/**
491 * @brief A fast and robust check if a given global coord is
492 * outside of a deformed element. For regular elements, this
493 * check is unnecessary
494 *
495 * @param coords Input Cartesian global coordinates
496 *
497 * @return 1 is inside of the element.
498 * 0 maybe inside
499 * -1 outside of the element
500 */
502{
503 // bounding box check
504 if (!MinMaxCheck(gloCoord))
505 {
506 return -1;
507 }
508
509 // regular element check
510 if (GetMetricInfo()->GetGtype() == eRegular)
511 {
512 return 0;
513 }
514
515 // All left check for straight edges/plane surfaces
516 return v_AllLeftCheck(gloCoord);
517}
518
519/**
520 * @brief Check if given global coord is within the BoundingBox
521 * of the element.
522 *
523 * @param coords Input Cartesian global coordinates
524 *
525 * @return True if within distance or False otherwise.
526 */
528{
529 // Validation checks
530 ASSERTL1(gloCoord.size() >= m_coordim,
531 "Expects number of global coordinates supplied to be greater than "
532 "or equal to the mesh dimension.");
533
534 std::array<NekDouble, 6> minMax = GetBoundingBox();
535 for (int i = 0; i < m_coordim; ++i)
536 {
537 if ((gloCoord[i] < minMax[i]) || (gloCoord[i] > minMax[i + 3]))
538 {
539 return false;
540 }
541 }
542 return true;
543}
544
545/**
546 * @brief Clamp local coords to be within standard regions [-1, 1]^dim.
547 *
548 * @param Lcoords Corresponding local coordinates
549 */
551{
552 // Validation checks
553 ASSERTL1(locCoord.size() >= GetShapeDim(),
554 "Expects local coordinates to be same or "
555 "larger than shape dimension.");
556
557 // If out of range clamp locCoord to be within [-1,1]^dim
558 // since any larger value will be very oscillatory if
559 // called by 'returnNearestElmt' option in
560 // ExpList::GetExpIndex
561 bool clamp = false;
562 for (int i = 0; i < GetShapeDim(); ++i)
563 {
564 if (!std::isfinite(locCoord[i]))
565 {
566 locCoord[i] = 0.;
567 clamp = true;
568 }
569 else if (locCoord[i] < -(1. + tol))
570 {
571 locCoord[i] = -(1. + tol);
572 clamp = true;
573 }
574 else if (locCoord[i] > (1. + tol))
575 {
576 locCoord[i] = 1. + tol;
577 clamp = true;
578 }
579 }
580 return clamp;
581}
582
583} // namespace Nektar::SpatialDomains
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.cpp:211
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:200
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:351
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition: Geometry.h:203
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:355
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:191
virtual void v_CalculateInverseIsoParam()
Definition: Geometry.cpp:240
virtual int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord)
Definition: Geometry.cpp:234
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:542
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:501
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:329
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:363
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:424
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:375
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:363
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
Definition: Geometry.cpp:169
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:135
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.cpp:219
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.cpp:249
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:384
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:550
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:281
int GetFid(int i) const
Get the ID of face i of this object.
Definition: Geometry.cpp:151
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:318
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
virtual ~Geometry()
Default destructor.
Definition: Geometry.cpp:70
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:187
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:308
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.cpp:159
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:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:433
static GeomFactorsVector m_regGeomFactorsManager
Definition: Geometry.h:180
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:411
int GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.h:397
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:371
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition: Geometry.cpp:527
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:296
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition: Geometry.cpp:284
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition: Geometry.h:299
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:88
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.cpp:179
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:340
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.cpp:227
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:307
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:260
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:143
static const NekDouble kGeomFactorsTol
static const NekDouble kFindDistanceMin
std::vector< GeomFactorsSharedPtr > GeomFactorsVector
A vector of GeomFactor pointers.
Definition: GeomFactors.h:62
bool GlobalIdEquality(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Definition: Geometry.cpp:126
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
bool SortByGlobalId(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Less than operator to sort Geometry objects by global id when sorting STL containers.
Definition: Geometry.cpp:120
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:59
@ eRegular
Geometry is straight-sided with constant geometric factors.
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:62
@ eNotFilled
Geometric information has not been generated.
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:61
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
double NekDouble