49 #ifndef __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_C2_FEM_HPP__
94 typedef struct Triangle<6> cell_topology_type;
98 template<EOperator opType>
100 template<
typename OutputViewType,
101 typename inputViewType>
102 KOKKOS_INLINE_FUNCTION
104 getValues( OutputViewType output,
105 const inputViewType input );
109 template<
typename DeviceType,
110 typename outputValueValueType,
class ...outputValueProperties,
111 typename inputPointValueType,
class ...inputPointProperties>
113 getValues(
const typename DeviceType::execution_space& space,
114 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
115 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
116 const EOperator operatorType);
121 template<
typename outputValueViewType,
122 typename inputPointViewType,
125 outputValueViewType _outputValues;
126 const inputPointViewType _inputPoints;
128 KOKKOS_INLINE_FUNCTION
129 Functor( outputValueViewType outputValues_,
130 inputPointViewType inputPoints_ )
131 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
133 KOKKOS_INLINE_FUNCTION
134 void operator()(
const ordinal_type pt)
const {
136 case OPERATOR_VALUE : {
137 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
138 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
146 case OPERATOR_MAX : {
147 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
148 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
153 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
154 opType != OPERATOR_GRAD &&
155 opType != OPERATOR_CURL &&
156 opType != OPERATOR_D1 &&
157 opType != OPERATOR_D2 &&
158 opType != OPERATOR_MAX,
159 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::Serial::getValues) operator is not supported");
167 template<
typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType =
double>
193 const EOperator operatorType = OPERATOR_VALUE )
const override {
194 #ifdef HAVE_INTREPID2_DEBUG
196 Intrepid2::getValues_HGRAD_Args(outputValues,
202 Impl::Basis_HGRAD_TRI_C2_FEM::
203 getValues<DeviceType>(space,
212 #ifdef HAVE_INTREPID2_DEBUG
214 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
217 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
223 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
229 #ifdef HAVE_INTREPID2_DEBUG
231 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
234 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
237 Kokkos::deep_copy(dofCoeffs, 1.0);
243 return "Intrepid2_HGRAD_TRI_C2_FEM";
255 BasisPtr<DeviceType,outputValueType,pointValueType>
260 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
263 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.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
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 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 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.