Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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. More...
 
 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 More...
 
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 >
 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.

54  {
55  }
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, and Nektar::LibUtilities::Points< NekDouble >::m_InterpManager.

71  :PointsBaseType(pkey)
72  {
73  m_InterpManager.RegisterCreator(PointsKey(0, eGaussGaussLegendre),
74  boost::bind(&GaussPoints::CreateMatrix, this, _1));
75  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMLegendre),
76  boost::bind(&GaussPoints::CreateMatrix, this, _1));
77  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauPLegendre),
78  boost::bind(&GaussPoints::CreateMatrix, this, _1));
79  m_InterpManager.RegisterCreator(PointsKey(0, eGaussLobattoLegendre),
80  boost::bind(&GaussPoints::CreateMatrix, this, _1));
81  m_InterpManager.RegisterCreator(PointsKey(0, eGaussGaussChebyshev),
82  boost::bind(&GaussPoints::CreateMatrix, this, _1));
83  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMChebyshev),
84  boost::bind(&GaussPoints::CreateMatrix, this, _1));
85  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauPChebyshev),
86  boost::bind(&GaussPoints::CreateMatrix, this, _1));
87  m_InterpManager.RegisterCreator(PointsKey(0, eGaussLobattoChebyshev),
88  boost::bind(&GaussPoints::CreateMatrix, this, _1));
89  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta1),
90  boost::bind(&GaussPoints::CreateMatrix, this, _1));
91  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
92  boost::bind(&GaussPoints::CreateMatrix, this, _1));
93  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
94  boost::bind(&GaussPoints::CreateMatrix, this, _1));
95  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
96  boost::bind(&GaussPoints::CreateMatrix, this, _1));
97  m_InterpManager.RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
98  boost::bind(&GaussPoints::CreateMatrix, this, _1));
99  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauKronrodMLegendre),
100  boost::bind(&GaussPoints::CreateMatrix, this, _1));
101  m_InterpManager.RegisterCreator(PointsKey(0, eGaussRadauKronrodMAlpha1Beta0),
102  boost::bind(&GaussPoints::CreateMatrix, this, _1));
103  m_InterpManager.RegisterCreator(PointsKey(0, eGaussLobattoKronrodLegendre),
104  boost::bind(&GaussPoints::CreateMatrix, this, _1));
105  m_InterpManager.RegisterCreator(PointsKey(0, eFourierEvenlySpaced),
106  boost::bind(&GaussPoints::CreateMatrix, this, _1));
107  m_InterpManager.RegisterCreator(PointsKey(0, ePolyEvenlySpaced),
108  boost::bind(&GaussPoints::CreateMatrix, this, _1));
109  m_InterpManager.RegisterCreator(PointsKey(0, eBoundaryLayerPoints),
110  boost::bind(&GaussPoints::CreateMatrix, this, _1));
111  m_InterpManager.RegisterCreator(PointsKey(0, eBoundaryLayerPointsRev),
112  boost::bind(&GaussPoints::CreateMatrix, this, _1));
113  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussGaussLegendre),
114  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
115  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMLegendre),
116  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
117  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauPLegendre),
118  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
119  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussLobattoLegendre),
120  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
121  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussGaussChebyshev),
122  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
123  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMChebyshev),
124  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
125  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauPChebyshev),
126  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
127  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussLobattoChebyshev),
128  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
129  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta1),
130  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
131  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
132  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
133  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
134  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
135  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
136  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
137  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
138  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
139  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauKronrodMLegendre),
140  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
141  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussRadauKronrodMAlpha1Beta0),
142  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
143  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eGaussLobattoKronrodLegendre),
144  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
145  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eFourierEvenlySpaced),
146  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
147  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, ePolyEvenlySpaced),
148  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
149  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eBoundaryLayerPoints),
150  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
151  m_GalerkinProjectionManager.RegisterCreator(PointsKey(0, eBoundaryLayerPointsRev),
152  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
153  }
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:375
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:62
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:54
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
1D Gauss-Kronrod-Legendre quadrature points
Definition: PointsType.h:61
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:63
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:65
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:64
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:66
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:55
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:53
1D power law distribution for boundary layer points
Definition: PointsType.h:69
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:60
1D power law distribution for boundary layer points
Definition: PointsType.h:68
boost::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
boost::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix(const PointsKey &pkey)
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:374
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:56
Points< NekDouble > PointsBaseType
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 ( )
privatevirtual

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

Definition at line 135 of file GaussPoints.cpp.

References ASSERTL0, Nektar::LibUtilities::Points< NekDouble >::CalculateDerivMatrix(), CellMLToNektar.pycml::copy(), 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.

