16 #ifndef __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
17 #define __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
61 typedef struct Triangle<6> cell_topology_type;
65 template<EOperator opType>
67 template<
typename OutputViewType,
68 typename inputViewType>
69 KOKKOS_INLINE_FUNCTION
71 getValues( OutputViewType output,
72 const inputViewType input );
76 template<
typename DeviceType,
77 typename outputValueValueType,
class ...outputValueProperties,
78 typename inputPointValueType,
class ...inputPointProperties>
80 getValues(
const typename DeviceType::execution_space& space,
81 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
82 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
83 const EOperator operatorType);
88 template<
typename outputValueViewType,
89 typename inputPointViewType,
92 outputValueViewType _outputValues;
93 const inputPointViewType _inputPoints;
95 KOKKOS_INLINE_FUNCTION
96 Functor( outputValueViewType outputValues_,
97 inputPointViewType inputPoints_ )
98 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
100 KOKKOS_INLINE_FUNCTION
101 void operator()(
const ordinal_type pt)
const {
103 case OPERATOR_VALUE : {
104 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
105 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
113 case OPERATOR_MAX : {
114 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
115 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
120 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
121 opType != OPERATOR_GRAD &&
122 opType != OPERATOR_CURL &&
123 opType != OPERATOR_D1 &&
124 opType != OPERATOR_D2 &&
125 opType != OPERATOR_MAX,
126 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::Serial::getValues) operator is not supported");
134 template<
typename DeviceType = void,
135 typename outputValueType = double,
136 typename pointValueType =
double>
160 const EOperator operatorType = OPERATOR_VALUE )
const override {
161 #ifdef HAVE_INTREPID2_DEBUG
163 Intrepid2::getValues_HGRAD_Args(outputValues,
169 Impl::Basis_HGRAD_TRI_C2_FEM::
170 getValues<DeviceType>(space,
178 ordinal_type& perThreadSpaceSize,
180 const EOperator operatorType = OPERATOR_VALUE)
const override;
182 KOKKOS_INLINE_FUNCTION
187 const EOperator operatorType,
188 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
189 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
190 const ordinal_type subcellDim = -1,
191 const ordinal_type subcellOrdinal = -1)
const override;
196 #ifdef HAVE_INTREPID2_DEBUG
198 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
199 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
201 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
202 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
204 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
205 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
207 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
213 #ifdef HAVE_INTREPID2_DEBUG
215 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
218 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
221 Kokkos::deep_copy(dofCoeffs, 1.0);
227 return "Intrepid2_HGRAD_TRI_C2_FEM";
239 BasisPtr<DeviceType,outputValueType,pointValueType>
244 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
247 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
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_TRI_C2_FEM.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.
See Intrepid2::Basis_HGRAD_TRI_C2_FEM.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle 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...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
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.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
See Intrepid2::Basis_HGRAD_TRI_C2_FEM.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TRI cells.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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 getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Basis_HGRAD_TRI_C2_FEM()
Constructor.
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.
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.