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