15 #ifndef __INTREPID2_HVOL_TET_CN_FEM_HPP__
16 #define __INTREPID2_HVOL_TET_CN_FEM_HPP__
21 #include "Teuchos_LAPACK.hpp"
52 template<EOperator opType>
54 template<
typename outputValueViewType,
55 typename inputPointViewType,
56 typename workViewType,
57 typename vinvViewType>
58 KOKKOS_INLINE_FUNCTION
60 getValues( outputValueViewType outputValues,
61 const inputPointViewType inputPoints,
63 const vinvViewType vinv );
66 KOKKOS_INLINE_FUNCTION
68 getWorkSizePerPoint(ordinal_type order) {
69 auto cardinality = getPnCardinality<3>(order);
76 return getDkCardinality<opType,3>()*cardinality;
81 template<
typename DeviceType, ordinal_type numPtsPerEval,
82 typename outputValueValueType,
class ...outputValueProperties,
83 typename inputPointValueType,
class ...inputPointProperties,
84 typename vinvValueType,
class ...vinvProperties>
86 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
87 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
88 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
89 const EOperator operatorType);
94 template<
typename outputValueViewType,
95 typename inputPointViewType,
96 typename vinvViewType,
97 typename workViewType,
99 ordinal_type numPtsEval>
101 outputValueViewType _outputValues;
102 const inputPointViewType _inputPoints;
103 const vinvViewType _vinv;
106 KOKKOS_INLINE_FUNCTION
107 Functor( outputValueViewType outputValues_,
108 inputPointViewType inputPoints_,
111 : _outputValues(outputValues_), _inputPoints(inputPoints_),
112 _vinv(vinv_), _work(work_) {}
114 KOKKOS_INLINE_FUNCTION
115 void operator()(
const size_type iter)
const {
119 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
120 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
122 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
124 auto vcprop = Kokkos::common_view_alloc_prop(_work);
125 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
128 case OPERATOR_VALUE : {
129 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_HVOL_TET_Cn_FEM::Functor) operator is not supported");
153 template<
typename DeviceType = void,
154 typename outputValueType = double,
155 typename pointValueType =
double>
157 :
public Basis<DeviceType,outputValueType,pointValueType> {
172 const EPointType pointType = POINTTYPE_EQUISPACED);
181 const EOperator operatorType = OPERATOR_VALUE)
const override {
182 #ifdef HAVE_INTREPID2_DEBUG
183 Intrepid2::getValues_HVOL_Args(outputValues,
190 Impl::Basis_HVOL_TET_Cn_FEM::
191 getValues<DeviceType,numPtsPerEval>( outputValues,
199 ordinal_type& perThreadSpaceSize,
201 const EOperator operatorType = OPERATOR_VALUE)
const override;
203 KOKKOS_INLINE_FUNCTION
208 const EOperator operatorType,
209 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
210 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
211 const ordinal_type subcellDim = -1,
212 const ordinal_type subcellOrdinal = -1)
const override;
218 #ifdef HAVE_INTREPID2_DEBUG
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
229 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
235 #ifdef HAVE_INTREPID2_DEBUG
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
240 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
243 Kokkos::deep_copy(dofCoeffs, 1.0);
249 Kokkos::deep_copy(vinv, this->
vinv_);
255 return "Intrepid2_HVOL_TET_Cn_FEM";
264 virtual HostBasisPtr<outputValueType,pointValueType>
273 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
274 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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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...
virtual const char * getName() const override
Returns basis name.
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.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual bool requireOrientation() const override
True if orientation is required.
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.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
See Intrepid2::Basis_HVOL_TET_Cn_FEM.
See Intrepid2::Basis_HVOL_TET_Cn_FEM.
Basis_HVOL_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Definition file for FEM basis functions of degree n for H(vol) functions on TET.
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...
See Intrepid2::Basis_HVOL_TET_Cn_FEM.
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Tetrahedron cell...
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.