Nektar++
Public Member Functions | Static Public 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)
 
const std::shared_ptr< NekMatrix< NekDouble > > GetI (const PointsKey &pkey)
 
const std::shared_ptr< NekMatrix< NekDouble > > GetI (const Array< OneD, const NekDouble > &x)
 
const std::shared_ptr< NekMatrix< NekDouble > > GetI (unsigned int numpoints, const Array< OneD, const NekDouble > &x)
 
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix (const PointsKey &pkey)
 
const std::shared_ptr< NekMatrix< NekDouble > > GetGalerkinProjection (const PointsKey &pkey)
 
 GaussPoints (const PointsKey &pkey)
 
- Public Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual ~Points ()
 
virtual void Initialize (void)
 
unsigned int GetPointsDim () const
 
unsigned int GetNumPoints () const
 
unsigned int 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
 
virtual const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
 
virtual const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y, const Array< OneD, const DataType > &z)
 

Static Public Member Functions

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

Private Member Functions

 GaussPoints ()
 These should not be called. All creation is done using the constructor requiring the key, declared above. More...
 
 GaussPoints (const GaussPoints &points)
 
void CalculatePoints ()
 
void CalculateWeights ()
 
void CalculateDerivMatrix ()
 
void CalculateInterpMatrix (unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
 
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix (const PointsKey &pkey)
 
NekDouble LagrangeInterpolant (NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
 functions used by the Kronrod points More...
 
NekDouble LagrangePoly (NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
 
NekDouble LagrangePolyDeriv (NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
 

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 Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual void CalculateBaryWeights ()
 This function calculates the barycentric weights used for enhanced interpolation speed. More...
 
 Points (const PointsKey &key)
 
- 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 47 of file GaussPoints.h.

Constructor & Destructor Documentation

◆ ~GaussPoints()

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

Definition at line 50 of file GaussPoints.h.

51  {
52  }

◆ GaussPoints() [1/3]

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

Definition at line 68 of file GaussPoints.h.

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

References CreateGPMatrix(), CreateMatrix(), Nektar::LibUtilities::eBoundaryLayerPoints, Nektar::LibUtilities::eBoundaryLayerPointsRev, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, 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 ( )
private

These should not be called. All creation is done using the constructor requiring the key, declared above.

◆ GaussPoints() [3/3]

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

Member Function Documentation

◆ CalculateDerivMatrix()

void Nektar::LibUtilities::GaussPoints::CalculateDerivMatrix ( )
privatevirtual

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

Definition at line 154 of file GaussPoints.cpp.

155  {
156  // Allocate the derivative matrix
158 
159  int numpoints = m_pointsKey.GetNumPoints();
160  int totpoints = m_pointsKey.GetTotNumPoints();
161  Array<OneD, NekDouble> dmtempSharedArray = Array<OneD, NekDouble>(totpoints*totpoints);
162  NekDouble *dmtemp = dmtempSharedArray.data();
163 
164  switch(m_pointsKey.GetPointsType())
165  {
166  case eGaussGaussLegendre:
167  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
168  break;
169 
171  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
172  break;
173 
175  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
176  break;
177 
179  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
180  break;
181 
183  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
184  break;
185 
187  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
188  break;
189 
191  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
192  break;
193 
195  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
196  break;
197 
199  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,1.0);
200  break;
201 
203  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,2.0);
204  break;
205 
207  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,1.0,0.0);
208  break;
209 
211  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,2.0,0.0);
212  break;
213 
218  {
219  for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
220  {
221  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
222  {
224  }
225  }
226  return;
227  }
228  break;
229 
230  default:
232  "Unknown Gauss quadrature point distribution requested");
233  }
234 
235  std::copy(dmtemp,dmtemp+totpoints*totpoints,m_derivmatrix[0]->begin());
236  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:390
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:396
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:387
unsigned int GetTotNumPoints() const
Definition: Points.h:180
PointsType GetPointsType() const
Definition: Points.h:112
unsigned int GetNumPoints() const
Definition: Points.h:107
def copy(self)
Definition: pycml.py:2663
double NekDouble
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:546
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:690
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:589
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:639

