Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | List of all members
Nektar::LibUtilities::PolyEPoints Class Reference

#include <PolyEPoints.h>

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

Public Types

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

Public Member Functions

virtual ~PolyEPoints ()
 
const boost::shared_ptr
< NekMatrix< NekDouble > > 
GetI (const PointsKey &pkey)
 
const boost::shared_ptr
< NekMatrix< NekDouble > > 
GetI (const Array< OneD, const NekDouble > &x)
 
const boost::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
 
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 boost::shared_ptr
< PointsBaseType
Create (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)
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
 Points (const PointsKey &key)
 
- Protected Attributes inherited from Nektar::LibUtilities::Points< NekDouble >
PointsKey m_pointsKey
 
Array< OneD, DataTypem_points [3]
 
Array< OneD, DataTypem_weights
 
MatrixSharedPtrType m_derivmatrix [3]
 
NekManager< PointsKey,
NekMatrix< DataType >
, PointsKey::opLess
m_InterpManager
 
NekManager< PointsKey,
NekMatrix< DataType >
, PointsKey::opLess
m_GalerkinProjectionManager
 

Detailed Description

Definition at line 55 of file PolyEPoints.h.

Member Typedef Documentation

Definition at line 58 of file PolyEPoints.h.

Constructor & Destructor Documentation

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

Definition at line 60 of file PolyEPoints.h.

61  {
62  }
Nektar::LibUtilities::PolyEPoints::PolyEPoints ( const PointsKey key)
inline

Definition at line 70 of file PolyEPoints.h.

70  :PointsBaseType(key)
71  {
72  }
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:58
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)
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:58
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  }
Points< NekDouble > PointsBaseType
Definition: PolyEPoints.h:58

Member Function Documentation

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

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

Definition at line 95 of file PolyEPoints.cpp.

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.

