Nektar++
Points.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Points.hpp
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: Header file of Points definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIB_UTILITIES_FOUNDATIONS_POINTS_H
36#define NEKTAR_LIB_UTILITIES_FOUNDATIONS_POINTS_H
37
44
46{
47
48/// Defines a specification for a set of points.
50{
51public:
52 // Used for looking up the creator. The creator for number of points
53 // can generate for any number, so we want the same creator called
54 // for all number.
55 struct opLess
56 {
58 const PointsKey &rhs) const;
59 };
60
61 /// Default constructor.
64 m_factor(NekConstants::kNekUnsetDouble)
65 {
66 }
67
68 /// Constructor defining the number and distribution of points.
69 PointsKey(const size_t &numpoints, const PointsType &pointstype,
71 : m_numpoints(numpoints), m_pointstype(pointstype), m_factor(factor)
72 {
73 }
74
75 /// Destructor.
76 virtual ~PointsKey()
77 {
78 }
79
80 /// Copy constructor.
81 PointsKey(const PointsKey &key) = default;
82
83 PointsKey &operator=(const PointsKey &key) = default;
84
85 inline size_t GetNumPoints() const
86 {
87 return m_numpoints;
88 }
89
91 {
92 return m_pointstype;
93 }
94
95 inline NekDouble GetFactor() const
96 {
97 return m_factor;
98 }
99
100 inline bool operator==(const PointsKey &key)
101 {
102
104 {
105 return (m_numpoints == key.m_numpoints &&
107 }
108
109 return false;
110 }
111
112 inline bool operator==(const PointsKey *y)
113 {
114 return (*this == *y);
115 }
116
117 inline bool operator!=(const PointsKey &y)
118 {
119 return (!(*this == y));
120 }
121
122 inline bool operator!=(const PointsKey *y)
123 {
124 return (!(*this == *y));
125 }
126
127 // If new points are added, this function must be modified
128 inline size_t GetPointsDim() const
129 {
130 size_t dimpoints = 1;
131
132 switch (m_pointstype)
133 {
134 case eNodalTriElec:
135 case eNodalTriFekete:
137 case eNodalTriSPI:
138 case eNodalQuadElec:
139 dimpoints = 2;
140 break;
141
142 case eNodalTetElec:
145 case eNodalPrismElec:
146 case eNodalHexElec:
147 dimpoints = 3;
148 break;
149
150 default:
151 break;
152 }
153
154 return dimpoints;
155 }
156
157 // If new points are added, this function must be modified
158 inline size_t GetTotNumPoints() const
159 {
160 size_t totpoints = m_numpoints;
161
162 switch (m_pointstype)
163 {
164 case eNodalTriElec:
165 case eNodalTriFekete:
167 totpoints = m_numpoints * (m_numpoints + 1) / 2;
168 break;
169 case eNodalTriSPI:
171 "This method cannot be implemented");
172 break;
173
174 case eNodalQuadElec:
175 totpoints = m_numpoints * m_numpoints;
176 break;
177
178 case eNodalTetElec:
180 totpoints =
181 m_numpoints * (m_numpoints + 1) * (m_numpoints + 2) / 6;
182 break;
183 case eNodalTetSPI:
185 "This method cannot be implemented");
186 break;
187
189 case eNodalPrismElec:
190 totpoints = m_numpoints * m_numpoints * (m_numpoints + 1) / 2;
191 break;
192 case eNodalPrismSPI:
194 "This method cannot be implemented");
195 break;
196
197 case eNodalHexElec:
198 totpoints = m_numpoints * m_numpoints * m_numpoints;
199 break;
200
201 default:
202 break;
203 }
204
205 return totpoints;
206 }
207
208 LIB_UTILITIES_EXPORT friend bool operator==(const PointsKey &lhs,
209 const PointsKey &rhs);
210 LIB_UTILITIES_EXPORT friend bool operator<(const PointsKey &lhs,
211 const PointsKey &rhs);
213 const PointsKey &lhs, const PointsKey &rhs) const;
214
215protected:
216 size_t m_numpoints; //!< number of the points (as appropriately
217 //!< defined for PointsType)
218 PointsType m_pointstype; //!< Type of Points
219 NekDouble m_factor; //!< optional factor
220private:
221};
222
224
226 const PointsKey &rhs);
227LIB_UTILITIES_EXPORT bool operator<(const PointsKey &lhs, const PointsKey &rhs);
228LIB_UTILITIES_EXPORT std::ostream &operator<<(std::ostream &os,
229 const PointsKey &rhs);
230
231typedef std::vector<PointsKey> PointsKeyVector;
232
233/// Stores a set of points of datatype DataT, defined by a PointKey.
234template <typename DataT> class Points
235{
236public:
237 typedef DataT DataType;
238 typedef std::shared_ptr<NekMatrix<DataType>> MatrixSharedPtrType;
239
240 virtual ~Points()
241 {
242 }
243
244 inline void Initialize(void)
245 {
246 v_Initialize();
247 }
248
249 inline size_t GetPointsDim() const
250 {
251 return m_pointsKey.GetPointsDim();
252 }
253
254 inline size_t GetNumPoints() const
255 {
256 return m_pointsKey.GetNumPoints();
257 }
258
259 inline size_t GetTotNumPoints() const
260 {
262 }
263
265 {
266 return m_pointsKey.GetPointsType();
267 }
268
269 inline const Array<OneD, const DataType> &GetZ() const
270 {
271 return m_points[0];
272 }
273
274 inline const Array<OneD, const DataType> &GetW() const
275 {
276 return m_weights;
277 }
278
281 {
282 z = m_points[0];
283 w = m_weights;
284 }
285
287 {
288 return m_bcweights;
289 }
290
292 {
293 x = m_points[0];
294 }
295
298 {
299 x = m_points[0];
300 y = m_points[1];
301 }
302
306 {
307 x = m_points[0];
308 y = m_points[1];
309 z = m_points[2];
310 }
311
312 inline const MatrixSharedPtrType &GetD(Direction dir = xDir) const
313 {
314 return m_derivmatrix[(int)dir];
315 }
316
318 {
319 return v_GetI(key);
320 }
321
323 {
324 return v_GetI(x);
325 }
326
327 const MatrixSharedPtrType GetI(size_t uint,
329 {
330 return v_GetI(uint, x);
331 }
332
335 {
336 return v_GetI(x, y);
337 }
338
342 {
343 return v_GetI(x, y, z);
344 }
345
347 {
348 return v_GetGalerkinProjection(pkey);
349 }
350
351protected:
352 /// Points type for this points distributions.
354 /// Storage for the point locations, allowing for up to a 3D points
355 /// storage.
357 /// Quadrature weights for the weights.
359 /// Barycentric weights.
361 /// Derivative matrices.
367
368 virtual void v_Initialize(void)
369 {
372 v_CalculateBaryWeights();
373 v_CalculateDerivMatrix();
374 }
375
376 virtual void v_CalculatePoints()
377 {
378 size_t pointsDim = GetPointsDim();
379 size_t totNumPoints = GetTotNumPoints();
380
381 for (size_t i = 0; i < pointsDim; ++i)
382 {
383 m_points[i] = Array<OneD, DataType>(totNumPoints);
384 }
385 }
386
387 virtual void v_CalculateWeights()
388 {
390 }
391
392 /**
393 * @brief This function calculates the barycentric weights used for
394 * enhanced interpolation speed.
395 *
396 * For the points distribution \f$ z_i \f$ with \f% 1\leq z_i \leq N
397 * \f$, the barycentric weights are computed as:
398 *
399 * \f[
400 * b_i=\prod_{\substack{1\leq j\leq N\\ i\neq j}} \frac{1}{z_i-z_j}
401 * \f]
402 */
403 virtual void v_CalculateBaryWeights()
404 {
405 const size_t totNumPoints = m_pointsKey.GetNumPoints();
406 m_bcweights = Array<OneD, DataType>(totNumPoints, 1.0);
407
409
410 for (size_t i = 0; i < totNumPoints; ++i)
411 {
412 for (size_t j = 0; j < totNumPoints; ++j)
413 {
414 if (i == j)
415 {
416 continue;
417 }
418
419 m_bcweights[i] *= (z[i] - z[j]);
420 }
421
422 m_bcweights[i] = 1.0 / m_bcweights[i];
423 }
424 }
425
426 virtual void v_CalculateDerivMatrix()
427 {
428 size_t totNumPoints = GetTotNumPoints();
429 for (size_t i = 0; i < m_pointsKey.GetPointsDim(); ++i)
430 {
431 m_derivmatrix[i] =
432 MemoryManager<NekMatrix<DataType>>::AllocateSharedPtr(
433 totNumPoints, totNumPoints);
434 }
435 }
436
437 Points(const PointsKey &key) : m_pointsKey(key)
438 {
439 }
440
441 virtual const MatrixSharedPtrType v_GetI(
442 [[maybe_unused]] const PointsKey &key)
443 {
444 NEKERROR(ErrorUtil::efatal, "Method not implemented ");
445 std::shared_ptr<NekMatrix<NekDouble>> returnval(
446 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
447 return returnval;
448 }
449
450 virtual const MatrixSharedPtrType v_GetI(
451 [[maybe_unused]] const Array<OneD, const DataType> &x)
452 {
453 NEKERROR(ErrorUtil::efatal, "Method not implemented");
454 std::shared_ptr<NekMatrix<NekDouble>> returnval(
455 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
456 return returnval;
457 }
458
459 virtual const MatrixSharedPtrType v_GetI(
460 size_t, [[maybe_unused]] const Array<OneD, const DataType> &x)
461 {
462 NEKERROR(ErrorUtil::efatal, "Method not implemented");
463 std::shared_ptr<NekMatrix<NekDouble>> returnval(
464 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
465 return returnval;
466 }
467
468 virtual const MatrixSharedPtrType v_GetI(
469 [[maybe_unused]] const Array<OneD, const DataType> &x,
470 [[maybe_unused]] const Array<OneD, const DataType> &y)
471 {
472 NEKERROR(ErrorUtil::efatal, "Method not implemented");
473 std::shared_ptr<NekMatrix<NekDouble>> returnval(
474 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
475 return returnval;
476 }
477
478 virtual const MatrixSharedPtrType v_GetI(
479 [[maybe_unused]] const Array<OneD, const DataType> &x,
480 [[maybe_unused]] const Array<OneD, const DataType> &y,
481 [[maybe_unused]] const Array<OneD, const DataType> &z)
482 {
483 NEKERROR(ErrorUtil::efatal, "Method not implemented");
484 std::shared_ptr<NekMatrix<NekDouble>> returnval(
485 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
486 return returnval;
487 }
488
489 virtual const MatrixSharedPtrType v_GetGalerkinProjection(
490 [[maybe_unused]] const PointsKey &pkey)
491 {
492 NEKERROR(ErrorUtil::efatal, "Method not implemented ");
493 std::shared_ptr<NekMatrix<NekDouble>> returnval(
494 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr());
495 return returnval;
496 }
497
498private:
499 Points(const Points &pts) = delete;
500 Points() = delete;
501};
502
503} // namespace Nektar::LibUtilities
504
505#endif // NEKTAR_LIB_UTILITIES_FOUNDATIONS_POINTS_H
#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 LIB_UTILITIES_EXPORT
std::shared_ptr< NekMatrix< NekDouble > > MatrixSharedPtrType
1D Array of constant elements with garbage collection and bounds checking.
Stores a set of points of datatype DataT, defined by a PointKey.
Definition: Points.h:235
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:356
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:364
size_t GetPointsDim() const
Definition: Points.h:249
void GetPoints(Array< OneD, const DataType > &x, Array< OneD, const DataType > &y) const
Definition: Points.h:296
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:366
virtual void v_CalculatePoints()
Definition: Points.h:376
const MatrixSharedPtrType GetI(size_t uint, const Array< OneD, const DataType > &x)
Definition: Points.h:327
size_t GetNumPoints() const
Definition: Points.h:254
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:353
void GetZW(Array< OneD, const DataType > &z, Array< OneD, const DataType > &w) const
Definition: Points.h:279
virtual void v_Initialize(void)
Definition: Points.h:368
void GetPoints(Array< OneD, const DataType > &x) const
Definition: Points.h:291
void GetPoints(Array< OneD, const DataType > &x, Array< OneD, const DataType > &y, Array< OneD, const DataType > &z) const
Definition: Points.h:303
PointsType GetPointsType() const
Definition: Points.h:264
const MatrixSharedPtrType GetI(const Array< OneD, const DataType > &x)
Definition: Points.h:322
const Array< OneD, const DataType > & GetW() const
Definition: Points.h:274
const Array< OneD, const NekDouble > & GetBaryWeights() const
Definition: Points.h:286
const MatrixSharedPtrType GetI(const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y, const Array< OneD, const DataType > &z)
Definition: Points.h:339
size_t GetTotNumPoints() const
Definition: Points.h:259
std::shared_ptr< NekMatrix< DataType > > MatrixSharedPtrType
Definition: Points.h:238
const MatrixSharedPtrType & GetD(Direction dir=xDir) const
Definition: Points.h:312
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:358
const MatrixSharedPtrType GetI(const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
Definition: Points.h:333
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:317
const Array< OneD, const DataType > & GetZ() const
Definition: Points.h:269
virtual void v_CalculateWeights()
Definition: Points.h:387
const MatrixSharedPtrType GetGalerkinProjection(const PointsKey &pkey)
Definition: Points.h:346
Array< OneD, DataType > m_bcweights
Barycentric weights.
Definition: Points.h:360
Defines a specification for a set of points.
Definition: Points.h:50
size_t m_numpoints
number of the points (as appropriately defined for PointsType)
Definition: Points.h:216
bool operator==(const PointsKey *y)
Definition: Points.h:112
size_t GetPointsDim() const
Definition: Points.h:128
PointsKey(const PointsKey &key)=default
Copy constructor.
PointsType GetPointsType() const
Definition: Points.h:90
PointsKey & operator=(const PointsKey &key)=default
friend bool operator<(const PointsKey &lhs, const PointsKey &rhs)
PointsType m_pointstype
Type of Points.
Definition: Points.h:218
PointsKey(const size_t &numpoints, const PointsType &pointstype, const NekDouble factor=NekConstants::kNekUnsetDouble)
Constructor defining the number and distribution of points.
Definition: Points.h:69
size_t GetTotNumPoints() const
Definition: Points.h:158
bool operator==(const PointsKey &key)
Definition: Points.h:100
size_t GetNumPoints() const
Definition: Points.h:85
NekDouble GetFactor() const
Definition: Points.h:95
NekDouble m_factor
optional factor
Definition: Points.h:219
virtual ~PointsKey()
Destructor.
Definition: Points.h:76
bool operator!=(const PointsKey &y)
Definition: Points.h:117
PointsKey(void)
Default constructor.
Definition: Points.h:62
bool operator!=(const PointsKey *y)
Definition: Points.h:122
bool operator==(const BasisKey &x, const BasisKey &y)
bool operator<(const BasisKey &lhs, const BasisKey &rhs)
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
@ eNodalPrismEvenlySpaced
3D Evenly-spaced points on a Prism
Definition: PointsType.h:86
@ eNodalTriFekete
2D Nodal Fekete Points on a Triangle
Definition: PointsType.h:82
@ eNodalPrismSPI
3D prism SPI
Definition: PointsType.h:92
@ eNodalTriElec
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:81
@ eNodalTriEvenlySpaced
2D Evenly-spaced points on a Triangle
Definition: PointsType.h:83
@ eNodalHexElec
3D GLL for hex
Definition: PointsType.h:94
@ eNodalQuadElec
2D GLL for quad
Definition: PointsType.h:93
@ eNodalTetEvenlySpaced
3D Evenly-spaced points on a Tetrahedron
Definition: PointsType.h:84
@ eNodalTetSPI
3D Nodal Symmetric positive internal tet (Whitherden, Vincent)
Definition: PointsType.h:90
@ eNodalPrismElec
3D electrostatically spaced points on a Prism
Definition: PointsType.h:87
@ eNodalTetElec
3D Nodal Electrostatic Points on a Tetrahedron
Definition: PointsType.h:85
@ eNodalTriSPI
2D Nodal Symmetric positive internal triangle (Whitherden, Vincent)
Definition: PointsType.h:88
static const PointsKey NullPointsKey(0, eNoPointsType)
static const NekDouble kNekUnsetDouble
static const NekDouble kNekZeroTol
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
double NekDouble
bool operator()(const PointsKey &lhs, const PointsKey &rhs) const