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

 ~PolyEPoints () override
 
 PolyEPoints (const PointsKey &key)
 
- 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< PointsBaseTypeCreate (const PointsKey &key)
 

Protected Member Functions

const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const PointsKey &pkey) override
 
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (const Array< OneD, const NekDouble > &x) override
 
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI (size_t numpoints, const Array< OneD, const NekDouble > &x) override
 
- 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

 PolyEPoints ()=delete
 
 PolyEPoints (const PolyEPoints &points)=delete
 
void v_CalculatePoints () final
 
void v_CalculateWeights () final
 
void v_CalculateDerivMatrix () final
 
void CalculateInterpMatrix (size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
 

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

Member Typedef Documentation

◆ PointsBaseType

Definition at line 47 of file PolyEPoints.h.

Constructor & Destructor Documentation

◆ ~PolyEPoints()

Nektar::LibUtilities::PolyEPoints::~PolyEPoints ( )
inlineoverride

Definition at line 49 of file PolyEPoints.h.

50 {
51 }

◆ PolyEPoints() [1/3]

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

Definition at line 56 of file PolyEPoints.h.

56 : PointsBaseType(key)
57 {
58 }
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:47

◆ PolyEPoints() [2/3]

Nektar::LibUtilities::PolyEPoints::PolyEPoints ( )
privatedelete

◆ PolyEPoints() [3/3]

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

Member Function Documentation

◆ CalculateInterpMatrix()

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

Definition at line 110 of file PolyEPoints.cpp.

113{
114 for (size_t i = 0; i < npts; ++i)
115 {
116 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
117 {
118 interp[i + j * npts] = Polylib::laginterp(
119 xpoints[i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
120 }
121 }
122}
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:356
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:353
size_t GetNumPoints() const
Definition: Points.h:85
double laginterp(double z, int j, const double *zj, int np)
Definition: Polylib.cpp:85

References Nektar::LibUtilities::PointsKey::GetNumPoints(), Polylib::laginterp(), 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 124 of file PolyEPoints.cpp.

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

◆ v_CalculateDerivMatrix()

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

Definition at line 95 of file PolyEPoints.cpp.

96{
97 // Allocate the derivative matrix.
98 PointsBaseType::v_CalculateDerivMatrix();
99
100 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
101 {
102 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
103 {
105 m_points[0][i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
106 }
107 }
108}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362
double laginterpderiv(double z, int k, const double *zj, int np)
Definition: Polylib.cpp:98

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

◆ v_CalculatePoints()

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

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

Definition at line 45 of file PolyEPoints.cpp.

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

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

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

Definition at line 65 of file PolyEPoints.cpp.

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

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Polylib::laginterp(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, Nektar::LibUtilities::Points< NekDouble >::m_weights, Nektar::LibUtilities::PointsManager(), Nektar::LibUtilities::Points< NekDouble >::v_CalculateWeights(), Nektar::UnitTests::w(), and Nektar::UnitTests::z().

◆ v_GetI() [1/3]

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

Delegate to function below.

Definition at line 149 of file PolyEPoints.cpp.

151{
152 size_t numpoints = 1;
153
154 /// Delegate to function below.
155 return GetI(numpoints, x);
156}
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:317

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

◆ v_GetI() [2/3]

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

Definition at line 135 of file PolyEPoints.cpp.

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

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 ( size_t  numpoints,
const Array< OneD, const NekDouble > &  x 
)
overrideprotected

Definition at line 158 of file PolyEPoints.cpp.

160{
161 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
162
163 CalculateInterpMatrix(numpoints, x, interp);
164
165 size_t np = GetTotNumPoints();
166 NekDouble *d = interp.data();
167 std::shared_ptr<NekMatrix<NekDouble>> returnval(
168 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
169 d));
170
171 return returnval;
172}
void CalculateInterpMatrix(size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
std::vector< double > d(NPUPPER *NPUPPER)

References CalculateInterpMatrix(), Nektar::UnitTests::d(), 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:168
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73

Definition at line 69 of file PolyEPoints.h.