96  {
97  // Allocate the derivative matrix.
99 
100  for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
101  {
102  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
103  {
105  }
106  }
107  }
MatrixSharedPtrType m_derivmatrix[3]
Definition: Points.h:351
Array< OneD, DataType > m_points[3]
Definition: Points.h:349
unsigned int GetNumPoints() const
Definition: Points.h:106
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
void Nektar::LibUtilities::PolyEPoints::CalculateInterpMatrix ( unsigned int  npts,
const Array< OneD, const NekDouble > &  xpoints,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 109 of file PolyEPoints.cpp.

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

Referenced by GetI().

110  {
111  for(unsigned int i=0;i<npts;++i)
112  {
113  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
114  {
115  interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
116  }
117  }
118  }
Array< OneD, DataType > m_points[3]
Definition: Points.h:349
static std::string npts
Definition: InputFld.cpp:43
unsigned int GetNumPoints() const
Definition: Points.h:106
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
void Nektar::LibUtilities::PolyEPoints::CalculatePoints ( )
privatevirtual

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

Definition at line 46 of file PolyEPoints.cpp.

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

47  {
48  // Allocate the storage for points
50 
51  unsigned int npts = m_pointsKey.GetNumPoints();
52  if(npts==1)
53  {
54  m_points[0][0] = 0.0;
55  }
56  else
57  {
58  NekDouble dx = 2.0/(NekDouble)(npts-1);
59  for(unsigned int i=0;i<npts;++i)
60  {
61  m_points[0][i] = -1.0 + i*dx;
62  }
63  }
64  }
Array< OneD, DataType > m_points[3]
Definition: Points.h:349
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
unsigned int GetNumPoints() const
Definition: Points.h:106
void Nektar::LibUtilities::PolyEPoints::CalculateWeights ( )
privatevirtual

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

Definition at line 66 of file PolyEPoints.cpp.

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, npts, and Nektar::LibUtilities::PointsManager().

67  {
68  // Allocate the storage for points
70 
71  unsigned int npts = m_pointsKey.GetNumPoints();
72  if(npts==1)
73  {
74  m_weights[0] = 2.0; //midpoint rule
75  }
76  else
77  {
78  PointsKey gaussKey(npts, eGaussLobattoLegendre);
79  boost::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
80  Array<OneD, const NekDouble> z;
81  Array<OneD, const NekDouble> w;
82 
83  ptr->GetZW(z,w);
84  for(unsigned int i=0; i<npts;++i)
85  {
86  m_weights[i] = 0.0;
87  for(unsigned j=0;j<npts;++j)
88  {
89  m_weights[i] += w[j]*LagrangePoly(z[j],i,npts,m_points[0]);
90  }
91  }
92  }
93  }
Array< OneD, DataType > m_points[3]
Definition: Points.h:349
static std::string npts
Definition: InputFld.cpp:43
PointsManagerT & PointsManager(void)
Array< OneD, DataType > m_weights
Definition: Points.h:350
unsigned int GetNumPoints() const
Definition: Points.h:106
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
boost::shared_ptr< PolyEPoints::PointsBaseType > Nektar::LibUtilities::PolyEPoints::Create ( const PointsKey key)
static

Definition at line 120 of file PolyEPoints.cpp.

121  {
122  boost::shared_ptr<PointsBaseType> returnval(MemoryManager<PolyEPoints>::AllocateSharedPtr(key));
123 
124  returnval->Initialize();
125 
126  return returnval;
127  }
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::PolyEPoints::GetI ( const PointsKey pkey)
virtual

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

Definition at line 129 of file PolyEPoints.cpp.

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

Referenced by GetI().

130  {
131  ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
132 
133  int numpoints = pkey.GetNumPoints();
134  Array<OneD, const NekDouble> xpoints;
135 
136  PointsManager()[pkey]->GetPoints(xpoints);
137 
138  return GetI(numpoints, xpoints);
139  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
const boost::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
PointsManagerT & PointsManager(void)
const boost::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 141 of file PolyEPoints.cpp.

References GetI().

142  {
143  int numpoints = 1;
144 
145  /// Delegate to function below.
146  return GetI(numpoints, x);
147  }
const boost::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
const boost::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 149 of file PolyEPoints.cpp.

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

150  {
151  Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
152 
153  CalculateInterpMatrix(numpoints, x, interp);
154 
155  unsigned int np = GetTotNumPoints();
156  NekDouble* d = interp.data();
157  boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np, d));
158 
159  return returnval;
160  }
unsigned int GetNumPoints() const
Definition: Points.h:246
double NekDouble
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
unsigned int GetTotNumPoints() const
Definition: Points.h:251
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 163 of file PolyEPoints.cpp.

References LagrangePoly(), and npts.

165  {
166  NekDouble sum = 0.0;
167 
168  for(int i=0;i<npts;++i)
169  {
170  sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
171  }
172  return sum;
173  }
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
NekDouble Nektar::LibUtilities::PolyEPoints::LagrangePoly ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 176 of file PolyEPoints.cpp.

References npts.

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

177  {
178  NekDouble h=1.0;
179 
180  for(int i=0;i<pt; ++i)
181  {
182  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
183  }
184 
185  for(int i=pt+1;i<npts;++i)
186  {
187  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
188  }
189 
190  return h;
191  }
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
NekDouble Nektar::LibUtilities::PolyEPoints::LagrangePolyDeriv ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 193 of file PolyEPoints.cpp.

References npts.

Referenced by CalculateDerivMatrix().

194  {
195  NekDouble h;
196  NekDouble y=0.0;
197 
198  for(int j=0;j<npts;++j)
199  {
200  if(j!=pt)
201  {
202  h=1.0;
203  for(int i=0;i<npts;++i)
204  {
205  if(i!=pt)
206  {
207  if(i!=j)
208  {
209  h = h*(x-xpts[i]);
210  }
211  h = h/(xpts[pt]-xpts[i]);
212  }
213  }
214  y = y + h;
215  }
216  }
217  return y;
218  }
static std::string npts
Definition: InputFld.cpp:43
double NekDouble