Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Static Private Attributes | List of all members
Nektar::LibUtilities::GaussPoints Class Reference

#include <GaussPoints.h>

Inheritance diagram for Nektar::LibUtilities::GaussPoints:
[legend]

Public Member Functions

 ~GaussPoints () override
 
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix (const PointsKey &pkey)
 
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix (const PointsKey &pkey)
 
 GaussPoints (const PointsKey &pkey)
 
- Public Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual ~Points ()
 
void Initialize (void)
 
size_t GetPointsDim () const
 
size_t GetNumPoints () const
 
size_t GetTotNumPoints () const
 
PointsType GetPointsType () const
 
const Array< OneD, const DataType > & GetZ () const
 
const Array< OneD, const DataType > & GetW () const
 
void GetZW (Array< OneD, const DataType > &z, Array< OneD, const DataType > &w) const
 
const Array< OneD, const NekDouble > & GetBaryWeights () const
 
void GetPoints (Array< OneD, const DataType > &x) const
 
void GetPoints (Array< OneD, const DataType > &x, Array< OneD, const DataType > &y) const
 
void GetPoints (Array< OneD, const DataType > &x, Array< OneD, const DataType > &y, Array< OneD, const DataType > &z) const
 
const MatrixSharedPtrTypeGetD (Direction dir=xDir) const
 
const MatrixSharedPtrType GetI (const PointsKey &key)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x)
 
const MatrixSharedPtrType GetI (size_t uint, const Array< OneD, const DataType > &x)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y, const Array< OneD, const DataType > &z)
 
const MatrixSharedPtrType GetGalerkinProjection (const PointsKey &pkey)
 

Static Public Member Functions

static std::shared_ptr< Points< NekDouble > > Create (const PointsKey &pkey)
 

Protected Member Functions

const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const PointsKey &pkey) override
 
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const Array< OneD, const NekDouble > &x) override
 
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (size_t numpoints, const Array< OneD, const NekDouble > &x) override
 
const std::shared_ptr< NekMatrix< NekDouble > > v_GetGalerkinProjection (const PointsKey &pkey) override
 
- Protected Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual void v_Initialize (void)
 

Private Member Functions

 GaussPoints ()=delete
 
 GaussPoints (const GaussPoints &points)=delete
 
void v_CalculatePoints () final
 
void v_CalculateWeights () final
 
void v_CalculateDerivMatrix () final
 