136  {
137  // Allocate the derivative matrix
139 
140  int numpoints = m_pointsKey.GetNumPoints();
141  int totpoints = m_pointsKey.GetTotNumPoints();
142  Array<OneD, NekDouble> dmtempSharedArray = Array<OneD, NekDouble>(totpoints*totpoints);
143  NekDouble *dmtemp = dmtempSharedArray.data();
144 
145  switch(m_pointsKey.GetPointsType())
146  {
147  case eGaussGaussLegendre:
148  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
149  break;
150 
152  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
153  break;
154 
156  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
157  break;
158 
160  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
161  break;
162 
164  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
165  break;
166 
168  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
169  break;
170 
172  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
173  break;
174 
176  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
177  break;
178 
180  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,1.0);
181  break;
182 
184  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,2.0);
185  break;
186 
188  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,1.0,0.0);
189  break;
190 
192  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,2.0,0.0);
193  break;
194 
199  {
200  for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
201  {
202  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
203  {
205  }
206  }
207  return;
208  }
209  break;
210 
211  default:
212  ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
213  }
214 
215  std::copy(dmtemp,dmtemp+totpoints*totpoints,m_derivmatrix[0]->begin());
216  }
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:972
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:62
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:54
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
MatrixSharedPtrType m_derivmatrix[3]
Definition: Points.h:373
Array< OneD, DataType > m_points[3]
Definition: Points.h:371
1D Gauss-Kronrod-Legendre quadrature points
Definition: PointsType.h:61
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:63
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:64
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:1058
unsigned int GetTotNumPoints() const
Definition: Points.h:179
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:55
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:53
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1158
double NekDouble
unsigned int GetNumPoints() const
Definition: Points.h:106
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:60
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1260
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
PointsType GetPointsType() const
Definition: Points.h:111
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:56
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().

410  {
411  int numpointsfrom = pkey.GetNumPoints();
412  int numpointsto = GetNumPoints();
413 
414  Array<OneD, const NekDouble> weightsfrom;
415 
416  weightsfrom = PointsManager()[pkey]->GetW();
417 
418  boost::shared_ptr< NekMatrix<NekDouble> > Interp = GetI(pkey);
419 
420  Array<OneD, NekDouble> GalProj(numpointsfrom*numpointsto);
421 
422  // set up inner product matrix and multiply by inverse of
423  // diagaonal mass matrix
424  for(int i = 0; i < numpointsto; ++i)
425  {
426  Vmath::Vmul(numpointsfrom,Interp->GetPtr().get() +i*numpointsfrom,1,
427  &weightsfrom[0],1,&GalProj[0] +i,numpointsto);
428  Vmath::Smul(numpointsfrom,1.0/m_weights[i],&GalProj[0]+i,numpointsto,
429  &GalProj[0]+i,numpointsto);
430  }
431 
432 
433  NekDouble* t = GalProj.data();
434  boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpointsto,numpointsfrom,t));
435 
436  return returnval;
437  }
const boost::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
unsigned int GetNumPoints() const
Definition: Points.h:268
PointsManagerT & PointsManager(void)
Array< OneD, DataType > m_weights
Definition: Points.h:372
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
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().

219  {
220  switch(m_pointsKey.GetPointsType())
221  {
222  case eGaussGaussLegendre:
223  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
224  break;
225 
227  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
228  break;
229 
231  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
232  break;
233 
235  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
236  break;
237 
239  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
240  break;
241 
243  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
244  break;
245 
247  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
248  break;
249 
251  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
252  break;
253 
255  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,1.0);
256  break;
257 
259  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,2.0);
260  break;
261 
263  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,1.0,0.0);
264  break;
265 
267  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,2.0,0.0);
268  break;
269 
274  {
275  for(unsigned int i=0;i<npts;++i)
276  {
277  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
278  {
279  interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
280  }
281  }
282  }
283  break;
284 
285  default:
286  ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
287  }
288  }
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:62
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:54
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1744
Array< OneD, DataType > m_points[3]
Definition: Points.h:371
1D Gauss-Kronrod-Legendre quadrature points
Definition: PointsType.h:61
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:63
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:64
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points.
Definition: Polylib.cpp:1806
unsigned int GetNumPoints() const
Definition: Points.h:268
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1684
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:55
static std::string npts
Definition: InputFld.cpp:43
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:53
unsigned int GetNumPoints() const
Definition: Points.h:106
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:60
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an.
Definition: Polylib.cpp:1626
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
PointsType GetPointsType() const
Definition: Points.h:111
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:56
void Nektar::LibUtilities::GaussPoints::CalculatePoints ( )
privatevirtual

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

Definition at line 51 of file GaussPoints.cpp.

