16 #ifndef __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
17 #define __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
38 template<EOperator opType>
40 template<
typename outputValueViewType,
41 typename inputPointViewType,
42 typename workViewType,
43 typename vinvViewType>
44 KOKKOS_INLINE_FUNCTION
46 getValues( outputValueViewType outputValues,
47 const inputPointViewType inputPoints,
49 const vinvViewType vinv,
50 const ordinal_type operatorDn = 0 );
53 template<
typename DeviceType, ordinal_type numPtsPerEval,
54 typename outputValueValueType,
class ...outputValueProperties,
55 typename inputPointValueType,
class ...inputPointProperties,
56 typename vinvValueType,
class ...vinvProperties>
58 getValues(
const typename DeviceType::execution_space& space,
59 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
60 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
61 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
62 const EOperator operatorType);
68 template<
typename outputValueViewType,
69 typename inputPointViewType,
70 typename vinvViewType,
71 typename workViewType,
73 ordinal_type numPtsEval>
75 outputValueViewType _outputValues;
76 const inputPointViewType _inputPoints;
77 const vinvViewType _vinv;
79 const ordinal_type _opDn;
81 KOKKOS_INLINE_FUNCTION
82 Functor( outputValueViewType outputValues_,
83 inputPointViewType inputPoints_,
86 const ordinal_type opDn_ = 0 )
87 : _outputValues(outputValues_), _inputPoints(inputPoints_),
88 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
90 KOKKOS_INLINE_FUNCTION
91 void operator()(
const size_type iter)
const {
95 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
96 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
98 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
100 auto vcprop = Kokkos::common_view_alloc_prop(_work);
101 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
104 case OPERATOR_VALUE : {
105 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
111 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
116 INTREPID2_TEST_FOR_ABORT(
true,
117 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
131 template<
typename DeviceType = void,
132 typename outputValueType = double,
133 typename pointValueType =
double>
135 :
public Basis<DeviceType,outputValueType,pointValueType> {
149 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
158 const EPointType pointType = POINTTYPE_EQUISPACED);
167 const EOperator operatorType = OPERATOR_VALUE )
const override {
168 #ifdef HAVE_INTREPID2_DEBUG
169 Intrepid2::getValues_HGRAD_Args(outputValues,
176 Impl::Basis_HGRAD_QUAD_Cn_FEM::
177 getValues<DeviceType,numPtsPerEval>(space,
186 ordinal_type& perThreadSpaceSize,
188 const EOperator operatorType = OPERATOR_VALUE)
const override;
190 KOKKOS_INLINE_FUNCTION
195 const EOperator operatorType,
196 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
197 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
198 const ordinal_type subcellDim = -1,
199 const ordinal_type subcellOrdinal = -1)
const override;
204 #ifdef HAVE_INTREPID2_DEBUG
206 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
209 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
210 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
212 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
215 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
221 #ifdef HAVE_INTREPID2_DEBUG
223 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
226 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
229 Kokkos::deep_copy(dofCoeffs, 1.0);
235 return "Intrepid2_HGRAD_QUAD_Cn_FEM";
244 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
245 getVandermondeInverse()
const {
250 getWorkSizePerPoint(
const EOperator operatorType)
const {
262 BasisPtr<DeviceType,outputValueType,pointValueType>
264 if(subCellDim == 1) {
265 return Teuchos::rcp(
new
269 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
272 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.
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.
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.
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 void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, 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...
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.