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

virtual ~FourierPoints ()
 
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)
 
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 MatrixSharedPtrType v_GetI (const PointsKey &pkey) override
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const NekDouble > &x) override
 
virtual const MatrixSharedPtrType 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

 FourierPoints ()
 Default constructor should not be called except by Create method. More...
 
 FourierPoints (const FourierPoints &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 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 47 of file FourierPoints.h.

Constructor & Destructor Documentation

◆ ~FourierPoints()

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

Definition at line 50 of file FourierPoints.h.

51  {
52  }

◆ FourierPoints() [1/3]

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

Definition at line 59 of file FourierPoints.h.

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

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

Default constructor should not be called except by Create method.

◆ FourierPoints() [3/3]

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

Copy constructor should not be called.

Member Function Documentation

◆ CalculateInterpMatrix()

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

Definition at line 176 of file FourierPoints.cpp.

179 {
180  const NekDouble h = 2.0 / m_pointsKey.GetNumPoints();
181  for (unsigned int i = 0; i < npts; ++i)
182  {
183  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
184  {
185  interp[i * m_pointsKey.GetNumPoints() + j] =
186  PeriodicSincFunction(M_PI * (xpoints[i] - m_points[0][j]), h);
187  }
188  }
189 }
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:375
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:372
unsigned int GetNumPoints() const
Definition: Points.h:104
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 120 of file FourierPoints.cpp.

121 {
122  std::shared_ptr<Points<NekDouble>> returnval(
124 
125  returnval->Initialize();
126 
127  return returnval;
128 }
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 130 of file FourierPoints.cpp.

132 {
133  int numpoints = pkey.GetNumPoints();
134  Array<OneD, const NekDouble> xpoints;
135 
136  PointsManager()[pkey]->GetPoints(xpoints);
137 
138  /// Delegate to function below.
139  return GetI(numpoints, xpoints);
140 }
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:336
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 191 of file FourierPoints.cpp.

193 {
194  // This formula is based upon a mapped version of
195  // the periodic sinc presented in Trefethen's "Spectral Methods
196  // in Matlab"
197 
198  NekDouble y = 1.0;
199 
200  if (fabs(x) > 1.0e-12)
201  {
202  y = sin(M_PI * x / (M_PI * h)) / ((2.0 / h) * tan(0.5 * x));
203  }
204 
205  return y;
206 }

Referenced by CalculateInterpMatrix().

◆ v_CalculateDerivMatrix()

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

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

Definition at line 90 of file FourierPoints.cpp.

91 {
92  //// Allocate the derivative matrix.
94 
95  unsigned int npts = m_pointsKey.GetNumPoints();
96 
97  for (unsigned int i = 1; i < npts; ++i)
98  {
99  m_derivmatrix[0]->SetValue(i, i, 0.0);
100  }
101 
102  for (unsigned int i = 1; i < npts; ++i)
103  {
104  m_derivmatrix[0]->SetValue(0, i,
105  -0.5 * M_PI * pow(-1.0, NekDouble(i)) *
106  cos(M_PI * i / npts) /
107  sin(M_PI * i / npts));
108  }
109 
110  for (unsigned int i = 1; i < npts; ++i)
111  {
112  for (unsigned int j = 0; j < npts; ++j)
113  {
114  m_derivmatrix[0]->SetValue(
115  i, j, (*m_derivmatrix[0])(0, (j - i + npts) % npts));
116  }
117  }
118 }
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:381
virtual void v_CalculateDerivMatrix()
Definition: Points.h:445

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

◆ v_CalculatePoints()

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

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

Definition at line 49 of file FourierPoints.cpp.

50 {
51  // Allocate the storage for points
53 
54  unsigned int npts = m_pointsKey.GetNumPoints();
55  ASSERTL0(!(npts % 2), "Fourier points need to be of even order");
56 
57  if (npts == 1)
58  {
59  m_points[0][0] = 0.0;
60  }
61  else
62  {
63  NekDouble dx = 2.0 / (NekDouble)(npts);
64  for (unsigned int i = 0; i < npts; ++i)
65  {
66  m_points[0][i] = -1.0 + i * dx;
67  }
68  }
69 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

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 71 of file FourierPoints.cpp.

72 {
73  // Allocate the storage for points
75 
76  unsigned int npts = m_pointsKey.GetNumPoints();
77  if (npts == 1)
78  {
79  m_weights[0] = 1.0; // midpoint rule
80  }
81  else
82  {
83  for (unsigned int i = 0; i < npts; ++i)
84  {
85  m_weights[i] = 1.0 / (NekDouble)npts;
86  }
87  }
88 }
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:377

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)
overrideprotectedvirtual

Delegate to function below.

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

Definition at line 151 of file FourierPoints.cpp.

153 {
154  int numpoints = 1;
155 
156  /// Delegate to function below.
157  return GetI(numpoints, x);
158 }

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

◆ v_GetI() [2/3]

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

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

Definition at line 142 of file FourierPoints.cpp.

144 {
145  ASSERTL0(pkey.GetPointsDim() == 1,
146  "Fourier Points can only interp to other 1d point distributions");
147 
148  return m_InterpManager[pkey];
149 }

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 ( unsigned int  numpoints,
const Array< OneD, const NekDouble > &  x 
)
overrideprotectedvirtual

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

Definition at line 160 of file FourierPoints.cpp.

162 {
163  Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
164 
165  CalculateInterpMatrix(numpoints, x, interp);
166 
167  NekDouble *t = interp.data();
168  unsigned int np = GetNumPoints();
169  std::shared_ptr<NekMatrix<NekDouble>> returnval(
170  MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
171  t));
172 
173  return returnval;
174 }
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
unsigned int GetNumPoints() const
Definition: Points.h:273

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:170

Definition at line 115 of file FourierPoints.h.