48 #ifndef __INTREPID2_HVOL_QUAD_CN_FEM_HPP__ 49 #define __INTREPID2_HVOL_QUAD_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 );
81 KOKKOS_INLINE_FUNCTION
83 getWorkSizePerPoint(ordinal_type order) {
return 3*getPnCardinality<1>(order); }
86 template<
typename DeviceType, ordinal_type numPtsPerEval,
87 typename outputValueValueType,
class ...outputValueProperties,
88 typename inputPointValueType,
class ...inputPointProperties,
89 typename vinvValueType,
class ...vinvProperties>
91 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
99 template<
typename outputValueViewType,
100 typename inputPointViewType,
101 typename vinvViewType,
102 typename workViewType,
104 ordinal_type numPtsEval>
106 outputValueViewType _outputValues;
107 const inputPointViewType _inputPoints;
108 const vinvViewType _vinv;
110 const ordinal_type _opDn;
112 KOKKOS_INLINE_FUNCTION
113 Functor( outputValueViewType outputValues_,
114 inputPointViewType inputPoints_,
117 const ordinal_type opDn_ = 0 )
118 : _outputValues(outputValues_), _inputPoints(inputPoints_),
119 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
121 KOKKOS_INLINE_FUNCTION
122 void operator()(
const size_type iter)
const {
126 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
127 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
129 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
131 auto vcprop = Kokkos::common_view_alloc_prop(_work);
132 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
135 case OPERATOR_VALUE : {
136 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
146 INTREPID2_TEST_FOR_ABORT(
true,
147 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Functor) operator is not supported");
162 template<
typename DeviceType = void,
163 typename outputValueType = double,
164 typename pointValueType =
double>
166 :
public Basis<DeviceType,outputValueType,pointValueType> {
175 const EPointType pointType = POINTTYPE_EQUISPACED);
185 getValues( OutputViewType outputValues,
186 const PointViewType inputPoints,
187 const EOperator operatorType = OPERATOR_VALUE )
const override {
188 #ifdef HAVE_INTREPID2_DEBUG 196 Impl::Basis_HVOL_QUAD_Cn_FEM::
197 getValues<DeviceType,numPtsPerEval>( outputValues,
205 getDofCoords( ScalarViewType dofCoords )
const override {
206 #ifdef HAVE_INTREPID2_DEBUG 208 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
211 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
217 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
222 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
223 #ifdef HAVE_INTREPID2_DEBUG 225 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
228 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
231 Kokkos::deep_copy(dofCoeffs, 1.0);
237 return "Intrepid2_HVOL_QUAD_Cn_FEM";
254 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by 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. ...
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 (...
Quadrilateral topology, 4 nodes.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type), allowing host access to input and output views, and ensuring host execution of basis evaluation.
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...
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)
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...
EPointType
Enumeration of types of point distributions in Intrepid.
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
virtual const char * getName() const override
Returns basis name.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.