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

#include <GaussPoints.h>

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

Public Member Functions

virtual ~GaussPoints ()
boost::shared_ptr< NekMatrix
< NekDouble > > 
CreateMatrix (const PointsKey &pkey)
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)
boost::shared_ptr< NekMatrix
< NekDouble > > 
CreateGPMatrix (const PointsKey &pkey)
const boost::shared_ptr
< NekMatrix< NekDouble > > 
GetGalerkinProjection (const PointsKey &pkey)
 GaussPoints (const PointsKey &pkey)
- 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)

Static Public Member Functions

static boost::shared_ptr
< Points< NekDouble > > 
Create (const PointsKey &pkey)

Private Member Functions

 GaussPoints ()
 These should not be called. All creation is done using the constructor requiring the key, declared above.
 GaussPoints (const GaussPoints &points)
void CalculatePoints ()
void CalculateWeights ()
void CalculateDerivMatrix ()
void CalculateInterpMatrix (unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
boost::shared_ptr< NekMatrix
< NekDouble > > 
CalculateGalerkinProjectionMatrix (const PointsKey &pkey)
NekDouble LagrangeInterpolant (NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
 functions used by the Kronrod points
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

- Public Types inherited from Nektar::LibUtilities::Points< NekDouble >
typedef NekDouble DataType
typedef boost::shared_ptr
< NekMatrix< DataType > > 
MatrixSharedPtrType
- 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 50 of file GaussPoints.h.

Constructor & Destructor Documentation

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

Definition at line 53 of file GaussPoints.h.

{
}
Nektar::LibUtilities::GaussPoints::GaussPoints ( const PointsKey pkey)
inline

Definition at line 71 of file GaussPoints.h.

References CreateGPMatrix(), CreateMatrix(), Nektar::LibUtilities::eBoundaryLayerPoints, Nektar::LibUtilities::eBoundaryLayerPointsRev, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::LibUtilities::Points< NekDouble >::m_GalerkinProjectionManager, Nektar::LibUtilities::Points< NekDouble >::m_InterpManager, and Nektar::LibUtilities::NekManager< KeyType, ValueT, opLessCreator >::RegisterCreator().

{
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
}
Nektar::LibUtilities::GaussPoints::GaussPoints ( )
private

These should not be called. All creation is done using the constructor requiring the key, declared above.

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

Member Function Documentation

void Nektar::LibUtilities::GaussPoints::CalculateDerivMatrix ( )
private

Definition at line 135 of file GaussPoints.cpp.

References ASSERTL0, Polylib::Dgj(), Polylib::Dglj(), Polylib::Dgrjm(), Polylib::Dgrjp(), Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Nektar::LibUtilities::PointsKey::GetTotNumPoints(), LagrangePolyDeriv(), Nektar::LibUtilities::Points< NekDouble >::m_derivmatrix, Nektar::LibUtilities::Points< NekDouble >::m_points, and Nektar::LibUtilities::Points< NekDouble >::m_pointsKey.

{
// Allocate the derivative matrix
int numpoints = m_pointsKey.GetNumPoints();
int totpoints = m_pointsKey.GetTotNumPoints();
Array<OneD, NekDouble> dmtempSharedArray = Array<OneD, NekDouble>(totpoints*totpoints);
NekDouble *dmtemp = dmtempSharedArray.data();
{
Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
break;
Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
break;
Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
break;
Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
break;
Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
break;
Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,1.0);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,2.0);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,1.0,0.0);
break;
Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,2.0,0.0);
break;
{
for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
{
for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
{
}
}
return;
}
break;
default:
ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
}
std::copy(dmtemp,dmtemp+totpoints*totpoints,m_derivmatrix[0]->begin());
}
boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CalculateGalerkinProjectionMatrix ( const PointsKey pkey)
private

Definition at line 409 of file GaussPoints.cpp.

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

Referenced by CreateGPMatrix().

{
int numpointsfrom = pkey.GetNumPoints();
int numpointsto = GetNumPoints();
Array<OneD, const NekDouble> weightsfrom;
weightsfrom = PointsManager()[pkey]->GetW();
boost::shared_ptr< NekMatrix<NekDouble> > Interp = GetI(pkey);
Array<OneD, NekDouble> GalProj(numpointsfrom*numpointsto);
// set up inner product matrix and multiply by inverse of
// diagaonal mass matrix
for(int i = 0; i < numpointsto; ++i)
{
Vmath::Vmul(numpointsfrom,Interp->GetPtr().get() +i*numpointsfrom,1,
&weightsfrom[0],1,&GalProj[0] +i,numpointsto);
Vmath::Smul(numpointsfrom,1.0/m_weights[i],&GalProj[0]+i,numpointsto,
&GalProj[0]+i,numpointsto);
}
NekDouble* t = GalProj.data();
boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpointsto,numpointsfrom,t));
return returnval;
}
void Nektar::LibUtilities::GaussPoints::CalculateInterpMatrix ( unsigned int  npts,
const Array< OneD, const NekDouble > &  xpoints,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 218 of file GaussPoints.cpp.

References ASSERTL0, Nektar::Array< OneD, DataType >::data(), Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Polylib::Imgj(), Polylib::Imglj(), Polylib::Imgrjm(), Polylib::Imgrjp(), LagrangePoly(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, and npts.

Referenced by GetI().

{
{
Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
break;
Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
break;
Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
break;
Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
break;
Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
break;
Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,1.0);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,2.0);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,1.0,0.0);
break;
Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,2.0,0.0);
break;
{
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]);
}
}
}
break;
default:
ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
}
}
void Nektar::LibUtilities::GaussPoints::CalculatePoints ( )
private

Definition at line 51 of file GaussPoints.cpp.

References ASSERTL0, CalculateWeights(), Nektar::Array< OneD, DataType >::data(), Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussKronrodLegendre, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauKronrodMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauKronrodMLegendre, Nektar::LibUtilities::eGaussRadauMAlpha0Beta1, Nektar::LibUtilities::eGaussRadauMAlpha0Beta2, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGaussRadauMChebyshev, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eGaussRadauPChebyshev, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::LibUtilities::PointsKey::GetPointsType(), Nektar::LibUtilities::Points< NekDouble >::m_points, Nektar::LibUtilities::Points< NekDouble >::m_pointsKey, Nektar::LibUtilities::Points< NekDouble >::m_weights, Polylib::zwgj(), Polylib::zwgk(), Polylib::zwglj(), Polylib::zwgrjm(), Polylib::zwgrjp(), Polylib::zwlk(), and Polylib::zwrk().

{
// Allocate the storage for points and weights
int numpoints = m_pointsKey.GetNumPoints();
{
Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
break;
Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
break;
Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,1.0);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,2.0);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
break;
Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,2.0,0.0);
break;
Polylib::zwgk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
break;
Polylib::zwlk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
break;
default:
ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
}
}
void Nektar::LibUtilities::GaussPoints::CalculateWeights ( )
private

Definition at line 130 of file GaussPoints.cpp.

Referenced by CalculatePoints().

{
// For Gauss Quadrature, this is done as part of the points computation
}
boost::shared_ptr< Points< NekDouble > > Nektar::LibUtilities::GaussPoints::Create ( const PointsKey pkey)
static

Definition at line 290 of file GaussPoints.cpp.

{
boost::shared_ptr< Points<NekDouble> > returnval(MemoryManager< GaussPoints >::AllocateSharedPtr(pkey));
returnval->Initialize();
return returnval;
}
boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CreateGPMatrix ( const PointsKey pkey)

Definition at line 401 of file GaussPoints.cpp.

References CalculateGalerkinProjectionMatrix().

Referenced by GaussPoints().

{
boost::shared_ptr< NekMatrix<NekDouble> > returnval = CalculateGalerkinProjectionMatrix(pkey);
// Delegate to function below
return returnval;
}
boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::CreateMatrix ( const PointsKey pkey)

Definition at line 299 of file GaussPoints.cpp.

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

Referenced by GaussPoints().

{
int numpoints = pkey.GetNumPoints();
Array<OneD, const NekDouble> xpoints;
PointsManager()[pkey]->GetPoints(xpoints);
// Delegate to function below
return GetI(numpoints, xpoints);
}
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetGalerkinProjection ( const PointsKey pkey)
virtual
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetI ( const PointsKey pkey)
virtual

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

Definition at line 311 of file GaussPoints.cpp.

References ASSERTL0, Nektar::LibUtilities::PointsKey::GetPointsDim(), and Nektar::LibUtilities::Points< NekDouble >::m_InterpManager.

Referenced by CalculateGalerkinProjectionMatrix(), CreateMatrix(), and GetI().

{
ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
return m_InterpManager[pkey];
}
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetI ( const Array< OneD, const NekDouble > &  x)
virtual

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

Definition at line 318 of file GaussPoints.cpp.

References GetI().

{
int numpoints = 1;
// Delegate to function below
return GetI(numpoints, x);
}
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetI ( unsigned int  numpoints,
const Array< OneD, const NekDouble > &  x 
)
virtual

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

Definition at line 326 of file GaussPoints.cpp.

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

{
Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
CalculateInterpMatrix(numpoints, x, interp);
NekDouble* t = interp.data();
unsigned int np = GetNumPoints();
boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np,t));
return returnval;
}
NekDouble Nektar::LibUtilities::GaussPoints::LagrangeInterpolant ( NekDouble  x,
int  npts,
const Array< OneD, const NekDouble > &  xpts,
const Array< OneD, const NekDouble > &  funcvals 
)
private

functions used by the Kronrod points

Definition at line 339 of file GaussPoints.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::GaussPoints::LagrangePoly ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 352 of file GaussPoints.cpp.

References npts.

Referenced by CalculateInterpMatrix(), 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::GaussPoints::LagrangePolyDeriv ( NekDouble  x,
int  pt,
int  npts,
const Array< OneD, const NekDouble > &  xpts 
)
private

Definition at line 369 of file GaussPoints.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;
}