Intrepid2
Classes | Public Types | Public Member Functions | Protected Attributes | List of all members
Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType > Class Template Reference

Defines cubature (integration) rules over Neumann boundaries for control volume method. More...

#include <Intrepid2_CubatureControlVolumeBoundary.hpp>

Inheritance diagram for Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >:
Intrepid2::Cubature< DeviceType, pointValueType, weightValueType >

Classes

struct  Functor
 

Public Types

typedef Cubature< DeviceType,
pointValueType,
weightValueType >
::PointViewType 
PointViewType
 
typedef Cubature< DeviceType,
pointValueType,
weightValueType >
::weightViewType 
weightViewType
 
- Public Types inherited from Intrepid2::Cubature< DeviceType, pointValueType, weightValueType >
using ExecSpaceType = typename DeviceType::execution_space
 
using PointViewType = Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, DeviceType >
 
using weightViewType = Kokkos::DynRankView< weightValueType, Kokkos::LayoutStride, DeviceType >
 
using PointViewTypeAllocatable = Kokkos::DynRankView< pointValueType, DeviceType >
 
using WeightViewTypeAllocatable = Kokkos::DynRankView< weightValueType, DeviceType >
 
using TensorPointDataType = TensorPoints< pointValueType, DeviceType >
 
using TensorWeightDataType = TensorData< weightValueType, DeviceType >
 

Public Member Functions

virtual void getCubature (PointViewType cubPoints, weightViewType cubWeights, PointViewType cellCoords) const override
 Returns cubature points and weights (return arrays must be pre-sized/pre-allocated). More...
 
virtual ordinal_type getNumPoints () const override
 Returns the number of cubature points.
 
virtual ordinal_type getDimension () const override
 Returns dimension of integration domain.
 
virtual const char * getName () const override
 Returns cubature name.
 
 CubatureControlVolumeBoundary (const shards::CellTopology cellTopology, const ordinal_type sideIndex)
 
- Public Member Functions inherited from Intrepid2::Cubature< DeviceType, pointValueType, weightValueType >
virtual TensorPointDataType allocateCubaturePoints () const
 Returns a points container appropriate for passing to getCubature(). More...
 
virtual TensorWeightDataType allocateCubatureWeights () const
 Returns a weight container appropriate for passing to getCubature(). More...
 
virtual void getCubature (PointViewType, weightViewType) const
 Returns cubature points and weights (return arrays must be pre-sized/pre-allocated). More...
 
virtual void getCubature (const TensorPointDataType &tensorCubPoints, const TensorWeightDataType &tensorCubWeights) const
 Returns tensor cubature points and weights. For non-tensor cubatures, the tensor structures are trivial, thin wrappers around the data returned by getCubature(). The provided containers should be pre-allocated through calls to allocateCubaturePoints() and allocateCubatureWeights(). More...
 
virtual ordinal_type getAccuracy () const
 Returns dimension of the integration domain.
 

Protected Attributes

shards::CellTopology primaryCellTopo_
 The topology of the primary cell side.
 
shards::CellTopology subcvCellTopo_
 The topology of the sub-control volume.
 
ordinal_type degree_
 The degree of the polynomials that are integrated exactly.
 
ordinal_type sideIndex_
 Index of cell side.
 
Kokkos::View< ordinal_type
**, Kokkos::HostSpace > 
boundarySidesHost_
 
Kokkos::View< ordinal_type
**, Kokkos::LayoutRight,
DeviceType > 
sideNodeMap_
 
Kokkos::DynRankView
< pointValueType, DeviceType > 
sidePoints_
 

Detailed Description

template<typename DeviceType = void, typename pointValueType = double, typename weightValueType = double>
class Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >

Defines cubature (integration) rules over Neumann boundaries for control volume method.

Integration on Neumann boundaries for the control volume method requires integration points defined on primary cell sides. These points are not equivalent to control volume points on lower dimensional topologies and therefore require a separate class to define them.

Definition at line 72 of file Intrepid2_CubatureControlVolumeBoundary.hpp.

Constructor & Destructor Documentation

template<typename DT , typename PT , typename WT >
Intrepid2::CubatureControlVolumeBoundary< DT, PT, WT >::CubatureControlVolumeBoundary ( const shards::CellTopology  cellTopology,
const ordinal_type  sideIndex 
)

brief Constructor.

Parameters
cellTopology[in] - The topology of the primary cell.
cellSide[in] - The index of the boundary side of the primary cell

Definition at line 57 of file Intrepid2_CubatureControlVolumeBoundaryDef.hpp.

References Intrepid2::CellTools< DeviceType >::mapToReferenceSubcell().

Member Function Documentation

template<typename DT , typename PT , typename WT >
void Intrepid2::CubatureControlVolumeBoundary< DT, PT, WT >::getCubature ( PointViewType  cubPoints,
weightViewType  cubWeights,
PointViewType  cellCoords 
) const
overridevirtual

Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).

Parameters
cubPoints[out] - Array containing the cubature points.
cubWeights[out] - Array of corresponding cubature weights.
cellCoords[in] - Array of cell coordinates

Reimplemented from Intrepid2::Cubature< DeviceType, pointValueType, weightValueType >.

Definition at line 179 of file Intrepid2_CubatureControlVolumeBoundaryDef.hpp.

References Intrepid2::FunctionSpaceTools< DeviceType >::computeEdgeMeasure(), Intrepid2::FunctionSpaceTools< DeviceType >::computeFaceMeasure(), Intrepid2::CellTools< DeviceType >::getSubcvCoords(), Intrepid2::CellTools< DeviceType >::setJacobian(), and Intrepid2::CellTools< DeviceType >::setJacobianDet().


The documentation for this class was generated from the following files: