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) {
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 {
116 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
117 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
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 );
130 Serial<opType>::getValues( output, input, work, _vinv );
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
139 Serial<opType>::getValues( output, input, work, _vinv );
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> {
159 using BasisBase = Basis<DeviceType,outputValueType,pointValueType>;
172 const EPointType pointType = POINTTYPE_EQUISPACED);
181 const EOperator operatorType = OPERATOR_VALUE)
const override {
182#ifdef HAVE_INTREPID2_DEBUG
190 Impl::Basis_HVOL_TET_Cn_FEM::
191 getValues<DeviceType,numPtsPerEval>( outputValues,
198 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
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";
273 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
274 EPointType pointType_;
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
KOKKOS_INLINE_FUNCTION ordinal_type getPnCardinality(ordinal_type n)
Returns cardinality of Polynomials of order n (P^n).
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns multiplicities of dx, dy, and dz based on the enumeration of the partial derivative,...
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....
Definition file for FEM basis functions of degree n for H(vol) functions on TET.
virtual const char * getName() const override
Returns basis name.
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...
virtual bool requireOrientation() const override
True if orientation is required.
Basis_HVOL_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Kokkos::DynRankView< scalarType, ExecutionSpace > vinv_
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
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
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
ordinal_type basisDegree_
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
ScalarTraits< PointScalar >::scalar_type scalarType
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
See Intrepid2::Basis_HVOL_TET_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HVOL_TET_Cn_FEM.