References ASSERTL0, Nektar::LibUtilities::Points< NekDouble >::CalculatePoints(), Nektar::LibUtilities::Points< NekDouble >::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().

52  {
53  // Allocate the storage for points and weights
56 
57  int numpoints = m_pointsKey.GetNumPoints();
58 
59  switch(m_pointsKey.GetPointsType())
60  {
62  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
63  break;
64 
66  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
67  break;
68 
70  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
71  break;
72 
74  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
75  break;
76 
78  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
79  break;
80 
82  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
83  break;
84 
86  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
87  break;
88 
90  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
91  break;
92 
94  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,1.0);
95  break;
96 
98  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,2.0);
99  break;
100 
102  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
103  break;
104 
106  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,2.0,0.0);
107  break;
108 
110  Polylib::zwgk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
111  break;
112 
114  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
115  break;
116 
118  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
119  break;
120 
122  Polylib::zwlk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
123  break;
124 
125  default:
126  ASSERTL0(false, "Unknown Gauss quadrature point distribution requested");
127  }
128  }
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:62
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:54
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition: Polylib.cpp:274
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition: Polylib.cpp:346
Array< OneD, DataType > m_points[3]
Definition: Points.h:371
1D Gauss-Kronrod-Legendre quadrature points
Definition: PointsType.h:61
void zwlk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:746
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:63
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:64
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:420
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:55
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:53
Array< OneD, DataType > m_weights
Definition: Points.h:372
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition: Polylib.cpp:198
unsigned int GetNumPoints() const
Definition: Points.h:106
void zwrk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Radau-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:554
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:60
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
PointsType GetPointsType() const
Definition: Points.h:111
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:56
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:142
void Nektar::LibUtilities::GaussPoints::CalculateWeights ( )
privatevirtual

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

Definition at line 130 of file GaussPoints.cpp.

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

Definition at line 290 of file GaussPoints.cpp.

291  {
292  boost::shared_ptr< Points<NekDouble> > returnval(MemoryManager< GaussPoints >::AllocateSharedPtr(pkey));
293 
294  returnval->Initialize();
295 
296  return returnval;
297  }
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().

402  {
403  boost::shared_ptr< NekMatrix<NekDouble> > returnval = CalculateGalerkinProjectionMatrix(pkey);
404 
405  // Delegate to function below
406  return returnval;
407  }
boost::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix(const PointsKey &pkey)
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().

300  {
301  int numpoints = pkey.GetNumPoints();
302  Array<OneD, const NekDouble> xpoints;
303 
304  PointsManager()[pkey]->GetPoints(xpoints);
305 
306  // Delegate to function below
307  return GetI(numpoints, xpoints);
308  }
const boost::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
PointsManagerT & PointsManager(void)
const boost::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::GaussPoints::GetGalerkinProjection ( const PointsKey pkey)
virtual

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

Definition at line 396 of file GaussPoints.cpp.

References Nektar::LibUtilities::Points< NekDouble >::m_GalerkinProjectionManager.

397  {
398  return m_GalerkinProjectionManager[pkey];
399  }
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:375
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().

312  {
313  ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
314 
315  return m_InterpManager[pkey];
316  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:374
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().

319  {
320  int numpoints = 1;
321 
322  // Delegate to function below
323  return GetI(numpoints, x);
324  }
const boost::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
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().

327  {
328  Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
329 
330  CalculateInterpMatrix(numpoints, x, interp);
331 
332  NekDouble* t = interp.data();
333  unsigned int np = GetNumPoints();
334  boost::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np,t));
335 
336  return returnval;
337  }
unsigned int GetNumPoints() const
Definition: Points.h:268
double NekDouble
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
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.

341  {
342  NekDouble sum = 0.0;
343 
344  for(int i=0;i<npts;++i)
345  {
346  sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
347  }
348  return sum;
349  }
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
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().

353  {
354  NekDouble h=1.0;
355 
356  for(int i=0;i<pt; ++i)
357  {
358  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
359  }
360 
361  for(int i=pt+1;i<npts;++i)
362  {
363  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
364  }
365 
366  return h;
367  }
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
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().

370  {
371  NekDouble h;
372  NekDouble y=0.0;
373 
374  for(int j=0;j<npts;++j)
375  {
376  if(j!=pt)
377  {
378  h=1.0;
379  for(int i=0;i<npts;++i)
380  {
381  if(i!=pt)
382  {
383  if(i!=j)
384  {
385  h = h*(x-xpts[i]);
386  }
387  h = h/(xpts[pt]-xpts[i]);
388  }
389  }
390  y = y + h;
391  }
392  }
393  return y;
394  }
static std::string npts
Definition: InputFld.cpp:43
double NekDouble