16#ifndef __INTREPID2_HGRAD_TET_C2_FEM_HPP__
17#define __INTREPID2_HGRAD_TET_C2_FEM_HPP__
79 template<EOperator opType>
81 template<
typename OutputViewType,
82 typename inputViewType>
83 KOKKOS_INLINE_FUNCTION
85 getValues( OutputViewType output,
86 const inputViewType input );
90 template<
typename DeviceType,
91 typename outputValueValueType,
class ...outputValueProperties,
92 typename inputPointValueType,
class ...inputPointProperties>
94 getValues(
const typename DeviceType::execution_space& space,
95 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
96 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
97 const EOperator operatorType);
102 template<
typename outputValueViewType,
103 typename inputPointViewType,
106 outputValueViewType _outputValues;
107 const inputPointViewType _inputPoints;
109 KOKKOS_INLINE_FUNCTION
110 Functor( outputValueViewType outputValues_,
111 inputPointViewType inputPoints_ )
112 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
114 KOKKOS_INLINE_FUNCTION
115 void operator()(
const ordinal_type pt)
const {
117 case OPERATOR_VALUE : {
118 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
119 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
126 case OPERATOR_MAX : {
127 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
128 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
133 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
134 opType != OPERATOR_GRAD &&
135 opType != OPERATOR_MAX,
136 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
144 template<
typename DeviceType = void,
145 typename outputValueType = double,
146 typename pointValueType =
double>
170 const EOperator operatorType = OPERATOR_VALUE )
const override {
171#ifdef HAVE_INTREPID2_DEBUG
179 Impl::Basis_HGRAD_TET_C2_FEM::
180 getValues<DeviceType>(space,
188 ordinal_type& perThreadSpaceSize,
190 const EOperator operatorType = OPERATOR_VALUE)
const override;
192 KOKKOS_INLINE_FUNCTION
197 const EOperator operatorType,
198 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
199 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
200 const ordinal_type subcellDim = -1,
201 const ordinal_type subcellOrdinal = -1)
const override;
206#ifdef HAVE_INTREPID2_DEBUG
208 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
211 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
217 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
224#ifdef HAVE_INTREPID2_DEBUG
226 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
229 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
232 Kokkos::deep_copy(dofCoeffs, 1.0);
238 return "Intrepid2_HGRAD_TET_C2_FEM";
251 if(subCellDim == 1) {
253 }
else if(subCellDim == 2) {
257 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Definition file for FEM basis functions of degree 2 for H(grad) functions on TET cells.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell.
Basis_HGRAD_TET_C2_FEM()
Constructor.
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 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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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 const char * getName() const override
Returns basis name.
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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
Tetrahedron topology, 10 nodes.