|
Intrepid2
|
Utilizes cubature (integration) rules contained in the library Polylib (Spencer Sherwin, Aeronautics, Imperial College London) within Intrepid. More...
#include <Intrepid2_CubaturePolylib.hpp>
Public Types | |
| typedef CubatureDirect< DeviceType, pointValueType, weightValueType >::PointViewType | PointViewType |
| typedef CubatureDirect< DeviceType, pointValueType, weightValueType >::weightViewType | weightViewType |
Public Types inherited from Intrepid2::CubatureDirect< DeviceType, pointValueType, weightValueType > | |
| 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 | |
| CubaturePolylib (const ordinal_type degree, const EPolyType polytype=POLYTYPE_GAUSS, const double alpha=0.0, const double beta=0.0) | |
| virtual const char * | getName () const override |
| Returns cubature name. | |
Public Member Functions inherited from Intrepid2::CubatureDirect< DeviceType, pointValueType, weightValueType > | |
| virtual void | getCubature (PointViewType cubPoints, weightViewType cubWeights) const override |
| Returns cubature points and weights (return arrays must be pre-sized/pre-allocated). | |
| virtual ordinal_type | getNumPoints () const override |
| Returns the number of cubature points. | |
| virtual ordinal_type | getDimension () const override |
| Returns dimension of integration domain. | |
| virtual ordinal_type | getAccuracy () const override |
| Returns max. degree of polynomials that are integrated exactly. The return vector has size 1. | |
| CubatureDirect (const CubatureDirect &b) | |
| CubatureDirect & | operator= (const CubatureDirect &b) |
| CubatureDirect (const ordinal_type degree, const ordinal_type dimension) | |
Public Member Functions inherited from Intrepid2::Cubature< DeviceType, pointValueType, weightValueType > | |
| virtual TensorPointDataType | allocateCubaturePoints () const |
| Returns a points container appropriate for passing to getCubature(). | |
| virtual TensorWeightDataType | allocateCubatureWeights () const |
| Returns a weight container appropriate for passing to getCubature(). | |
| virtual void | getCubature (PointViewType, weightViewType, PointViewType) const |
| Returns cubature points and weights on physical cells (return arrays must be pre-sized/pre-allocated). | |
| 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(). | |
Additional Inherited Members | |
Protected Member Functions inherited from Intrepid2::CubatureDirect< DeviceType, pointValueType, weightValueType > | |
| template<typename cubPointValueType , class ... cubPointProperties, typename cubWeightValueType , class ... cubWeightProperties> | |
| void | getCubatureFromData (Kokkos::DynRankView< cubPointValueType, cubPointProperties... > cubPoints, Kokkos::DynRankView< cubWeightValueType, cubWeightProperties... > cubWeights, const CubatureData cubData) const |
| Returns cubature points and weights. | |
Protected Attributes inherited from Intrepid2::CubatureDirect< DeviceType, pointValueType, weightValueType > | |
| ordinal_type | degree_ |
| The degree of polynomials that are integrated exactly by this cubature rule. | |
| ordinal_type | dimension_ |
| Dimension of integration domain. | |
| CubatureData | cubatureData_ |
| Cubature data on device. | |
Utilizes cubature (integration) rules contained in the library Polylib (Spencer Sherwin, Aeronautics, Imperial College London) within Intrepid.
They are based on zeros of Jacobi polynomials, e.g. Legendre (alpha=beta=0, default), Chebyshev (alpha=beta=-0.5), etc. They are given on the interval [-1,1] and are optimal with respect to the following requirements, yielding 4 subclasses:
Definition at line 43 of file Intrepid2_CubaturePolylib.hpp.
| typedef CubatureDirect<DeviceType,pointValueType,weightValueType>::PointViewType Intrepid2::CubaturePolylib< DeviceType, pointValueType, weightValueType >::PointViewType |
Definition at line 46 of file Intrepid2_CubaturePolylib.hpp.
| typedef CubatureDirect<DeviceType,pointValueType,weightValueType>::weightViewType Intrepid2::CubaturePolylib< DeviceType, pointValueType, weightValueType >::weightViewType |
Definition at line 47 of file Intrepid2_CubaturePolylib.hpp.
| Intrepid2::CubaturePolylib< DT, PT, WT >::CubaturePolylib | ( | const ordinal_type | degree, |
| const EPolyType | polytype = POLYTYPE_GAUSS, |
||
| const double | alpha = 0.0, |
||
| const double | beta = 0.0 |
||
| ) |
Definition at line 19 of file Intrepid2_CubaturePolylibDef.hpp.
|
inlineoverridevirtual |
Returns cubature name.
Reimplemented from Intrepid2::CubatureDirect< DeviceType, pointValueType, weightValueType >.
Definition at line 58 of file Intrepid2_CubaturePolylib.hpp.