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