References Nektar::LibUtilities::Points< NekDouble >::CalculateDerivMatrix(), CellMLToNektar.pycml::copy(), Polylib::Dgj(), Polylib::Dglj(), Polylib::Dgrjm(), Polylib::Dgrjp(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Nektar::LibUtilities::PointsKey::GetTotNumPoints(), LagrangePolyDeriv(), Nektar::LibUtilities::Points< NekDouble >::m_derivmatrix, Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, and NEKERROR.

◆ CalculateGalerkinProjectionMatrix()

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

Definition at line 430 of file GaussPoints.cpp.

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

References 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 ( unsigned int  npts,
const Array< OneD, const NekDouble > &  xpoints,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 238 of file GaussPoints.cpp.

239  {
240  switch(m_pointsKey.GetPointsType())
241  {
242  case eGaussGaussLegendre:
243  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
244  break;
245 
247  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
248  break;
249 
251  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
252  break;
253 
255  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
256  break;
257 
259  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
260  break;
261 
263  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
264  break;
265 
267  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
268  break;
269 
271  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
272  break;
273 
275  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,1.0);
276  break;
277 
279  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,2.0);
280  break;
281 
283  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,1.0,0.0);
284  break;
285 
287  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,2.0,0.0);
288  break;
289 
294  {
295  for(unsigned int i=0;i<npts;++i)
296  {
297  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
298  {
299  interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
300  }
301  }
302  }
303  break;
304 
305  default:
307  "Unknown Gauss quadrature point distribution requested");
308  }
309  }
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
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:946
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:887
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:916
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:977

References Nektar::Array< OneD, DataType >::data(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, 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(), LagrangePoly(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, and NEKERROR.

Referenced by GetI().

◆ CalculatePoints()

void Nektar::LibUtilities::GaussPoints::CalculatePoints ( )
privatevirtual

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

Definition at line 69 of file GaussPoints.cpp.

70  {
71  // Allocate the storage for points and weights
74 
75  int numpoints = m_pointsKey.GetNumPoints();
76 
77  switch(m_pointsKey.GetPointsType())
78  {
80  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
81  break;
82 
84  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
85  break;
86 
88  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
89  break;
90 
92  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
93  break;
94 
96  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
97  break;
98 
100  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
101  break;
102 
104  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
105  break;
106 
108  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
109  break;
110 
112  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,1.0);
113  break;
114 
116  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,2.0);
117  break;
118 
120  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
121  break;
122 
124  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,2.0,0.0);
125  break;
126 
128  Polylib::zwgk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
129  break;
130 
132  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
133  break;
134 
136  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
137  break;
138 
140  Polylib::zwlk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
141  break;
142 
143  default:
145  "Unknown Gauss quadrature point distribution requested");
146  }
147  }
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:119
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:158
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:195
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:91
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:240
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:313
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:420

