Nektar++
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Attributes | List of all members
Nektar::LibUtilities::PolyEPoints Class Reference

#include <PolyEPoints.h>

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

Public Types

typedef Points< NekDoublePointsBaseType
 
- Public Types inherited from Nektar::LibUtilities::Points< NekDouble >
typedef NekDouble DataType
 
typedef std::shared_ptr< NekMatrix< DataType > > MatrixSharedPtrType
 

Public Member Functions

virtual ~PolyEPoints ()
 
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)
 
 PolyEPoints (const PointsKey &key)
 
- 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)
 
virtual const MatrixSharedPtrType GetGalerkinProjection (const PointsKey &pkey)
 

Static Public Member Functions

static std::shared_ptr< PointsBaseTypeCreate (const PointsKey &key)
 

Private Member Functions

 PolyEPoints ()
 Default constructor should not be called except by Create method. More...
 
 PolyEPoints (const PolyEPoints &points)
 Copy constructor should not be called. More...
 
void CalculatePoints ()
 
void CalculateWeights ()
 
void CalculateDerivMatrix ()
 
void CalculateInterpMatrix (unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
 
NekDouble LagrangeInterpolant (NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
 
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

- 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 53 of file PolyEPoints.h.

Member Typedef Documentation

◆ PointsBaseType

Definition at line 56 of file PolyEPoints.h.

Constructor & Destructor Documentation

◆ ~PolyEPoints()

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

Definition at line 58 of file PolyEPoints.h.

59  {
60  }

◆ PolyEPoints() [1/3]

Nektar::LibUtilities::PolyEPoints::PolyEPoints ( const PointsKey key)
inline

Definition at line 68 of file PolyEPoints.h.

68  :PointsBaseType(key)
69  {
70  }
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:56

◆ PolyEPoints() [2/3]

Nektar::LibUtilities::PolyEPoints::PolyEPoints ( )
inlineprivate

Default constructor should not be called except by Create method.

Definition at line 76 of file PolyEPoints.h.

77  {
78  }
static const PointsKey NullPointsKey(0, eNoPointsType)

◆ PolyEPoints() [3/3]

Nektar::LibUtilities::PolyEPoints::PolyEPoints ( const PolyEPoints points)
inlineprivate

Copy constructor should not be called.

Definition at line 81 of file PolyEPoints.h.

81  :PointsBaseType(points.m_pointsKey)
82  {
83  }

Member Function Documentation

◆ CalculateDerivMatrix()

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

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

Definition at line 98 of file PolyEPoints.cpp.

99  {
100  // Allocate the derivative matrix.
102 
103  for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
104  {
105  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
106  {
108  }
109  }
110  }
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 GetNumPoints() const
Definition: Points.h:107
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)

References Nektar::LibUtilities::Points< NekDouble >::CalculateDerivMatrix(), Nektar::LibUtilities::PointsKey::GetNumPoints(), LagrangePolyDeriv(), Nektar::LibUtilities::Points< NekDouble >::m_derivmatrix, Nektar::LibUtilities::Points< NekDouble >::m_points, and Nektar::LibUtilities::Points< NekDouble >::m_pointsKey.

◆ CalculateInterpMatrix()

void Nektar::LibUtilities::PolyEPoints::CalculateInterpMatrix ( unsigned int  npts,
const Array< OneD, const NekDouble > &  xpoints,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 112 of file PolyEPoints.cpp.

113  {
114  for(unsigned int i=0;i<npts;++i)
115  {
116  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
117  {
118  interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
119  }
120  }
121  }
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)

References Nektar::LibUtilities::PointsKey::GetNumPoints(), LagrangePoly(), Nektar::LibUtilities::Points< NekDouble >::m_points, and Nektar::LibUtilities::Points< NekDouble >::m_pointsKey.

Referenced by GetI().

◆ CalculatePoints()

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

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

Definition at line 49 of file PolyEPoints.cpp.

50  {
51  // Allocate the storage for points
53 
54  unsigned int npts = m_pointsKey.GetNumPoints();
55  if(npts==1)
56  {
57  m_points[0][0] = 0.0;
58  }
59  else
60  {
61  NekDouble dx = 2.0/(NekDouble)(npts-1);
62  for(unsigned int i=0;i<npts;++i)
63  {
64  m_points[0][i] = -1.0 + i*dx;
65  }
66  }
67  }
double NekDouble

References Nektar::LibUtilities::Points< NekDouble >::CalculatePoints(), Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::m_points, and Nektar::LibUtilities::Points< NekDouble >::m_pointsKey.

◆ CalculateWeights()

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

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

Definition at line 69 of file PolyEPoints.cpp.

70  {
71  // Allocate the storage for points
73 
74  unsigned int npts = m_pointsKey.GetNumPoints();
75  if(npts==1)
76  {
77  m_weights[0] = 2.0; //midpoint rule
78  }
79  else
80  {
81  PointsKey gaussKey(npts, eGaussLobattoLegendre);
82  std::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
83  Array<OneD, const NekDouble> z;
84  Array<OneD, const NekDouble> w;
85 
86  ptr->GetZW(z,w);
87  for(unsigned int i=0; i<npts;++i)
88  {
89  m_weights[i] = 0.0;
90  for(unsigned j=0;j<npts;++j)
91  {
92  m_weights[i] += w[j]*LagrangePoly(z[j],i,npts,m_points[0]);
93  }
94  }
95  }
96  }
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:392
PointsManagerT & PointsManager(void)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

References Nektar::LibUtilities::Points< NekDouble >::CalculateWeights(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), LagrangePoly(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, Nektar::LibUtilities::Points< NekDouble >::m_weights, and Nektar::LibUtilities::PointsManager().

◆ Create()

std::shared_ptr< PolyEPoints::PointsBaseType > Nektar::LibUtilities::PolyEPoints::Create ( const PointsKey key)
static

Definition at line 123 of file PolyEPoints.cpp.

124  {
125  std::shared_ptr<PointsBaseType> returnval(MemoryManager<PolyEPoints>::AllocateSharedPtr(key));
126 
127  returnval->Initialize();
128 
129  return returnval;
130  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ GetI() [1/3]

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

Delegate to function below.

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

Definition at line 144 of file PolyEPoints.cpp.

145  {
146  int numpoints = 1;
147 
148  /// Delegate to function below.
149  return GetI(numpoints, x);
150  }
const std::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)

References GetI().

◆ GetI() [2/3]

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

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

Definition at line 132 of file PolyEPoints.cpp.

133  {
134  ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
135 
136  int numpoints = pkey.GetNumPoints();
137  Array<OneD, const NekDouble> xpoints;
138 
139  PointsManager()[pkey]->GetPoints(xpoints);
140 
141  return GetI(numpoints, xpoints);
142  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216

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

Referenced by GetI().

◆ GetI() [3/3]

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

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

Definition at line 152 of file PolyEPoints.cpp.

153  {
154  Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
155 
156  CalculateInterpMatrix(numpoints, x, interp);
157 
158  unsigned int np = GetTotNumPoints();
159  NekDouble* d = interp.data();
160  std::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np, d));
161 
162  return returnval;
163  }
unsigned int GetNumPoints() const
Definition: Points.h:273
unsigned int GetTotNumPoints() const
Definition: Points.h:278
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)

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

