15 #ifndef __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
16 #define __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
34 template<EOperator opType>
36 template<
typename outputValueViewType,
37 typename inputPointViewType,
38 typename workViewType,
39 typename vinvViewType>
40 KOKKOS_INLINE_FUNCTION
42 getValues( outputValueViewType outputValues,
43 const inputPointViewType inputPoints,
45 const vinvViewType vinv,
46 const ordinal_type operatorDn = 0 );
48 KOKKOS_INLINE_FUNCTION
50 getWorkSizePerPoint(ordinal_type order) {
return 3*getPnCardinality<1>(order); }
53 template<
typename DeviceType, ordinal_type numPtsPerEval,
54 typename outputValueValueType,
class ...outputValueProperties,
55 typename inputPointValueType,
class ...inputPointProperties,
56 typename vinvValueType,
class ...vinvProperties>
58 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
59 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
60 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
61 const EOperator operatorType);
66 template<
typename outputValueViewType,
67 typename inputPointViewType,
68 typename vinvViewType,
69 typename workViewType,
71 ordinal_type numPtsEval>
73 outputValueViewType _outputValues;
74 const inputPointViewType _inputPoints;
75 const vinvViewType _vinv;
77 const ordinal_type _opDn;
79 KOKKOS_INLINE_FUNCTION
80 Functor( outputValueViewType outputValues_,
81 inputPointViewType inputPoints_,
84 const ordinal_type opDn_ = 0 )
85 : _outputValues(outputValues_), _inputPoints(inputPoints_),
86 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
88 KOKKOS_INLINE_FUNCTION
89 void operator()(
const size_type iter)
const {
93 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
94 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
96 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
98 auto vcprop = Kokkos::common_view_alloc_prop(_work);
99 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
102 case OPERATOR_VALUE : {
103 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
108 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
113 INTREPID2_TEST_FOR_ABORT(
true,
114 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Functor) operator is not supported");
129 template<
typename DeviceType = void,
130 typename outputValueType = double,
131 typename pointValueType =
double>
133 :
public Basis<DeviceType,outputValueType,pointValueType> {
148 const EPointType pointType = POINTTYPE_EQUISPACED);
156 const EOperator operatorType = OPERATOR_VALUE )
const override {
157 #ifdef HAVE_INTREPID2_DEBUG
158 Intrepid2::getValues_HVOL_Args(outputValues,
165 Impl::Basis_HVOL_QUAD_Cn_FEM::
166 getValues<DeviceType,numPtsPerEval>( outputValues,
174 ordinal_type& perThreadSpaceSize,
176 const EOperator operatorType = OPERATOR_VALUE)
const override;
178 KOKKOS_INLINE_FUNCTION
183 const EOperator operatorType,
184 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
185 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
186 const ordinal_type subcellDim = -1,
187 const ordinal_type subcellOrdinal = -1)
const override;
193 #ifdef HAVE_INTREPID2_DEBUG
195 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
196 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
198 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
199 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
201 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
202 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
204 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
210 #ifdef HAVE_INTREPID2_DEBUG
212 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
215 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
218 Kokkos::deep_copy(dofCoeffs, 1.0);
224 return "Intrepid2_HVOL_QUAD_Cn_FEM";
233 virtual HostBasisPtr<outputValueType,pointValueType>
241 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
242 EPointType pointType_;
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Definition file for FEM basis functions of degree n for H(vol) functions on QUAD. ...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements ...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
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< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.