Nektar++
Public Types | Public Member Functions | Static Public Member Functions | Protected 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 ()
 
 PolyEPoints (const PointsKey &key)
 
- Public Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual ~Points ()
 
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
 
const MatrixSharedPtrType GetI (const PointsKey &key)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x)
 
const MatrixSharedPtrType GetI (unsigned int 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< PointsBaseTypeCreate (const PointsKey &key)
 

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 (unsigned int numpoints, const Array< OneD, const NekDouble > &x) override
 
- Protected Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual void v_Initialize (void)
 
virtual void v_CalculateBaryWeights ()
 This function calculates the barycentric weights used for enhanced interpolation speed. More...
 
 Points (const PointsKey &key)
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y, const Array< OneD, const DataType > &z)
 
virtual const MatrixSharedPtrType v_GetGalerkinProjection (const PointsKey &pkey)
 

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...
 
virtual void v_CalculatePoints () override
 
virtual void v_CalculateWeights () override
 
virtual void v_CalculateDerivMatrix () override
 
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 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 51 of file PolyEPoints.h.

Member Typedef Documentation

◆ PointsBaseType

Definition at line 54 of file PolyEPoints.h.

Constructor & Destructor Documentation

◆ ~PolyEPoints()

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

Definition at line 56 of file PolyEPoints.h.

57  {
58  }

◆ PolyEPoints() [1/3]

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

Definition at line 63 of file PolyEPoints.h.

63  : PointsBaseType(key)
64  {
65  }
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:54

◆ PolyEPoints() [2/3]

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

Default constructor should not be called except by Create method.

Definition at line 79 of file PolyEPoints.h.

80  {
81  }
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 84 of file PolyEPoints.h.

84  : PointsBaseType(points.m_pointsKey)
85  {
86  }

Member Function Documentation

◆ 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.

115 {
116  for (unsigned int i = 0; i < npts; ++i)
117  {
118  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
119  {
120  interp[i + j * npts] = LagrangePoly(
121  xpoints[i], j, m_pointsKey.GetNumPoints(), m_points[0]);
122  }
123  }
124 }
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:375
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:372
unsigned int GetNumPoints() const
Definition: Points.h:104
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 v_GetI().

◆ Create()

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

Definition at line 126 of file PolyEPoints.cpp.

128 {
129  std::shared_ptr<PointsBaseType> returnval(
131 
132  returnval->Initialize();
133 
134  return returnval;
135 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

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

179 {
180  NekDouble sum = 0.0;
181 
182  for (int i = 0; i < npts; ++i)
183  {
184  sum += funcvals[i] * LagrangePoly(x, i, npts, xpts);
185  }
186  return sum;
187 }
double NekDouble

References LagrangePoly().

◆ LagrangePoly()

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

Definition at line 189 of file PolyEPoints.cpp.

191 {
192  NekDouble h = 1.0;
193 
194  for (int i = 0; i < pt; ++i)
195  {
196  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
197  }
198 
199  for (int i = pt + 1; i < npts; ++i)
200  {
201  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
202  }
203 
204  return h;
205 }

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

◆ LagrangePolyDeriv()

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

Definition at line 207 of file PolyEPoints.cpp.

209 {
210  NekDouble h;
211  NekDouble y = 0.0;
212 
213  for (int j = 0; j < npts; ++j)
214  {
215  if (j != pt)
216  {
217  h = 1.0;
218  for (int i = 0; i < npts; ++i)
219  {
220  if (i != pt)
221  {
222  if (i != j)
223  {
224  h = h * (x - xpts[i]);
225  }
226  h = h / (xpts[pt] - xpts[i]);
227  }
228  }
229  y = y + h;
230  }
231  }
232  return y;
233 }

Referenced by v_CalculateDerivMatrix().

◆ v_CalculateDerivMatrix()

void Nektar::LibUtilities::PolyEPoints::v_CalculateDerivMatrix ( )
overrideprivatevirtual

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

Definition at line 97 of file PolyEPoints.cpp.

98 {
99  // Allocate the derivative matrix.
101 
102  for (unsigned int i = 0; i < m_pointsKey.GetNumPoints(); ++i)
103  {
104  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
105  {
106  (*m_derivmatrix[0])(i, j) = LagrangePolyDeriv(
107  m_points[0][i], j, m_pointsKey.GetNumPoints(), m_points[0]);
108  }
109  }
110 }
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:381
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)

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

◆ v_CalculatePoints()

void Nektar::LibUtilities::PolyEPoints::v_CalculatePoints ( )
overrideprivatevirtual

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

Definition at line 48 of file PolyEPoints.cpp.

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

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

◆ v_CalculateWeights()

void Nektar::LibUtilities::PolyEPoints::v_CalculateWeights ( )
overrideprivatevirtual

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

Definition at line 68 of file PolyEPoints.cpp.

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

References 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, Nektar::LibUtilities::PointsManager(), and Nektar::LibUtilities::Points< NekDouble >::v_CalculateWeights().

◆ v_GetI() [1/3]

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

Delegate to function below.

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

Definition at line 151 of file PolyEPoints.cpp.

153 {
154  int numpoints = 1;
155 
156  /// Delegate to function below.
157  return GetI(numpoints, x);
158 }
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:336

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

◆ v_GetI() [2/3]

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

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

Definition at line 137 of file PolyEPoints.cpp.

139 {
140  ASSERTL0(pkey.GetPointsDim() == 1,
141  "Gauss Points can only interp to other 1d point distributions");
142 
143  int numpoints = pkey.GetNumPoints();
144  Array<OneD, const NekDouble> xpoints;
145 
146  PointsManager()[pkey]->GetPoints(xpoints);
147 
148  return GetI(numpoints, xpoints);
149 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

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

◆ v_GetI() [3/3]

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

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

Definition at line 160 of file PolyEPoints.cpp.

162 {
163  Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
164 
165  CalculateInterpMatrix(numpoints, x, interp);
166 
167  unsigned int np = GetTotNumPoints();
168  NekDouble *d = interp.data();
169  std::shared_ptr<NekMatrix<NekDouble>> returnval(
170  MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
171  d));
172 
173  return returnval;
174 }
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().

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:170
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75

Definition at line 76 of file PolyEPoints.h.