◆ LagrangeInterpolant()

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

Definition at line 166 of file PolyEPoints.cpp.

168  {
169  NekDouble sum = 0.0;
170 
171  for(int i=0;i<npts;++i)
172  {
173  sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
174  }
175  return sum;
176  }

References LagrangePoly().

◆ LagrangePoly()

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

Definition at line 179 of file PolyEPoints.cpp.

180  {
181  NekDouble h=1.0;
182 
183  for(int i=0;i<pt; ++i)
184  {
185  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
186  }
187 
188  for(int i=pt+1;i<npts;++i)
189  {
190  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
191  }
192 
193  return h;
194  }

Referenced by CalculateInterpMatrix(), CalculateWeights(), and LagrangeInterpolant().

◆ LagrangePolyDeriv()

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

Definition at line 196 of file PolyEPoints.cpp.

197  {
198  NekDouble h;
199  NekDouble y=0.0;
200 
201  for(int j=0;j<npts;++j)
202  {
203  if(j!=pt)
204  {
205  h=1.0;
206  for(int i=0;i<npts;++i)
207  {
208  if(i!=pt)
209  {
210  if(i!=j)
211  {
212  h = h*(x-xpts[i]);
213  }
214  h = h/(xpts[pt]-xpts[i]);
215  }
216  }
217  y = y + h;
218  }
219  }
220  return y;
221  }

Referenced by CalculateDerivMatrix().

Member Data Documentation

◆ initPointsManager

bool Nektar::LibUtilities::PolyEPoints::initPointsManager
staticprivate
Initial value:
= {
}
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
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:64

Definition at line 73 of file PolyEPoints.h.