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