References Nektar::LibUtilities::Points< NekDouble >::CalculatePoints(), Nektar::LibUtilities::Points< NekDouble >::CalculateWeights(), Nektar::Array< OneD, DataType >::data(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, 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, Polylib::zwgj(), Polylib::zwgk(), Polylib::zwglj(), Polylib::zwgrjm(), Polylib::zwgrjp(), Polylib::zwlk(), and Polylib::zwrk().

◆ CalculateWeights()

void Nektar::LibUtilities::GaussPoints::CalculateWeights ( )
privatevirtual

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

Definition at line 149 of file GaussPoints.cpp.

150  {
151  // For Gauss Quadrature, this is done as part of the points computation
152  }

◆ Create()

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

Definition at line 311 of file GaussPoints.cpp.

312  {
313  std::shared_ptr< Points<NekDouble> > returnval(MemoryManager< GaussPoints >::AllocateSharedPtr(pkey));
314 
315  returnval->Initialize();
316 
317  return returnval;
318  }
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 422 of file GaussPoints.cpp.

423  {
424  std::shared_ptr< NekMatrix<NekDouble> > returnval = CalculateGalerkinProjectionMatrix(pkey);
425 
426  // Delegate to function below
427  return returnval;
428  }
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 320 of file GaussPoints.cpp.

321  {
322  int numpoints = pkey.GetNumPoints();
323  Array<OneD, const NekDouble> xpoints;
324 
325  PointsManager()[pkey]->GetPoints(xpoints);
326 
327  // Delegate to function below
328  return GetI(numpoints, xpoints);
329  }

References GetI(), Nektar::LibUtilities::PointsKey::GetNumPoints(), and Nektar::LibUtilities::PointsManager().

Referenced by GaussPoints().

◆ GetGalerkinProjection()

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

◆ GetI() [1/3]

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

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

Definition at line 339 of file GaussPoints.cpp.

340  {
341  int numpoints = 1;
342 
343  // Delegate to function below
344  return GetI(numpoints, x);
345  }

References GetI().

◆ GetI() [2/3]

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

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

Definition at line 332 of file GaussPoints.cpp.

333  {
334  ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
335 
336  return m_InterpManager[pkey];
337  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216

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

Referenced by CalculateGalerkinProjectionMatrix(), CreateMatrix(), and GetI().

◆ GetI() [3/3]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetI ( unsigned int  numpoints,
const Array< OneD, const NekDouble > &  x 
)
virtual

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

Definition at line 347 of file GaussPoints.cpp.

348  {
349  Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
350 
351  CalculateInterpMatrix(numpoints, x, interp);
352 
353  NekDouble* t = interp.data();
354  unsigned int np = GetNumPoints();
355  std::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np,t));
356 
357  return returnval;
358  }
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)

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

◆ LagrangeInterpolant()

NekDouble Nektar::LibUtilities::GaussPoints::LagrangeInterpolant ( NekDouble  x,
int  npts,
const Array< OneD, const NekDouble > &  xpts,
const Array< OneD, const NekDouble > &  funcvals 
)
private

functions used by the Kronrod points

Definition at line 360 of file GaussPoints.cpp.

362  {
363  NekDouble sum = 0.0;
364 
365  for(int i=0;i<npts;++i)
366  {
367  sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
368  }
369  return sum;
370  }

References LagrangePoly().

◆ LagrangePoly()

NekDouble Nektar::LibUtilities::GaussPoints::LagrangePoly ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 373 of file GaussPoints.cpp.

374  {
375  NekDouble h=1.0;
376 
377  for(int i=0;i<pt; ++i)
378  {
379  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
380  }
381 
382  for(int i=pt+1;i<npts;++i)
383  {
384  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
385  }
386 
387  return h;
388  }

Referenced by CalculateInterpMatrix(), and LagrangeInterpolant().

◆ LagrangePolyDeriv()

NekDouble Nektar::LibUtilities::GaussPoints::LagrangePolyDeriv ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 390 of file GaussPoints.cpp.

391  {
392  NekDouble h;
393  NekDouble y=0.0;
394 
395  for(int j=0;j<npts;++j)
396  {
397  if(j!=pt)
398  {
399  h=1.0;
400  for(int i=0;i<npts;++i)
401  {
402  if(i!=pt)
403  {
404  if(i!=j)
405  {
406  h = h*(x-xpts[i]);
407  }
408  h = h/(xpts[pt]-xpts[i]);
409  }
410  }
411  y = y + h;
412  }
413  }
414  return y;
415  }

Referenced by CalculateDerivMatrix().

Member Data Documentation

◆ initPointsManager

bool Nektar::LibUtilities::GaussPoints::initPointsManager
staticprivate
Initial value:
= {
}
static std::shared_ptr< Points< NekDouble > > Create(const PointsKey &pkey)
bool RegisterCreator(const KeyType &key, const CreateFuncType &createFunc)
Register the given function and associate it with the key. The return value is just to facilitate cal...
Definition: NekManager.hpp:166

Definition at line 155 of file GaussPoints.h.