49 #ifndef __INTREPID2_HGRAD_HEX_CN_FEM_HPP__ 50 #define __INTREPID2_HGRAD_HEX_CN_FEM_HPP__ 67 template<EOperator opType>
69 template<
typename outputValueViewType,
70 typename inputPointViewType,
71 typename workViewType,
72 typename vinvViewType>
73 KOKKOS_INLINE_FUNCTION
75 getValues( outputValueViewType outputValues,
76 const inputPointViewType inputPoints,
78 const vinvViewType vinv,
79 const ordinal_type operatorDn = 0 );
82 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
83 typename outputValueValueType,
class ...outputValueProperties,
84 typename inputPointValueType,
class ...inputPointProperties,
85 typename vinvValueType,
class ...vinvProperties>
87 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
88 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
89 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
95 template<
typename outputValueViewType,
96 typename inputPointViewType,
97 typename vinvViewType,
98 typename workViewType,
100 ordinal_type numPtsEval>
102 outputValueViewType _outputValues;
103 const inputPointViewType _inputPoints;
104 const vinvViewType _vinv;
106 const ordinal_type _opDn;
108 KOKKOS_INLINE_FUNCTION
109 Functor( outputValueViewType outputValues_,
110 inputPointViewType inputPoints_,
113 const ordinal_type opDn_ = 0 )
114 : _outputValues(outputValues_), _inputPoints(inputPoints_),
115 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
117 KOKKOS_INLINE_FUNCTION
118 void operator()(
const size_type iter)
const {
122 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
123 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
125 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
127 auto vcprop = Kokkos::common_view_alloc_prop(_work);
128 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
131 case OPERATOR_VALUE : {
132 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
143 INTREPID2_TEST_FOR_ABORT(
true,
144 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::Functor) operator is not supported");
164 template<
typename ExecSpaceType = void,
165 typename outputValueType = double,
166 typename pointValueType =
double>
168 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
180 Kokkos::DynRankView<typename scalarViewType::value_type,ExecSpaceType>
vinv_;
187 const EPointType pointType = POINTTYPE_EQUISPACED);
195 const EOperator operatorType = OPERATOR_VALUE )
const {
196 #ifdef HAVE_INTREPID2_DEBUG 204 Impl::Basis_HGRAD_HEX_Cn_FEM::
205 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
215 #ifdef HAVE_INTREPID2_DEBUG 217 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
232 #ifdef HAVE_INTREPID2_DEBUG 234 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
237 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
240 Kokkos::deep_copy(dofCoeffs, 1.0);
246 return "Intrepid2_HGRAD_HEX_Cn_FEM";
255 Kokkos::DynRankView<typename scalarViewType::const_value_type,ExecSpaceType>
256 getVandermondeInverse() {
261 getWorkSizePerPoint(
const EOperator operatorType) {
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
Basis_HGRAD_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Definition file for basis function of degree n for H(grad) functions on HEX cells.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
void getValues_HGRAD_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 HGRAD-conforming FEM basis...
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
EPointType
Enumeration of types of point distributions in Intrepid.
Hexahedron topology, 8 nodes.
Kokkos::DynRankView< typename scalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
virtual const char * getName() const
Returns basis name.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.