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 
43 namespace Nektar
44 {
45 namespace 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)
58 {
59 }
60 
61 /**
62  * @brief Constructor when supplied a coordinate dimension.
63  */
64 Geometry::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)
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 
124 bool SortByGlobalId(const std::shared_ptr<Geometry> &lhs,
125  const std::shared_ptr<Geometry> &rhs)
126 {
127  return lhs->GetGlobalID() < rhs->GetGlobalID();
128 }
129 
130 bool 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  */
139 int 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  */
147 int 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  */
155 int 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.");
200  return StdRegions::eForwards;
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 
240 /**
241  * @copydoc Geometry::GetXmap()
242  */
244 {
245  return m_xmap;
246 }
247 
248 /**
249  * @copydoc Geometry::ContainsPoint(
250  * const Array<OneD, const NekDouble> &, Array<OneD, NekDouble> &,
251  * NekDouble, NekDouble&)
252  */
254  Array<OneD, NekDouble> &locCoord, NekDouble tol,
255  NekDouble &dist)
256 {
257  // Convert to the local (xi) coordinates.
258  dist = GetLocCoords(gloCoord, locCoord);
259  if (dist <= tol + NekConstants::kNekMachineEpsilon)
260  {
261  return true;
262  }
264  m_xmap->LocCoordToLocCollapsed(locCoord, eta);
265  if (ClampLocCoords(eta, tol))
266  {
267  m_xmap->LocCollapsedToLocCoord(eta, locCoord);
268  return false;
269  }
270  else
271  {
272  return true;
273  }
274 }
275 
278 {
279  boost::ignore_unused(xs, xi);
281  "This function has not been defined for this geometry");
282  return false;
283 }
284 
285 /**
286  * @copydoc Geometry::GetVertexEdgeMap()
287  */
288 int Geometry::v_GetVertexEdgeMap(const int i, const int j) const
289 {
290  boost::ignore_unused(i, j);
292  "This function has not been defined for this geometry");
293  return 0;
294 }
295 
296 /**
297  * @copydoc Geometry::GetVertexFaceMap()
298  */
299 int Geometry::v_GetVertexFaceMap(const int i, const int j) const
300 {
301  boost::ignore_unused(i, j);
303  "This function has not been defined for this geometry");
304  return 0;
305 }
306 
307 /**
308  * @copydoc Geometry::GetEdgeFaceMap()
309  */
310 int Geometry::v_GetEdgeFaceMap(const int i, const int j) const
311 {
312  boost::ignore_unused(i, j);
314  "This function has not been defined for this geometry");
315  return 0;
316 }
317 
318 /**
319  * @copydoc Geometry::GetEdgeNormalToFaceVert()
320  */
321 int Geometry::v_GetEdgeNormalToFaceVert(const int i, const int j) const
322 {
323  boost::ignore_unused(i, j);
325  "This function has not been defined for this geometry");
326  return 0;
327 }
328 
329 /**
330  * @copydoc Geometry::GetDir()
331  */
332 int Geometry::v_GetDir(const int i, const int j) const
333 {
334  boost::ignore_unused(i, j);
336  "This function has not been defined for this geometry");
337  return 0;
338 }
339 
340 /**
341  * @copydoc Geometry::GetCoord()
342  */
344  const Array<OneD, const NekDouble> &Lcoord)
345 {
346  boost::ignore_unused(i, Lcoord);
348  "This function is only valid for expansion type geometries");
349  return 0.0;
350 }
351 
352 /**
353  * @copydoc Geometry::GetLocCoords()
354  */
356  Array<OneD, NekDouble> &Lcoords)
357 {
358  boost::ignore_unused(coords, Lcoords);
360  "This function is only valid for expansion type geometries");
361  return 0.0;
362 }
363 
364 /**
365  * @copydoc Geometry::FillGeom()
366  */
368 {
370  "This function is only valid for expansion type geometries");
371 }
372 
373 /**
374  * @copydoc Geometry::Reset()
375  */
376 void Geometry::v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
377 {
378  boost::ignore_unused(curvedEdges, curvedFaces);
379 
380  // Reset state
383 
384  // Junk geometric factors
386 }
387 
389 {
391  "This function is only valid for expansion type geometries");
392 }
393 
394 /**
395  * @brief Generates the bounding box for the element.
396  *
397  * For regular elements, the vertices are sufficient to define the extent of
398  * the bounding box. For non-regular elements, the extremes of the quadrature
399  * point coordinates are used. A 10% margin is added around this computed
400  * region to account for convex hull elements where the true extent of the
401  * element may extend slightly beyond the quadrature points.
402  */
403 std::array<NekDouble, 6> Geometry::GetBoundingBox()
404 {
405  if (m_boundingBox.size() == 6)
406  {
407  return {{m_boundingBox[0], m_boundingBox[1], m_boundingBox[2],
409  }
410  // NekDouble minx, miny, minz, maxx, maxy, maxz;
411  Array<OneD, NekDouble> min(3), max(3);
412 
413  // Always get vertexes min/max
415  Array<OneD, NekDouble> x(3, 0.0);
416  p->GetCoords(x[0], x[1], x[2]);
417  for (int j = 0; j < 3; ++j)
418  {
419  min[j] = x[j];
420  max[j] = x[j];
421  }
422  for (int i = 1; i < GetNumVerts(); ++i)
423  {
424  p = GetVertex(i);
425  p->GetCoords(x[0], x[1], x[2]);
426  for (int j = 0; j < 3; ++j)
427  {
428  min[j] = (x[j] < min[j] ? x[j] : min[j]);
429  max[j] = (x[j] > max[j] ? x[j] : max[j]);
430  }
431  }
432  // If element is deformed loop over quadrature points
433  if (GetGeomFactors()->GetGtype() != eRegular)
434  {
435  const int nq = GetXmap()->GetTotPoints();
437  for (int j = 0; j < 3; ++j)
438  {
439  x[j] = Array<OneD, NekDouble>(nq, 0.0);
440  }
441  for (int j = 0; j < GetCoordim(); ++j)
442  {
443  GetXmap()->BwdTrans(m_coeffs[j], x[j]);
444  }
445  for (int j = 0; j < 3; ++j)
446  {
447  for (int i = 0; i < nq; ++i)
448  {
449  min[j] = (x[j][i] < min[j] ? x[j][i] : min[j]);
450  max[j] = (x[j][i] > max[j] ? x[j][i] : max[j]);
451  }
452  }
453  }
454  // Add 10% margin to bounding box, in order to
455  // return the nearest element
456  for (int j = 0; j < 3; ++j)
457  {
458  const NekDouble len = max[j] - min[j];
459  min[j] -= (0.1 * len + NekConstants::kGeomFactorsTol);
460  max[j] += (0.1 * len + NekConstants::kGeomFactorsTol);
461  }
462 
463  // save bounding box
465  for (int j = 0; j < 3; ++j)
466  {
467  m_boundingBox[j] = min[j];
468  m_boundingBox[j + 3] = max[j];
469  }
470  // Return bounding box
471  return {{min[0], min[1], min[2], max[0], max[1], max[2]}};
472 }
473 
475 {
476  m_boundingBox = {};
477 }
478 
479 /**
480  * @brief Check if given global coord is within the BoundingBox
481  * of the element.
482  *
483  * @param coords Input Cartesian global coordinates
484  *
485  * @return True if within distance or False otherwise.
486  */
488 {
489  // Validation checks
490  ASSERTL1(gloCoord.size() >= m_coordim,
491  "Expects number of global coordinates supplied to be greater than "
492  "or equal to the mesh dimension.");
493 
494  std::array<NekDouble, 6> minMax = GetBoundingBox();
495  for (int i = 0; i < m_coordim; ++i)
496  {
497  if ((gloCoord[i] < minMax[i]) || (gloCoord[i] > minMax[i + 3]))
498  {
499  return false;
500  }
501  }
502  return true;
503 }
504 
505 /**
506  * @brief Clamp local coords to be within standard regions [-1, 1]^dim.
507  *
508  * @param Lcoords Corresponding local coordinates
509  */
511 {
512  // Validation checks
513  ASSERTL1(locCoord.size() >= GetShapeDim(),
514  "Expects local coordinates to be same or "
515  "larger than shape dimension.");
516 
517  // If out of range clamp locCoord to be within [-1,1]^dim
518  // since any larger value will be very oscillatory if
519  // called by 'returnNearestElmt' option in
520  // ExpList::GetExpIndex
521  bool clamp = false;
522  for (int i = 0; i < GetShapeDim(); ++i)
523  {
524  if (!std::isfinite(locCoord[i]))
525  {
526  locCoord[i] = 0.;
527  clamp = true;
528  }
529  else if (locCoord[i] < -(1. + tol))
530  {
531  locCoord[i] = -(1. + tol);
532  clamp = true;
533  }
534  else if (locCoord[i] > (1. + tol))
535  {
536  locCoord[i] = 1. + tol;
537  clamp = true;
538  }
539  }
540  return clamp;
541 }
542 
543 } // namespace SpatialDomains
544 } // 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:343
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition: Geometry.h:204
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:351
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:192
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:538
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:321
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:359
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:420
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:367
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:355
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:243
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:376
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:510
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:277
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:310
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:202
virtual ~Geometry()
Default destructor.
Definition: Geometry.cpp:74
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:188
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:190
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:429
static GeomFactorsVector m_regGeomFactorsManager
Definition: Geometry.h:181
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:403
int GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.h:393
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:367
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition: Geometry.cpp:487
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:288
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition: Geometry.cpp:276
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition: Geometry.h:295
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:186
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:184
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:332
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:299
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
Definition: Geometry.cpp:253
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 kNekMachineEpsilon
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:59
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
@ eNotFilled
Geometric information has not been generated.
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:63
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble