Nektar++
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)
 
virtual void v_CalculatePoints ()
 
virtual void v_CalculateWeights ()
 

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. More...
 
Array< OneD, DataTypem_points [3]
 Storage for the point locations, allowing for up to a 3D points storage. More...
 
Array< OneD, DataTypem_weights
 Quadrature weights for the weights. More...
 
Array< OneD, DataTypem_bcweights
 Barycentric weights. More...
 
MatrixSharedPtrType m_derivmatrix [3]
 Derivative matrices. More...
 
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_GalerkinProjectionManager.RegisterCreator(
122 PointsKey(0, eGaussGaussLegendre),
123 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
124 m_GalerkinProjectionManager.RegisterCreator(
125 PointsKey(0, eGaussRadauMLegendre),
126 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
127 m_GalerkinProjectionManager.RegisterCreator(
128 PointsKey(0, eGaussRadauPLegendre),
129 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
130 m_GalerkinProjectionManager.RegisterCreator(
131 PointsKey(0, eGaussLobattoLegendre),
132 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
133 m_GalerkinProjectionManager.RegisterCreator(
134 PointsKey(0, eGaussGaussChebyshev),
135 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
136 m_GalerkinProjectionManager.RegisterCreator(
137 PointsKey(0, eGaussRadauMChebyshev),
138 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
139 m_GalerkinProjectionManager.RegisterCreator(
140 PointsKey(0, eGaussRadauPChebyshev),
141 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
142 m_GalerkinProjectionManager.RegisterCreator(
143 PointsKey(0, eGaussLobattoChebyshev),
144 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
145 m_GalerkinProjectionManager.RegisterCreator(
146 PointsKey(0, eGaussRadauMAlpha0Beta1),
147 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
148 m_GalerkinProjectionManager.RegisterCreator(
149 PointsKey(0, eGaussRadauMAlpha0Beta2),
150 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
151 m_GalerkinProjectionManager.RegisterCreator(
152 PointsKey(0, eGaussRadauMAlpha1Beta0),
153 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
154 m_GalerkinProjectionManager.RegisterCreator(
155 PointsKey(0, eGaussRadauMAlpha2Beta0),
156 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
157 m_GalerkinProjectionManager.RegisterCreator(
158 PointsKey(0, eGaussKronrodLegendre),
159 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
160 m_GalerkinProjectionManager.RegisterCreator(
161 PointsKey(0, eGaussRadauKronrodMLegendre),
162 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
163 m_GalerkinProjectionManager.RegisterCreator(
164 PointsKey(0, eGaussRadauKronrodMAlpha1Beta0),
165 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
166 m_GalerkinProjectionManager.RegisterCreator(
167 PointsKey(0, eGaussLobattoKronrodLegendre),
168 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
169 m_GalerkinProjectionManager.RegisterCreator(
170 PointsKey(0, eFourierEvenlySpaced),
171 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
172 m_GalerkinProjectionManager.RegisterCreator(
173 PointsKey(0, ePolyEvenlySpaced),
174 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
175 m_GalerkinProjectionManager.RegisterCreator(
176 PointsKey(0, eBoundaryLayerPoints),
177 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
178 m_GalerkinProjectionManager.RegisterCreator(
179 PointsKey(0, eBoundaryLayerPointsRev),
180 std::bind(&GaussPoints::CreateGPMatrix, this, pl::_1));
181 }
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
@ 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
@ 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::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 425 of file GaussPoints.cpp.

427{
428 size_t numpointsfrom = pkey.GetNumPoints();
429 size_t numpointsto = GetNumPoints();
430
431 Array<OneD, const NekDouble> weightsfrom;
432
433 weightsfrom = PointsManager()[pkey]->GetW();
434
435 std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
436
437 Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
438
439 // set up inner product matrix and multiply by inverse of
440 // diagaonal mass matrix
441 for (size_t i = 0; i < numpointsto; ++i)
442 {
443 Vmath::Vmul(numpointsfrom, Interp->GetPtr().get() + i * numpointsfrom,
444 1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
445 Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
446 numpointsto, &GalProj[0] + i, numpointsto);
447 }
448
449 NekDouble *t = GalProj.data();
450 std::shared_ptr<NekMatrix<NekDouble>> returnval(
451 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
452 numpointsto, numpointsfrom, t));
453
454 return returnval;
455}
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)
double NekDouble
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

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(), and Vmath::Vmul().

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 264 of file GaussPoints.cpp.

267{
268 switch (m_pointsKey.GetPointsType())
269 {
271 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
272 GetNumPoints(), npts, 0.0, 0.0);
273 break;
274
276 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
277 GetNumPoints(), npts, 0.0, 0.0);
278 break;
279
281 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
282 GetNumPoints(), npts, 0.0, 0.0);
283 break;
284
286 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
287 GetNumPoints(), npts, 0.0, 0.0);
288 break;
289
291 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
292 GetNumPoints(), npts, -0.5, -0.5);
293 break;
294
296 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
297 GetNumPoints(), npts, -0.5, -0.5);
298 break;
299
301 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
302 GetNumPoints(), npts, -0.5, -0.5);
303 break;
304
306 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
307 GetNumPoints(), npts, -0.5, -0.5);
308 break;
309
311 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
312 GetNumPoints(), npts, 0.0, 1.0);
313 break;
314
315 case eGaussRadauMAlpha0Beta2:
316 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
317 GetNumPoints(), npts, 0.0, 2.0);
318 break;
319
320 case eGaussRadauMAlpha1Beta0:
321 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
322 GetNumPoints(), npts, 1.0, 0.0);
323 break;
324
325 case eGaussRadauMAlpha2Beta0:
326 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
327 GetNumPoints(), npts, 2.0, 0.0);
328 break;
329
330 case eGaussKronrodLegendre:
334 {
335 for (size_t i = 0; i < npts; ++i)
336 {
337 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
338 {
339 interp[i + j * npts] =
340 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
342 }
343 }
344 }
345 break;
346
347 default:
349 "Unknown Gauss quadrature point distribution requested");
350 }
351}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
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::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 353 of file GaussPoints.cpp.

354{
355 std::shared_ptr<Points<NekDouble>> returnval(
357
358 returnval->Initialize();
359
360 return returnval;
361}
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 415 of file GaussPoints.cpp.

417{
418 std::shared_ptr<NekMatrix<NekDouble>> returnval =
420
421 // Delegate to function below
422 return returnval;
423}
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 363 of file GaussPoints.cpp.

365{
366 size_t numpoints = pkey.GetNumPoints();
367 Array<OneD, const NekDouble> xpoints;
368
369 PointsManager()[pkey]->GetPoints(xpoints);
370
371 // Delegate to function below
372 return GetI(numpoints, xpoints);
373}

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 176 of file GaussPoints.cpp.

177{
178 // Allocate the derivative matrix
179 PointsBaseType::v_CalculateDerivMatrix();
180
181 size_t numpoints = m_pointsKey.GetNumPoints();
182 size_t totpoints = m_pointsKey.GetTotNumPoints();
183 Array<OneD, NekDouble> dmtempSharedArray =
184 Array<OneD, NekDouble>(totpoints * totpoints);
185 NekDouble *dmtemp = dmtempSharedArray.data();
186
187 switch (m_pointsKey.GetPointsType())
188 {
190 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
191 break;
192
194 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
195 break;
196
198 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
199 break;
200
202 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
203 break;
204
206 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
207 break;
208
210 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
211 break;
212
214 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
215 break;
216
218 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
219 break;
220
222 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
223 break;
224
225 case eGaussRadauMAlpha0Beta2:
226 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
227 break;
228
229 case eGaussRadauMAlpha1Beta0:
230 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
231 break;
232
233 case eGaussRadauMAlpha2Beta0:
234 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
235 break;
236
237 case eGaussKronrodLegendre:
241 {
242 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
243 {
244 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
245 {
247 m_points[0][i], j, &m_points[0][0],
249 }
250 }
251 return;
252 }
253 break;
254
255 default:
257 "Unknown Gauss quadrature point distribution requested");
258 }
259
260 std::copy(dmtemp, dmtemp + totpoints * totpoints,
261 m_derivmatrix[0]->begin());
262}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362
size_t GetTotNumPoints() const
Definition: Points.h:158
def copy(self)
Definition: pycml.py:2663
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 CellMLToNektar.pycml::copy(), Polylib::Dgj(), Polylib::Dglj(), Polylib::Dgrjm(), Polylib::Dgrjp(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, 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 75 of file GaussPoints.cpp.

76{
77 // Allocate the storage for points and weights
80
81 size_t numpoints = m_pointsKey.GetNumPoints();
82
83 switch (m_pointsKey.GetPointsType())
84 {
86 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
87 0.0);
88 break;
89
91 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
92 0.0, 0.0);
93 break;
94
96 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
97 0.0, 0.0);
98 break;
99
101 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
102 0.0);
103 break;
104
106 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
107 -0.5);
108 break;
109
111 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
112 -0.5, -0.5);
113 break;
114
116 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
117 -0.5, -0.5);
118 break;
119
121 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
122 -0.5, -0.5);
123 break;
124
126 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
127 0.0, 1.0);
128 break;
129
130 case eGaussRadauMAlpha0Beta2:
131 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
132 0.0, 2.0);
133 break;
134
135 case eGaussRadauMAlpha1Beta0:
136 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
137 1.0, 0.0);
138 break;
139
140 case eGaussRadauMAlpha2Beta0:
141 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
142 2.0, 0.0);
143 break;
144
145 case eGaussKronrodLegendre:
146 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
147 0.0);
148 break;
149
151 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
152 0.0);
153 break;
154
156 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
157 0.0);
158 break;
159
161 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
162 0.0);
163 break;
164
165 default:
167 "Unknown Gauss quadrature point distribution requested");
168 }
169}
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::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 171 of file GaussPoints.cpp.

172{
173 // For Gauss Quadrature, this is done as part of the points computation
174}

◆ 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 384 of file GaussPoints.cpp.

386{
387 size_t numpoints = 1;
388
389 // Delegate to function below
390 return GetI(numpoints, x);
391}

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 375 of file GaussPoints.cpp.

377{
378 ASSERTL0(pkey.GetPointsDim() == 1,
379 "Gauss Points can only interp to other 1d point distributions");
380
381 return m_InterpManager[pkey];
382}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

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 393 of file GaussPoints.cpp.

395{
396 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
397
398 CalculateInterpMatrix(numpoints, x, interp);
399
400 NekDouble *t = interp.data();
401 size_t np = GetNumPoints();
402 std::shared_ptr<NekMatrix<NekDouble>> returnval(
403 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
404 t));
405
406 return returnval;
407}
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 194 of file GaussPoints.h.