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

#include <FourierPoints.h>

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

Public Member Functions

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

 FourierPoints ()=delete
 
 FourierPoints (const FourierPoints &points)=delete
 
void v_CalculatePoints () override
 
void v_CalculateWeights () override
 
void v_CalculateDerivMatrix () override
 
void CalculateInterpMatrix (size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
 
NekDouble PeriodicSincFunction (const NekDouble x, const NekDouble h)
 

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 42 of file FourierPoints.h.

Constructor & Destructor Documentation

◆ ~FourierPoints()

Nektar::LibUtilities::FourierPoints::~FourierPoints ( )
inlineoverride

Definition at line 45 of file FourierPoints.h.

46 {
47 }

◆ FourierPoints() [1/3]

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

Definition at line 54 of file FourierPoints.h.

54 : PointsBaseType(key)
55 {
56 namespace pl = std::placeholders;
57 m_InterpManager.RegisterCreator(
58 PointsKey(0, eGaussGaussLegendre),
59 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
60 m_InterpManager.RegisterCreator(
61 PointsKey(0, eGaussRadauMLegendre),
62 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
63 m_InterpManager.RegisterCreator(
64 PointsKey(0, eGaussRadauPLegendre),
65 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
66 m_InterpManager.RegisterCreator(
67 PointsKey(0, eGaussLobattoLegendre),
68 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
69 m_InterpManager.RegisterCreator(
70 PointsKey(0, eGaussGaussChebyshev),
71 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
72 m_InterpManager.RegisterCreator(
73 PointsKey(0, eGaussRadauMChebyshev),
74 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
75 m_InterpManager.RegisterCreator(
76 PointsKey(0, eGaussRadauPChebyshev),
77 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
78 m_InterpManager.RegisterCreator(
79 PointsKey(0, eGaussLobattoChebyshev),
80 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
81 m_InterpManager.RegisterCreator(
82 PointsKey(0, eGaussRadauMAlpha0Beta1),
83 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
84 m_InterpManager.RegisterCreator(
85 PointsKey(0, eGaussRadauMAlpha0Beta2),
86 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
87 m_InterpManager.RegisterCreator(
88 PointsKey(0, eGaussRadauMAlpha1Beta0),
89 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
90 m_InterpManager.RegisterCreator(
91 PointsKey(0, eGaussRadauMAlpha2Beta0),
92 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
93 m_InterpManager.RegisterCreator(
94 PointsKey(0, ePolyEvenlySpaced),
95 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
96 m_InterpManager.RegisterCreator(
97 PointsKey(0, eFourierEvenlySpaced),
98 std::bind(&FourierPoints::CreateMatrix, this, pl::_1));
99 }
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:364
Points< NekDouble > PointsBaseType
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:47
@ eGaussRadauMAlpha0Beta1
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:58
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:57
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eGaussRadauPChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:55
@ eGaussRadauMChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:53
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:52
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:49

References CreateMatrix(), Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::ePolyEvenlySpaced, and Nektar::LibUtilities::Points< NekDouble >::m_InterpManager.

◆ FourierPoints() [2/3]

Nektar::LibUtilities::FourierPoints::FourierPoints ( )
privatedelete

◆ FourierPoints() [3/3]

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

Member Function Documentation

◆ CalculateInterpMatrix()

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

Definition at line 156 of file FourierPoints.cpp.

159{
160 const NekDouble h = 2.0 / m_pointsKey.GetNumPoints();
161 for (size_t i = 0; i < npts; ++i)
162 {
163 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
164 {
165 interp[i * m_pointsKey.GetNumPoints() + j] =
166 PeriodicSincFunction(M_PI * (xpoints[i] - m_points[0][j]), h);
167 }
168 }
169}
NekDouble PeriodicSincFunction(const NekDouble x, const NekDouble h)
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 NekDouble

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

Referenced by v_GetI().

◆ Create()

std::shared_ptr< Points< NekDouble > > Nektar::LibUtilities::FourierPoints::Create ( const PointsKey key)
static

Definition at line 100 of file FourierPoints.cpp.

101{
102 std::shared_ptr<Points<NekDouble>> returnval(
104
105 returnval->Initialize();
106
107 return returnval;
108}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ CreateMatrix()

std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::FourierPoints::CreateMatrix ( const PointsKey pkey)

Delegate to function below.

Definition at line 110 of file FourierPoints.cpp.

112{
113 size_t numpoints = pkey.GetNumPoints();
114 Array<OneD, const NekDouble> xpoints;
115
116 PointsManager()[pkey]->GetPoints(xpoints);
117
118 /// Delegate to function below.
119 return GetI(numpoints, xpoints);
120}
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:317
PointsManagerT & PointsManager(void)

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

Referenced by FourierPoints().

◆ PeriodicSincFunction()

NekDouble Nektar::LibUtilities::FourierPoints::PeriodicSincFunction ( const NekDouble  x,
const NekDouble  h 
)
private

Definition at line 171 of file FourierPoints.cpp.

173{
174 // This formula is based upon a mapped version of
175 // the periodic sinc presented in Trefethen's "Spectral Methods
176 // in Matlab"
177
178 NekDouble y = 1.0;
179
180 if (fabs(x) > 1.0e-12)
181 {
182 y = sin(M_PI * x / (M_PI * h)) / ((2.0 / h) * tan(0.5 * x));
183 }
184
185 return y;
186}

Referenced by CalculateInterpMatrix().

◆ v_CalculateDerivMatrix()

void Nektar::LibUtilities::FourierPoints::v_CalculateDerivMatrix ( )
overrideprivate

Definition at line 70 of file FourierPoints.cpp.

71{
72 //// Allocate the derivative matrix.
73 Points<NekDouble>::v_CalculateDerivMatrix();
74
75 size_t npts = m_pointsKey.GetNumPoints();
76
77 for (size_t i = 1; i < npts; ++i)
78 {
79 m_derivmatrix[0]->SetValue(i, i, 0.0);
80 }
81
82 for (size_t i = 1; i < npts; ++i)
83 {
84 m_derivmatrix[0]->SetValue(0, i,
85 -0.5 * M_PI * pow(-1.0, NekDouble(i)) *
86 cos(M_PI * i / npts) /
87 sin(M_PI * i / npts));
88 }
89
90 for (size_t i = 1; i < npts; ++i)
91 {
92 for (size_t j = 0; j < npts; ++j)
93 {
94 m_derivmatrix[0]->SetValue(
95 i, j, (*m_derivmatrix[0])(0, (j - i + npts) % npts));
96 }
97 }
98}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362

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

◆ v_CalculatePoints()

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

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

Definition at line 43 of file FourierPoints.cpp.

44{
45 // Allocate the storage for points
47
48 size_t npts = m_pointsKey.GetNumPoints();
49 ASSERTL0(!(npts % 2), "Fourier points need to be of even order");
50
51 NekDouble dx = 2.0 / (NekDouble)(npts);
52 for (size_t i = 0; i < npts; ++i)
53 {
54 m_points[0][i] = -1.0 + i * dx;
55 }
56}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, 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::FourierPoints::v_CalculateWeights ( )
overrideprivatevirtual

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

Definition at line 58 of file FourierPoints.cpp.

59{
60 // Allocate the storage for points
62
63 size_t npts = m_pointsKey.GetNumPoints();
64 for (size_t i = 0; i < npts; ++i)
65 {
66 m_weights[i] = 2.0 / (NekDouble)npts;
67 }
68}
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:358

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

◆ v_GetI() [1/3]

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

Delegate to function below.

Definition at line 131 of file FourierPoints.cpp.

133{
134 size_t numpoints = 1;
135
136 /// Delegate to function below.
137 return GetI(numpoints, x);
138}

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

◆ v_GetI() [2/3]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::FourierPoints::v_GetI ( const PointsKey pkey)
finalprotected

Definition at line 122 of file FourierPoints.cpp.

124{
125 ASSERTL0(pkey.GetPointsDim() == 1,
126 "Fourier Points can only interp to other 1d point distributions");
127
128 return m_InterpManager[pkey];
129}

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::FourierPoints::v_GetI ( size_t  numpoints,
const Array< OneD, const NekDouble > &  x 
)
finalprotected

Definition at line 140 of file FourierPoints.cpp.

142{
143 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
144
145 CalculateInterpMatrix(numpoints, x, interp);
146
147 NekDouble *t = interp.data();
148 size_t np = GetNumPoints();
149 std::shared_ptr<NekMatrix<NekDouble>> returnval(
150 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
151 t));
152
153 return returnval;
154}
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::FourierPoints::initPointsManager
staticprivate
Initial value:
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
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

Definition at line 110 of file FourierPoints.h.