49 #ifndef __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
71 template<EOperator opType>
73 template<
typename outputValueViewType,
74 typename inputPointViewType,
75 typename workViewType,
76 typename vinvViewType>
77 KOKKOS_INLINE_FUNCTION
79 getValues( outputValueViewType outputValues,
80 const inputPointViewType inputPoints,
82 const vinvViewType vinv,
83 const ordinal_type operatorDn = 0 );
86 template<
typename DeviceType, ordinal_type numPtsPerEval,
87 typename outputValueValueType,
class ...outputValueProperties,
88 typename inputPointValueType,
class ...inputPointProperties,
89 typename vinvValueType,
class ...vinvProperties>
91 getValues(
const typename DeviceType::execution_space& space,
92 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
93 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
94 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
95 const EOperator operatorType);
101 template<
typename outputValueViewType,
102 typename inputPointViewType,
103 typename vinvViewType,
104 typename workViewType,
106 ordinal_type numPtsEval>
108 outputValueViewType _outputValues;
109 const inputPointViewType _inputPoints;
110 const vinvViewType _vinv;
112 const ordinal_type _opDn;
114 KOKKOS_INLINE_FUNCTION
115 Functor( outputValueViewType outputValues_,
116 inputPointViewType inputPoints_,
119 const ordinal_type opDn_ = 0 )
120 : _outputValues(outputValues_), _inputPoints(inputPoints_),
121 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
123 KOKKOS_INLINE_FUNCTION
124 void operator()(
const size_type iter)
const {
128 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
129 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
131 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
133 auto vcprop = Kokkos::common_view_alloc_prop(_work);
134 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
137 case OPERATOR_VALUE : {
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
144 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
149 INTREPID2_TEST_FOR_ABORT(
true,
150 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
164 template<
typename DeviceType = void,
165 typename outputValueType = double,
166 typename pointValueType =
double>
168 :
public Basis<DeviceType,outputValueType,pointValueType> {
182 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
191 const EPointType pointType = POINTTYPE_EQUISPACED);
200 const EOperator operatorType = OPERATOR_VALUE )
const override {
201 #ifdef HAVE_INTREPID2_DEBUG
202 Intrepid2::getValues_HGRAD_Args(outputValues,
209 Impl::Basis_HGRAD_QUAD_Cn_FEM::
210 getValues<DeviceType,numPtsPerEval>(space,
220 #ifdef HAVE_INTREPID2_DEBUG
222 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
225 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
231 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
237 #ifdef HAVE_INTREPID2_DEBUG
239 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
242 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
245 Kokkos::deep_copy(dofCoeffs, 1.0);
251 return "Intrepid2_HGRAD_QUAD_Cn_FEM";
260 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
261 getVandermondeInverse()
const {
266 getWorkSizePerPoint(
const EOperator operatorType)
const {
278 BasisPtr<DeviceType,outputValueType,pointValueType>
280 if(subCellDim == 1) {
281 return Teuchos::rcp(
new
285 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
288 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
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...
ordinal_type getCardinality() const
Returns cardinality of the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree n for H(grad) functions on QUAD cells...
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM work is a rank 1 view having the same value_type of inputPoint...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.