49 #ifndef INTREPID_CUBATURE_CONTROLVOLUMEDEF_HPP 
   50 #define INTREPID_CUBATURE_CONTROLVOLUMEDEF_HPP 
   54 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
   58   primaryCellTopo_ = cellTopology;
 
   61   const CellTopologyData &myCellData =
 
   62          (primaryCellTopo_->getDimension() > 2) ? *shards::getCellTopologyData<shards::Hexahedron<8> >() :
 
   63                                                   *shards::getCellTopologyData<shards::Quadrilateral<4> >();
 
   64   subCVCellTopo_ = Teuchos::rcp(
new shards::CellTopology(&myCellData));
 
   69   numPoints_ = primaryCellTopo_->getNodeCount();
 
   71   cubDimension_ = primaryCellTopo_->getDimension();
 
   75 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
   77                                                                        ArrayWeight& cubWeights)
 const 
   79     TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
 
   80                       ">>> ERROR (CubatureControlVolume): Cubature defined in physical space calling method for reference space cubature.");
 
   83 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
   85                                                                        ArrayWeight& cubWeights,
 
   86                                                                        ArrayPoint& cellCoords)
 const 
   89   int numCells         = cellCoords.dimension(0);
 
   90   int numNodesPerCell  = cellCoords.dimension(1);
 
   91   int spaceDim         = cellCoords.dimension(2);
 
   92   int numNodesPerSubCV = subCVCellTopo_->getNodeCount();
 
  100   int subcvCubDegree = 2;
 
  101   Teuchos::RCP<Intrepid::Cubature<double,Intrepid::FieldContainer<double>  > > subCVCubature;
 
  102   subCVCubature = subCVCubFactory.
create(*(subCVCellTopo_), subcvCubDegree);
 
  104   int subcvCubDim       = subCVCubature -> getDimension();
 
  105   int numSubcvCubPoints = subCVCubature -> getNumPoints();
 
  111   subCVCubature -> getCubature(subcvCubPoints, subcvCubWeights);
 
  114   for (
int icell = 0; icell < numCells; icell++){
 
  119      for (
int isubcv = 0; isubcv < numNodesPerCell; isubcv++){
 
  120        for (
int idim = 0; idim < spaceDim; idim++){
 
  121           for (
int inode = 0; inode < numNodesPerSubCV; inode++){
 
  122               subCVCenter(isubcv,0,idim) += subCVCoords(icell,isubcv,inode,idim)/numNodesPerSubCV;
 
  123               cellCVCoords(isubcv,inode,idim) = subCVCoords(icell,isubcv,inode,idim);
 
  125           cubPoints(icell,isubcv,idim) = subCVCenter(isubcv,0,idim);
 
  136      for (
int inode = 0; inode < numNodesPerCell; inode++){
 
  138          for (
int ipt = 0; ipt < numSubcvCubPoints; ipt++){
 
  139             vol += subcvCubWeights(ipt)*subCVJacobDet(inode,ipt);
 
  141          cubWeights(icell,inode) = vol;
 
  149 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
  154 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
  156   return cubDimension_;
 
  159 template<
class Scalar, 
class ArrayPo
int, 
class ArrayWeight>
 
  161   accuracy.assign(1,degree_);
 
void getCubature(ArrayPoint &cubPoints, ArrayWeight &cubWeights) const 
Returns cubature points and weights Method for reference space cubature - throws an exception...
CubatureControlVolume(const Teuchos::RCP< const shards::CellTopology > &cellTopology)
int getDimension() const 
Returns dimension of integration domain. 
int getNumPoints() const 
Returns the number of cubature points. 
void getAccuracy(std::vector< int > &accuracy) const 
Returns max. degree of polynomials that are integrated exactly. The return vector has size 1...
A factory class that generates specific instances of cubatures. 
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > °ree)
Factory method.