Sundance::QuadratureFamilyBase Class Reference

Inheritance diagram for Sundance::QuadratureFamilyBase:

Sundance::QuadratureFamilyStub Playa::Handleable< QuadratureFamilyStub > Playa::Printable Sundance::Noncopyable Sundance::ClosedNewtonCotes Sundance::FeketeQuadrature Sundance::GaussianQuadrature Sundance::GaussLobattoQuadrature Sundance::PolygonQuadrature Sundance::ReducedQuadrature Sundance::SurfQuadrature Sundance::TrapesoidQuadrature

List of all members.

Public Member Functions

 QuadratureFamilyBase (int order)
virtual ~QuadratureFamilyBase ()
virtual int getNumPoints (const CellType &cellType) const
virtual void getPoints (const CellType &cellType, Array< Point > &quadPoints, Array< double > &quadWeights) const
virtual void getAdaptedWeights (const CellType &cellType, int cellDim, int celLID, int facetIndex, const Mesh &mesh, const ParametrizedCurve &globalCurve, Array< Point > &quadPoints, Array< double > &quadWeights, bool &isCut) const

Protected Member Functions

virtual void getLineRule (Array< Point > &quadPoints, Array< double > &quadWeights) const
virtual void getTriangleRule (Array< Point > &quadPoints, Array< double > &quadWeights) const
virtual void getQuadRule (Array< Point > &quadPoints, Array< double > &quadWeights) const
virtual void getTetRule (Array< Point > &quadPoints, Array< double > &quadWeights) const
virtual void getBrickRule (Array< Point > &quadPoints, Array< double > &quadWeights) const


Detailed Description

QuadratureFamilyBase extends QuadratureFamilyStub to provide an interface for getting quadrature points for a given cell type.

Definition at line 51 of file SundanceQuadratureFamilyBase.hpp.


Constructor & Destructor Documentation

Sundance::QuadratureFamilyBase::QuadratureFamilyBase ( int  order  )  [inline]

Definition at line 55 of file SundanceQuadratureFamilyBase.hpp.

virtual Sundance::QuadratureFamilyBase::~QuadratureFamilyBase (  )  [inline, virtual]

Definition at line 58 of file SundanceQuadratureFamilyBase.hpp.


Member Function Documentation

void QuadratureFamilyBase::getAdaptedWeights ( const CellType cellType,
int  cellDim,
int  celLID,
int  facetIndex,
const Mesh mesh,
const ParametrizedCurve &  globalCurve,
Array< Point > &  quadPoints,
Array< double > &  quadWeights,
bool &  isCut 
) const [virtual]

This method is used for the Adaptive Cell Integration, which returns special quadrature weights for cells which are cut by the curve
The quadPoints and the quadWeights are set to the default values. (quadrature without the curve , no curve) quadWeights will be changed depending on the curve, if the curve cuts the cell isCut should be set by this method, if it is true then the quadWeights will be used for quadrature of the cell

Reimplemented in Sundance::FeketeQuadrature, Sundance::GaussianQuadrature, Sundance::GaussLobattoQuadrature, Sundance::ReducedQuadrature, and Sundance::TrapesoidQuadrature.

Definition at line 73 of file SundanceQuadratureFamilyBase.cpp.

References SUNDANCE_ERROR, and Sundance::QuadratureFamilyStub::toXML().

Referenced by Sundance::QuadratureFamily::getAdaptedWeights().

void QuadratureFamilyBase::getBrickRule ( Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [protected, virtual]

void QuadratureFamilyBase::getLineRule ( Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [protected, virtual]

virtual int Sundance::QuadratureFamilyBase::getNumPoints ( const CellType cellType  )  const [inline, virtual]

Gets number of points associated with a particular cell type: WARNING: this is slow. Call it once and store the result. TODO: make it pure virtual and override with queries in the derived classes, making them supply the information.

Reimplemented in Sundance::ReducedQuadrature.

Definition at line 64 of file SundanceQuadratureFamilyBase.hpp.

References getPoints().

Referenced by Sundance::QuadratureFamily::getNumPoints().

void QuadratureFamilyBase::getPoints ( const CellType cellType,
Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [virtual]

void QuadratureFamilyBase::getQuadRule ( Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [protected, virtual]

void QuadratureFamilyBase::getTetRule ( Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [protected, virtual]

compute a rule for the reference tet cell

Reimplemented in Sundance::GaussianQuadrature, Sundance::GaussLobattoQuadrature, and Sundance::TrapesoidQuadrature.

Definition at line 119 of file SundanceQuadratureFamilyBase.cpp.

References SUNDANCE_ERROR, and Sundance::QuadratureFamilyStub::toXML().

Referenced by getPoints().

void QuadratureFamilyBase::getTriangleRule ( Array< Point > &  quadPoints,
Array< double > &  quadWeights 
) const [protected, virtual]

Site Contact