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

virtual ~GaussPoints ()
 
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

virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const PointsKey &pkey) override
 
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const Array< OneD, const NekDouble > &x) override
 
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (size_t numpoints, const Array< OneD, const NekDouble > &x) override
 
virtual 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
 
virtual void v_CalculatePoints () override final
 
virtual void v_CalculateWeights () override final
 
virtual void v_CalculateDerivMatrix () override 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 44 of file GaussPoints.h.

Constructor & Destructor Documentation

◆ ~GaussPoints()

virtual Nektar::LibUtilities::GaussPoints::~GaussPoints ( )
inlinevirtual

Definition at line 47 of file GaussPoints.h.

48 {
49 }

◆ GaussPoints() [1/3]

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

Definition at line 60 of file GaussPoints.h.

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

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

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

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

269{
270 switch (m_pointsKey.GetPointsType())
271 {
273 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
274 GetNumPoints(), npts, 0.0, 0.0);
275 break;
276
278 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
279 GetNumPoints(), npts, 0.0, 0.0);
280 break;
281
283 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
284 GetNumPoints(), npts, 0.0, 0.0);
285 break;
286
288 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
289 GetNumPoints(), npts, 0.0, 0.0);
290 break;
291
293 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
294 GetNumPoints(), npts, -0.5, -0.5);
295 break;
296
298 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
299 GetNumPoints(), npts, -0.5, -0.5);
300 break;
301
303 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
304 GetNumPoints(), npts, -0.5, -0.5);
305 break;
306
308 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
309 GetNumPoints(), npts, -0.5, -0.5);
310 break;
311
313 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
314 GetNumPoints(), npts, 0.0, 1.0);
315 break;
316
317 case eGaussRadauMAlpha0Beta2:
318 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
319 GetNumPoints(), npts, 0.0, 2.0);
320 break;
321
322 case eGaussRadauMAlpha1Beta0:
323 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
324 GetNumPoints(), npts, 1.0, 0.0);
325 break;
326
327 case eGaussRadauMAlpha2Beta0:
328 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
329 GetNumPoints(), npts, 2.0, 0.0);
330 break;
331
332 case eGaussKronrodLegendre:
336 {
337 for (size_t i = 0; i < npts; ++i)
338 {
339 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
340 {
341 interp[i + j * npts] =
342 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
344 }
345 }
346 }
347 break;
348
349 default:
351 "Unknown Gauss quadrature point distribution requested");
352 }
353}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:361
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:358
PointsType GetPointsType() const
Definition: Points.h:95
size_t GetNumPoints() const
Definition: Points.h:90
double laginterp(double z, int j, const double *zj, int np)
Definition: Polylib.cpp:87
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:1074
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:1012
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:1043
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:1105

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

356{
357 std::shared_ptr<Points<NekDouble>> returnval(
359
360 returnval->Initialize();
361
362 return returnval;
363}
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 417 of file GaussPoints.cpp.

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

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

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 ( )
finaloverrideprivatevirtual

Definition at line 178 of file GaussPoints.cpp.

179{
180 // Allocate the derivative matrix
181 PointsBaseType::v_CalculateDerivMatrix();
182
183 size_t numpoints = m_pointsKey.GetNumPoints();
184 size_t totpoints = m_pointsKey.GetTotNumPoints();
185 Array<OneD, NekDouble> dmtempSharedArray =
186 Array<OneD, NekDouble>(totpoints * totpoints);
187 NekDouble *dmtemp = dmtempSharedArray.data();
188
189 switch (m_pointsKey.GetPointsType())
190 {
192 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
193 break;
194
196 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
197 break;
198
200 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
201 break;
202
204 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
205 break;
206
208 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
209 break;
210
212 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
213 break;
214
216 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
217 break;
218
220 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
221 break;
222
224 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
225 break;
226
227 case eGaussRadauMAlpha0Beta2:
228 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
229 break;
230
231 case eGaussRadauMAlpha1Beta0:
232 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
233 break;
234
235 case eGaussRadauMAlpha2Beta0:
236 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
237 break;
238
239 case eGaussKronrodLegendre:
243 {
244 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
245 {
246 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
247 {
249 m_points[0][i], j, &m_points[0][0],
251 }
252 }
253 return;
254 }
255 break;
256
257 default:
259 "Unknown Gauss quadrature point distribution requested");
260 }
261
262 std::copy(dmtemp, dmtemp + totpoints * totpoints,
263 m_derivmatrix[0]->begin());
264}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:367
size_t GetTotNumPoints() const
Definition: Points.h:163
def copy(self)
Definition: pycml.py:2663
double laginterpderiv(double z, int k, const double *zj, int np)
Definition: Polylib.cpp:100
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:647
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:807
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:694
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:750

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 ( )
finaloverrideprivatevirtual

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

Definition at line 77 of file GaussPoints.cpp.

78{
79 // Allocate the storage for points and weights
82
83 size_t numpoints = m_pointsKey.GetNumPoints();
84
85 switch (m_pointsKey.GetPointsType())
86 {
88 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
89 0.0);
90 break;
91
93 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
94 0.0, 0.0);
95 break;
96
98 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
99 0.0, 0.0);
100 break;
101
103 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
104 0.0);
105 break;
106
108 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
109 -0.5);
110 break;
111
113 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
114 -0.5, -0.5);
115 break;
116
118 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
119 -0.5, -0.5);
120 break;
121
123 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
124 -0.5, -0.5);
125 break;
126
128 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
129 0.0, 1.0);
130 break;
131
132 case eGaussRadauMAlpha0Beta2:
133 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
134 0.0, 2.0);
135 break;
136
137 case eGaussRadauMAlpha1Beta0:
138 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
139 1.0, 0.0);
140 break;
141
142 case eGaussRadauMAlpha2Beta0:
143 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
144 2.0, 0.0);
145 break;
146
147 case eGaussKronrodLegendre:
148 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
149 0.0);
150 break;
151
153 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
154 0.0);
155 break;
156
158 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
159 0.0);
160 break;
161
163 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
164 0.0);
165 break;
166
167 default:
169 "Unknown Gauss quadrature point distribution requested");
170 }
171}
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:223
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:262
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:154
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:311
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:383
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:489

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 ( )
finaloverrideprivatevirtual

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

Definition at line 173 of file GaussPoints.cpp.

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

◆ v_GetGalerkinProjection()

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

◆ v_GetI() [1/3]

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

Definition at line 386 of file GaussPoints.cpp.

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

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

◆ v_GetI() [2/3]

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

Definition at line 377 of file GaussPoints.cpp.

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

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 
)
overrideprotectedvirtual

Definition at line 395 of file GaussPoints.cpp.

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