48#ifndef __INTREPID2_HVOL_HEX_CN_FEM_HPP__
49#define __INTREPID2_HVOL_HEX_CN_FEM_HPP__
62 typedef struct Hexahedron<8> cell_topology_type;
66 template<EOperator opType>
68 template<
typename outputValueViewType,
69 typename inputPointViewType,
70 typename workViewType,
71 typename vinvViewType>
72 KOKKOS_INLINE_FUNCTION
74 getValues( outputValueViewType outputValues,
75 const inputPointViewType inputPoints,
77 const vinvViewType vinv,
78 const ordinal_type operatorDn = 0 );
80 KOKKOS_INLINE_FUNCTION
85 template<
typename DeviceType, ordinal_type numPtsPerEval,
86 typename outputValueValueType,
class ...outputValueProperties,
87 typename inputPointValueType,
class ...inputPointProperties,
88 typename vinvValueType,
class ...vinvProperties>
90 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
91 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
92 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
93 const EOperator operatorType );
98 template<
typename outputValueViewType,
99 typename inputPointViewType,
100 typename vinvViewType,
101 typename workViewType,
103 ordinal_type numPtsEval>
105 outputValueViewType _outputValues;
106 const inputPointViewType _inputPoints;
107 const vinvViewType _vinv;
109 const ordinal_type _opDn;
111 KOKKOS_INLINE_FUNCTION
112 Functor( outputValueViewType outputValues_,
113 inputPointViewType inputPoints_,
116 const ordinal_type opDn_ = 0 )
117 : _outputValues(outputValues_), _inputPoints(inputPoints_),
118 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
120 KOKKOS_INLINE_FUNCTION
121 void operator()(
const size_type iter)
const {
122 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
123 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
125 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
126 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
128 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
130 auto vcprop = Kokkos::common_view_alloc_prop(_work);
131 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
136 Serial<opType>::getValues( output, input, work, _vinv );
140 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
141 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
145 INTREPID2_TEST_FOR_ABORT(
true,
146 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::Functor) operator is not supported");
167 template<
typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType =
double>
171 :
public Basis<DeviceType,outputValueType,pointValueType> {
180 const EPointType pointType = POINTTYPE_EQUISPACED);
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE )
const override {
193#ifdef HAVE_INTREPID2_DEBUG
201 Impl::Basis_HVOL_HEX_Cn_FEM::
202 getValues<DeviceType,numPtsPerEval>( outputValues,
212#ifdef HAVE_INTREPID2_DEBUG
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
217 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
223 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
229#ifdef HAVE_INTREPID2_DEBUG
231 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
234 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
237 Kokkos::deep_copy(dofCoeffs, 1.0);
243 return "Intrepid2_HVOL_HEX_Cn_FEM";
260 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
261 EPointType pointType_;
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
KOKKOS_INLINE_FUNCTION ordinal_type getPnCardinality(ordinal_type n)
Returns cardinality of Polynomials of order n (P^n).
void getValues_HVOL_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on HEX cells.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
virtual const char * getName() const override
Returns basis name.
Basis_HVOL_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecutionSpace > vinv_
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HVOL_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HVOL_HEX_Cn_FEM.