49#ifndef INTREPID_CUBATURE_CONTROLVOLUME_HPP
50#define INTREPID_CUBATURE_CONTROLVOLUME_HPP
53#include "Teuchos_Assert.hpp"
54#include "Shards_CellTopology.hpp"
66 template<
class Scalar,
class ArrayPo
int,
class ArrayWeight>
83 ArrayWeight& cubWeights)
const;
93 ArrayWeight& cubWeights,
94 ArrayPoint& cellCoords)
const;
107 void getAccuracy(std::vector<int> & accuracy)
const;
139#include "Intrepid_CubatureControlVolumeDef.hpp"
Header file for the Intrepid::Cubature class.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Defines cubature (integration) rules over control volumes.
int getNumPoints() const
Returns the number of cubature points.
int getDimension() const
Returns dimension of integration domain.
void getAccuracy(std::vector< int > &accuracy) const
Returns max. degree of polynomials that are integrated exactly. The return vector has size 1.
void getCubature(ArrayPoint &cubPoints, ArrayWeight &cubWeights) const
Returns cubature points and weights Method for reference space cubature - throws an exception.
Teuchos::RCP< const shards::CellTopology > subCVCellTopo_
The topology of the sub-control volume.
int degree_
The degree of the polynomials that are integrated exactly.
Teuchos::RCP< const shards::CellTopology > primaryCellTopo_
The topology of the primary cell.
int cubDimension_
Dimension of integration domain.
int numPoints_
The number of cubature points.
Defines the base class for cubature (integration) rules in Intrepid.