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.
 PolyEPoints (const PolyEPoints &points)
 Copy constructor should not be called.
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 >
virtual void CalculatePoints ()
virtual void CalculateWeights ()
virtual void CalculateDerivMatrix ()
 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.

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

Definition at line 70 of file PolyEPoints.h.

{
}
Nektar::LibUtilities::PolyEPoints::PolyEPoints ( )
inlineprivate

Default constructor should not be called except by Create method.

Definition at line 76 of file PolyEPoints.h.

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

Copy constructor should not be called.

Definition at line 81 of file PolyEPoints.h.

:PointsBaseType(points.m_pointsKey)
{
}

Member Function Documentation

void Nektar::LibUtilities::PolyEPoints::CalculateDerivMatrix ( )
private
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().

{
for(unsigned int i=0;i<npts;++i)
{
for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
{
interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
}
}
}
void Nektar::LibUtilities::PolyEPoints::CalculatePoints ( )
private

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.

{
// Allocate the storage for points
unsigned int npts = m_pointsKey.GetNumPoints();
if(npts==1)
{
m_points[0][0] = 0.0;
}
else
{
NekDouble dx = 2.0/(NekDouble)(npts-1);
for(unsigned int i=0;i<npts;++i)
{
m_points[0][i] = -1.0 + i*dx;
}
}
}
void Nektar::LibUtilities::PolyEPoints::CalculateWeights ( )
private

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

{
// Allocate the storage for points
unsigned int npts = m_pointsKey.GetNumPoints();
if(npts==1)
{
m_weights[0] = 2.0; //midpoint rule
}
else
{
PointsKey gaussKey(npts, eGaussLobattoLegendre);
boost::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
Array<OneD, const NekDouble> z;
Array<OneD, const NekDouble> w;
ptr->GetZW(z,w);
for(unsigned int i=0; i<npts;++i)
{
m_weights[i] = 0.0;
for(unsigned j=0;j<npts;++j)
{
m_weights[i] += w[j]*LagrangePoly(z[j],i,npts,m_points[0]);
}
}
}
}
boost::shared_ptr< PolyEPoints::PointsBaseType > Nektar::LibUtilities::PolyEPoints::Create ( const PointsKey key)
static

Definition at line 120 of file PolyEPoints.cpp.

{
boost::shared_ptr<PointsBaseType> returnval(MemoryManager<PolyEPoints>::AllocateSharedPtr(key));
returnval->Initialize();
return returnval;
}
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().

{
ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
int numpoints = pkey.GetNumPoints();
Array<OneD, const NekDouble> xpoints;
PointsManager()[pkey]->GetPoints(xpoints);
return GetI(numpoints, xpoints);
}
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().

{
int numpoints = 1;
/// Delegate to function below.
return GetI(numpoints, x);
}
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().

{
Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
CalculateInterpMatrix(numpoints, x, interp);
unsigned int np = GetTotNumPoints();
NekDouble* d = interp.data();
boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np, d));
return returnval;
}
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.

{
NekDouble sum = 0.0;
for(int i=0;i<npts;++i)
{
sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
}
return sum;
}
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().

{
NekDouble h=1.0;
for(int i=0;i<pt; ++i)
{
h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
}
for(int i=pt+1;i<npts;++i)
{
h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
}
return h;
}
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().

{
NekDouble y=0.0;
for(int j=0;j<npts;++j)
{
if(j!=pt)
{
h=1.0;
for(int i=0;i<npts;++i)
{
if(i!=pt)
{
if(i!=j)
{
h = h*(x-xpts[i]);
}
h = h/(xpts[pt]-xpts[i]);
}
}
y = y + h;
}
}
return y;
}