16 #ifndef __INTREPID2_HGRAD_TET_COMP12_FEM_HPP__
17 #define __INTREPID2_HGRAD_TET_COMP12_FEM_HPP__
79 template<
typename po
intValueType>
80 KOKKOS_INLINE_FUNCTION
83 const pointValueType y,
84 const pointValueType z );
89 template<EOperator opType>
91 template<
typename outputValueViewType,
92 typename inputPointViewType>
93 KOKKOS_INLINE_FUNCTION
95 getValues( outputValueViewType outputValues,
96 const inputPointViewType inputPoints );
100 template<
typename DeviceType,
101 typename outputValueValueType,
class ...outputValueProperties,
102 typename inputPointValueType,
class ...inputPointProperties>
104 getValues(
const typename DeviceType::execution_space& space,
105 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
106 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
107 const EOperator operatorType );
112 template<
typename outputValueViewType,
113 typename inputPointViewType,
116 outputValueViewType _outputValues;
117 const inputPointViewType _inputPoints;
119 KOKKOS_INLINE_FUNCTION
120 Functor( outputValueViewType outputValues_,
121 inputPointViewType inputPoints_ )
122 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
124 KOKKOS_INLINE_FUNCTION
125 void operator()(
const ordinal_type pt)
const {
127 case OPERATOR_VALUE : {
128 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
129 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
134 case OPERATOR_MAX : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
136 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
141 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
142 opType != OPERATOR_GRAD &&
143 opType != OPERATOR_MAX,
144 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::Functor::operator() operator is not supported");
152 template<
typename DeviceType = void,
153 typename outputValueType = double,
154 typename pointValueType =
double>
191 const EOperator operatorType = OPERATOR_VALUE )
const override {
192 #ifdef HAVE_INTREPID2_DEBUG
193 Intrepid2::getValues_HGRAD_Args(outputValues,
199 Impl::Basis_HGRAD_TET_COMP12_FEM::
200 getValues<DeviceType>(space,
208 ordinal_type& perThreadSpaceSize,
210 const EOperator operatorType = OPERATOR_VALUE)
const override;
212 KOKKOS_INLINE_FUNCTION
217 const EOperator operatorType,
218 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
219 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
220 const ordinal_type subcellDim = -1,
221 const ordinal_type subcellOrdinal = -1)
const override;
226 #ifdef HAVE_INTREPID2_DEBUG
228 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) rank = 2 required for dofCoords array");
231 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
234 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
237 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
243 return "Intrepid2_HGRAD_TET_COMP12_FEM";
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
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.
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
FEM basis evaluation on a reference Tetrahedron cell.
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...
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Basis_HGRAD_TET_COMP12_FEM()
Constructor.
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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
static KOKKOS_INLINE_FUNCTION ordinal_type getLocalSubTet(const pointValueType x, const pointValueType y, const pointValueType z)
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Definition file for the composite H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell with 1...
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.