void CalculateInterpMatrix (size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
 
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix (const PointsKey &pkey)
 

Static Private Attributes

static bool initPointsManager []
 

Additional Inherited Members

- Public Types inherited from Nektar::LibUtilities::Points< NekDouble >
typedef NekDouble DataType
 
typedef std::shared_ptr< NekMatrix< DataType > > MatrixSharedPtrType
 
- Protected Attributes inherited from Nektar::LibUtilities::Points< NekDouble >
PointsKey m_pointsKey
 Points type for this points distributions.
 
Array< OneD, DataTypem_points [3]
 Storage for the point locations, allowing for up to a 3D points storage.
 
Array< OneD, DataTypem_weights
 Quadrature weights for the weights.
 
Array< OneD, DataTypem_bcweights
 Barycentric weights.
 
MatrixSharedPtrType m_derivmatrix [3]
 Derivative matrices.
 
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLessm_InterpManager
 
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLessm_GalerkinProjectionManager
 

Detailed Description

Definition at line 42 of file GaussPoints.h.

Constructor & Destructor Documentation

◆ ~GaussPoints()

Nektar::LibUtilities::GaussPoints::~GaussPoints ( )
inlineoverride

Definition at line 45 of file GaussPoints.h.

46 {
47 }

◆ GaussPoints() [1/3]

Nektar::LibUtilities::GaussPoints::GaussPoints ( const PointsKey pkey)
inline

Definition at line 58 of file GaussPoints.h.

58 : PointsBaseType(pkey)
59 {
60 namespace pl = std::placeholders;
61 m_InterpManager.RegisterCreator(
62 PointsKey(0, eGaussGaussLegendre),
63 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
64 m_InterpManager.RegisterCreator(
65 PointsKey(0, eGaussRadauMLegendre),
66 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
67 m_InterpManager.RegisterCreator(
68 PointsKey(0, eGaussRadauPLegendre),
69 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
70 m_InterpManager.RegisterCreator(
71 PointsKey(0, eGaussLobattoLegendre),
72 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
73 m_InterpManager.RegisterCreator(
74 PointsKey(0, eGaussGaussChebyshev),
75 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
76 m_InterpManager.RegisterCreator(
77 PointsKey(0, eGaussRadauMChebyshev),
78 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
79 m_InterpManager.RegisterCreator(
80 PointsKey(0, eGaussRadauPChebyshev),
81 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
82 m_InterpManager.RegisterCreator(
83 PointsKey(0, eGaussLobattoChebyshev),
84 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
85 m_InterpManager.RegisterCreator(
86 PointsKey(0, eGaussRadauMAlpha0Beta1),
87 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
88 m_InterpManager.RegisterCreator(
89 PointsKey(0, eGaussRadauMAlpha0Beta2),
90 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
91 m_InterpManager.RegisterCreator(
92 PointsKey(0, eGaussRadauMAlpha1Beta0),
93 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
94 m_InterpManager.RegisterCreator(
95 PointsKey(0, eGaussRadauMAlpha2Beta0),
96 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
97 m_InterpManager.RegisterCreator(
98 PointsKey(0, eGaussKronrodLegendre),
99 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
100 m_InterpManager.RegisterCreator(
101 PointsKey(0, eGaussRadauKronrodMLegendre),
102 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
103 m_InterpManager.RegisterCreator(
104 PointsKey(0, eGaussRadauKronrodMAlpha1Beta0),
105 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
106 m_InterpManager.RegisterCreator(
107 PointsKey(0, eGaussLobattoKronrodLegendre),
108 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
109 m_InterpManager.RegisterCreator(
110 PointsKey(0, eFourierEvenlySpaced),
111 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
112 m_InterpManager.RegisterCreator(
113 PointsKey(0, ePolyEvenlySpaced),
114 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
115 m_InterpManager.RegisterCreator(
116 PointsKey(0, eBoundaryLayerPoints),
117 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
118 m_InterpManager.RegisterCreator(
119 PointsKey(0, eBoundaryLayerPointsRev),
120 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
121 m_InterpManager.RegisterCreator(
122 PointsKey(0, eGaussLegendreWithMP),
123 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
124 m_InterpManager.RegisterCreator(
125 PointsKey(0, eGaussLegendreWithM),
126 std::bind(&GaussPoints::CreateMatrix, this, pl::_1));
127 m_GalerkinProjectionManager.RegisterCreator(
128 PointsKey(0, eGaussGaussLegendre),
129 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
130 m_GalerkinProjectionManager.RegisterCreator(
131 PointsKey(0, eGaussRadauMLegendre),
132 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
133 m_GalerkinProjectionManager.RegisterCreator(
134 PointsKey(0, eGaussRadauPLegendre),
135 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
136 m_GalerkinProjectionManager.RegisterCreator(
137 PointsKey(0, eGaussLobattoLegendre),
138 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
139 m_GalerkinProjectionManager.RegisterCreator(
140 PointsKey(0, eGaussGaussChebyshev),
141 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
142 m_GalerkinProjectionManager.RegisterCreator(
143 PointsKey(0, eGaussRadauMChebyshev),
144 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
145 m_GalerkinProjectionManager.RegisterCreator(
146 PointsKey(0, eGaussRadauPChebyshev),
147 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
148 m_GalerkinProjectionManager.RegisterCreator(
149 PointsKey(0, eGaussLobattoChebyshev),
150 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
151 m_GalerkinProjectionManager.RegisterCreator(
152 PointsKey(0, eGaussRadauMAlpha0Beta1),
153 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
154 m_GalerkinProjectionManager.RegisterCreator(
155 PointsKey(0, eGaussRadauMAlpha0Beta2),
156 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
157 m_GalerkinProjectionManager.RegisterCreator(
158 PointsKey(0, eGaussRadauMAlpha1Beta0),
159 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
160 m_GalerkinProjectionManager.RegisterCreator(
161 PointsKey(0, eGaussRadauMAlpha2Beta0),
162 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
163 m_GalerkinProjectionManager.RegisterCreator(
164 PointsKey(0, eGaussKronrodLegendre),
165 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
166 m_GalerkinProjectionManager.RegisterCreator(
167 PointsKey(0, eGaussRadauKronrodMLegendre),
168 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
169 m_GalerkinProjectionManager.RegisterCreator(
170 PointsKey(0, eGaussRadauKronrodMAlpha1Beta0),
171 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
172 m_GalerkinProjectionManager.RegisterCreator(
173 PointsKey(0, eGaussLobattoKronrodLegendre),
174 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
175 m_GalerkinProjectionManager.RegisterCreator(
176 PointsKey(0, eFourierEvenlySpaced),
177 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
178 m_GalerkinProjectionManager.RegisterCreator(
179 PointsKey(0, ePolyEvenlySpaced),
180 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
181 m_GalerkinProjectionManager.RegisterCreator(
182 PointsKey(0, eBoundaryLayerPoints),
183 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
184 m_GalerkinProjectionManager.RegisterCreator(
185 PointsKey(0, eBoundaryLayerPointsRev),
186 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
187 m_GalerkinProjectionManager.RegisterCreator(
188 PointsKey(0, eGaussLegendreWithMP),
189 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
190 m_GalerkinProjectionManager.RegisterCreator(
191 PointsKey(0, eGaussLegendreWithM),
192 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
193 }
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix(const PointsKey &pkey)
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition Points.h:364
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition Points.h:366
Points< NekDouble > PointsBaseType
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition PointsType.h:47
@ eGaussLegendreWithMP
1D Gauss-Legendre quadrature points with additional x=-1 and x=1 end points
Definition PointsType.h:95
@ eBoundaryLayerPoints
1D power law distribution for boundary layer points
Definition PointsType.h:77
@ eGaussRadauMAlpha0Beta1
Gauss Radau pinned at x=-1,.
Definition PointsType.h:58
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition PointsType.h:67
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:57
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition PointsType.h:74
@ eGaussRadauPChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition PointsType.h:55
@ eGaussRadauMChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition PointsType.h:53
@ eGaussRadauKronrodMAlpha1Beta0
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition PointsType.h:69
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
@ eGaussLegendreWithM
1D Gauss-Legendre quadrature points with additional x=-1 point
Definition PointsType.h:97
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition PointsType.h:52
@ eBoundaryLayerPointsRev
1D power law distribution for boundary layer points
Definition PointsType.h:79
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition PointsType.h:73
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition PointsType.h:72
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition PointsType.h:46
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition PointsType.h:49

References CreateGPMatrix(), CreateMatrix(), Nektar::LibUtilities::eBoundaryLayerPoints, Nektar::LibUtilities::eBoundaryLayerPointsRev, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLegendreWithM, Nektar::LibUtilities::eGaussLegendreWithMP, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::LibUtilities::Points< NekDouble >::m_GalerkinProjectionManager, and Nektar::LibUtilities::Points< NekDouble >::m_InterpManager.

◆ GaussPoints() [2/3]

Nektar::LibUtilities::GaussPoints::GaussPoints ( )
privatedelete

◆ GaussPoints() [3/3]

Nektar::LibUtilities::GaussPoints::GaussPoints ( const GaussPoints points)
privatedelete

Member Function Documentation

◆ CalculateGalerkinProjectionMatrix()

std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CalculateGalerkinProjectionMatrix ( const PointsKey pkey)
private

Definition at line 466 of file GaussPoints.cpp.

468{
469 size_t numpointsfrom = pkey.GetNumPoints();
470 size_t numpointsto = GetNumPoints();
471
472 Array<OneD, const NekDouble> weightsfrom;
473
474 weightsfrom = PointsManager()[pkey]->GetW();
475
476 std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
477
478 Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
479
480 // set up inner product matrix and multiply by inverse of
481 // diagaonal mass matrix
482 for (size_t i = 0; i < numpointsto; ++i)
483 {
484 Vmath::Vmul(numpointsfrom, Interp->GetPtr().data() + i * numpointsfrom,
485 1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
486 if (m_weights[i] == 0.0)
487 {
488 Vmath::Zero(numpointsto, &GalProj[0] + i, numpointsto);
489 }
490 else
491 {
492 Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
493 numpointsto, &GalProj[0] + i, numpointsto);
494 }
495 }
496
497 NekDouble *t = GalProj.data();
498 std::shared_ptr<NekMatrix<NekDouble>> returnval(
499 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
500 numpointsto, numpointsfrom, t));
501
502 return returnval;
503}
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition Points.h:358
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition Points.h:317
PointsManagerT & PointsManager(void)
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition Vmath.hpp:72
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273

References Nektar::LibUtilities::Points< NekDouble >::GetI(), Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::m_weights, Nektar::LibUtilities::PointsManager(), Vmath::Smul(), Vmath::Vmul(), and Vmath::Zero().

Referenced by CreateGPMatrix().

◆ CalculateInterpMatrix()

void Nektar::LibUtilities::GaussPoints::CalculateInterpMatrix ( size_t  npts,
const Array< OneD, const NekDouble > &  xpoints,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 283 of file GaussPoints.cpp.

286{
287 switch (m_pointsKey.GetPointsType())
288 {
290 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
291 GetNumPoints(), npts, 0.0, 0.0);
292 break;
293
295 // first / last line of interp are zeros:
296 for (int i = 0; i < interp.size(); ++i)
297 {
298 interp[i] = 0.0;
299 }
300 // only interpolate from interior points
301 Polylib::Imgj(interp.data() + npts, m_points[0].data() + 1,
302 xpoints.data(), GetNumPoints() - 2, npts, 0.0, 0.0);
303 break;
304
306 // first / last line of interp are zeros:
307 for (int i = 0; i < interp.size(); ++i)
308 {
309 interp[i] = 0.0;
310 }
311 // only interpolate from interior points
312 Polylib::Imgj(interp.data() + npts, m_points[0].data() + 1,
313 xpoints.data(), GetNumPoints() - 1, npts, 0.0, 0.0);
314 break;
315
317 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
318 GetNumPoints(), npts, 0.0, 0.0);
319 break;
320
322 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
323 GetNumPoints(), npts, 0.0, 0.0);
324 break;
325
327 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
328 GetNumPoints(), npts, 0.0, 0.0);
329 break;
330
332 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
333 GetNumPoints(), npts, -0.5, -0.5);
334 break;
335
337 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
338 GetNumPoints(), npts, -0.5, -0.5);
339 break;
340
342 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
343 GetNumPoints(), npts, -0.5, -0.5);
344 break;
345
347 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
348 GetNumPoints(), npts, -0.5, -0.5);
349 break;
350
352 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
353 GetNumPoints(), npts, 0.0, 1.0);
354 break;
355
356 case eGaussRadauMAlpha0Beta2:
357 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
358 GetNumPoints(), npts, 0.0, 2.0);
359 break;
360
361 case eGaussRadauMAlpha1Beta0:
362 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
363 GetNumPoints(), npts, 1.0, 0.0);
364 break;
365
366 case eGaussRadauMAlpha2Beta0:
367 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
368 GetNumPoints(), npts, 2.0, 0.0);
369 break;
370
371 case eGaussKronrodLegendre:
375 {
376 for (size_t i = 0; i < npts; ++i)
377 {
378 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
379 {
380 interp[i + j * npts] =
381 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
383 }
384 }
385 }
386 break;
387
388 default:
390 "Unknown Gauss quadrature point distribution requested");
391 }
392}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition Points.h:356
PointsKey m_pointsKey
Points type for this points distributions.
Definition Points.h:353
PointsType GetPointsType() const
Definition Points.h:90
size_t GetNumPoints() const
Definition Points.h:85
double laginterp(double z, int j, const double *zj, int np)
Definition Polylib.cpp:85
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at ...
Definition Polylib.cpp:1116
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm.
Definition Polylib.cpp:1054
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at...
Definition Polylib.cpp:1085
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.
Definition Polylib.cpp:1147

References Nektar::Array< OneD, DataType >::data(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLegendreWithM, Nektar::LibUtilities::eGaussLegendreWithMP, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Polylib::Imgj(), Polylib::Imglj(), Polylib::Imgrjm(), Polylib::Imgrjp(), Polylib::laginterp(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, and NEKERROR.

Referenced by v_GetI().

◆ Create()

std::shared_ptr< Points< NekDouble > > Nektar::LibUtilities::GaussPoints::Create ( const PointsKey pkey)
static

Definition at line 394 of file GaussPoints.cpp.

395{
396 std::shared_ptr<Points<NekDouble>> returnval(
398
399 returnval->Initialize();
400
401 return returnval;
402}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ CreateGPMatrix()

std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CreateGPMatrix ( const PointsKey pkey)

Definition at line 456 of file GaussPoints.cpp.

458{
459 std::shared_ptr<NekMatrix<NekDouble>> returnval =
461
462 // Delegate to function below
463 return returnval;
464}
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix(const PointsKey &pkey)

References CalculateGalerkinProjectionMatrix().

Referenced by GaussPoints().

◆ CreateMatrix()

std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CreateMatrix ( const PointsKey pkey)

Definition at line 404 of file GaussPoints.cpp.

406{
407 size_t numpoints = pkey.GetNumPoints();
408 Array<OneD, const NekDouble> xpoints;
409
410 PointsManager()[pkey]->GetPoints(xpoints);
411
412 // Delegate to function below
413 return GetI(numpoints, xpoints);
414}

References Nektar::LibUtilities::Points< NekDouble >::GetI(), Nektar::LibUtilities::PointsKey::GetNumPoints(), and Nektar::LibUtilities::PointsManager().

Referenced by GaussPoints().

◆ v_CalculateDerivMatrix()

void Nektar::LibUtilities::GaussPoints::v_CalculateDerivMatrix ( )
finalprivate

Definition at line 193 of file GaussPoints.cpp.

194{
195 // Allocate the derivative matrix
196 PointsBaseType::v_CalculateDerivMatrix();
197
198 size_t numpoints = m_pointsKey.GetNumPoints();
199 size_t totpoints = m_pointsKey.GetTotNumPoints();
200 Array<OneD, NekDouble> dmtempSharedArray =
201 Array<OneD, NekDouble>(totpoints * totpoints);
202 NekDouble *dmtemp = dmtempSharedArray.data();
203
204 switch (m_pointsKey.GetPointsType())
205 {
207 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
208 break;
209
211 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
212 break;
213
215 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
216 break;
217
219 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
220 break;
221
223 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
224 break;
225
227 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
228 break;
229
231 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
232 break;
233
235 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
236 break;
237
239 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
240 break;
241
242 case eGaussRadauMAlpha0Beta2:
243 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
244 break;
245
246 case eGaussRadauMAlpha1Beta0:
247 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
248 break;
249
250 case eGaussRadauMAlpha2Beta0:
251 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
252 break;
253
256 case eGaussKronrodLegendre:
260 {
261 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
262 {
263 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
264 {
266 m_points[0][i], j, &m_points[0][0],
268 }
269 }
270 return;
271 }
272 break;
273
274 default:
276 "Unknown Gauss quadrature point distribution requested");
277 }
278
279 std::copy(dmtemp, dmtemp + totpoints * totpoints,
280 m_derivmatrix[0]->begin());
281}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition Points.h:362
size_t GetTotNumPoints() const
Definition Points.h:158
double laginterpderiv(double z, int k, const double *zj, int np)
Definition Polylib.cpp:98
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.
Definition Polylib.cpp:653
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
Definition Polylib.cpp:833
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a z...
Definition Polylib.cpp:704
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.
Definition Polylib.cpp:768

References Polylib::Dgj(), Polylib::Dglj(), Polylib::Dgrjm(), Polylib::Dgrjp(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLegendreWithM, Nektar::LibUtilities::eGaussLegendreWithMP, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Nektar::LibUtilities::PointsKey::GetTotNumPoints(), Polylib::laginterpderiv(), Nektar::LibUtilities::Points< NekDouble >::m_derivmatrix, Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, and NEKERROR.

◆ v_CalculatePoints()

void Nektar::LibUtilities::GaussPoints::v_CalculatePoints ( )
finalprivatevirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 79 of file GaussPoints.cpp.

80{
81 // Allocate the storage for points and weights
84
85 size_t numpoints = m_pointsKey.GetNumPoints();
86
87 switch (m_pointsKey.GetPointsType())
88 {
90 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
91 0.0);
92 break;
94 Polylib::zwgj(m_points[0].data() + 1, m_weights.data() + 1,
95 numpoints - 2, 0.0, 0.0);
96 m_points[0][0] = -1.0;
97 m_points[0][numpoints - 1] = 1.0;
98 m_weights[0] = 0.0;
99 m_weights[numpoints - 1] = 0.0;
100 break;
102 Polylib::zwgj(m_points[0].data() + 1, m_weights.data() + 1,
103 numpoints - 1, 0.0, 0.0);
104 m_points[0][0] = -1.0;
105 m_weights[0] = 0.0;
106 break;
108 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
109 0.0, 0.0);
110 break;
111
113 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
114 0.0, 0.0);
115 break;
116
118 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
119 0.0);
120 break;
121
123 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
124 -0.5);
125 break;
126
128 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
129 -0.5, -0.5);
130 break;
131
133 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
134 -0.5, -0.5);
135 break;
136
138 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
139 -0.5, -0.5);
140 break;
141
143 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
144 0.0, 1.0);
145 break;
146
147 case eGaussRadauMAlpha0Beta2:
148 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
149 0.0, 2.0);
150 break;
151
152 case eGaussRadauMAlpha1Beta0:
153 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
154 1.0, 0.0);
155 break;
156
157 case eGaussRadauMAlpha2Beta0:
158 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
159 2.0, 0.0);
160 break;
161
162 case eGaussKronrodLegendre:
163 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
164 0.0);
165 break;
166
168 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
169 0.0);
170 break;
171
173 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
174 0.0);
175 break;
176
178 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
179 0.0);
180 break;
181
182 default:
184 "Unknown Gauss quadrature point distribution requested");
185 }
186}
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition Polylib.cpp:182
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition Polylib.cpp:225
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition Polylib.cpp:266
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition Polylib.cpp:152
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition Polylib.cpp:317
void zwrk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Radau-Kronrod-Jacobi zeros and weights.
Definition Polylib.cpp:389
void zwlk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
Definition Polylib.cpp:495

References Nektar::Array< OneD, DataType >::data(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLegendreWithM, Nektar::LibUtilities::eGaussLegendreWithMP, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, Nektar::LibUtilities::Points< NekDouble >::m_weights, NEKERROR, Nektar::LibUtilities::Points< NekDouble >::v_CalculatePoints(), Nektar::LibUtilities::Points< NekDouble >::v_CalculateWeights(), Polylib::zwgj(), Polylib::zwgk(), Polylib::zwglj(), Polylib::zwgrjm(), Polylib::zwgrjp(), Polylib::zwlk(), and Polylib::zwrk().

◆ v_CalculateWeights()

void Nektar::LibUtilities::GaussPoints::v_CalculateWeights ( )
finalprivatevirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 188 of file GaussPoints.cpp.

189{
190 // For Gauss Quadrature, this is done as part of the points computation
191}

◆ v_GetGalerkinProjection()

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::v_GetGalerkinProjection ( const PointsKey pkey)
overrideprotected

◆ v_GetI() [1/3]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::v_GetI ( const Array< OneD, const NekDouble > &  x)
overrideprotected

Definition at line 425 of file GaussPoints.cpp.

427{
428 size_t numpoints = 1;
429
430 // Delegate to function below
431 return GetI(numpoints, x);
432}

References Nektar::LibUtilities::Points< NekDouble >::GetI().

◆ v_GetI() [2/3]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::v_GetI ( const PointsKey pkey)
overrideprotected

Definition at line 416 of file GaussPoints.cpp.

418{
419 ASSERTL0(pkey.GetPointsDim() == 1,
420 "Gauss Points can only interp to other 1d point distributions");
421
422 return m_InterpManager[pkey];
423}
#define ASSERTL0(condition, msg)

References ASSERTL0, Nektar::LibUtilities::PointsKey::GetPointsDim(), and Nektar::LibUtilities::Points< NekDouble >::m_InterpManager.

◆ v_GetI() [3/3]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::v_GetI ( size_t  numpoints,
const Array< OneD, const NekDouble > &  x 
)
overrideprotected

Definition at line 434 of file GaussPoints.cpp.

436{
437 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
438
439 CalculateInterpMatrix(numpoints, x, interp);
440
441 NekDouble *t = interp.data();
442 size_t np = GetNumPoints();
443 std::shared_ptr<NekMatrix<NekDouble>> returnval(
444 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
445 t));
446
447 return returnval;
448}
void CalculateInterpMatrix(size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)

References CalculateInterpMatrix(), and Nektar::LibUtilities::Points< NekDouble >::GetNumPoints().

Member Data Documentation

◆ initPointsManager

bool Nektar::LibUtilities::GaussPoints::initPointsManager
staticprivate

Definition at line 206 of file GaussPoints.h.