49 #ifndef __INTREPID2_HGRAD_TET_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_TET_C2_FEM_HPP__
112 template<EOperator opType>
114 template<
typename OutputViewType,
115 typename inputViewType>
116 KOKKOS_INLINE_FUNCTION
118 getValues( OutputViewType output,
119 const inputViewType input );
123 template<
typename DeviceType,
124 typename outputValueValueType,
class ...outputValueProperties,
125 typename inputPointValueType,
class ...inputPointProperties>
127 getValues(
const typename DeviceType::execution_space& space,
128 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130 const EOperator operatorType);
135 template<
typename outputValueViewType,
136 typename inputPointViewType,
139 outputValueViewType _outputValues;
140 const inputPointViewType _inputPoints;
142 KOKKOS_INLINE_FUNCTION
143 Functor( outputValueViewType outputValues_,
144 inputPointViewType inputPoints_ )
145 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
147 KOKKOS_INLINE_FUNCTION
148 void operator()(
const ordinal_type pt)
const {
150 case OPERATOR_VALUE : {
151 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
152 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
159 case OPERATOR_MAX : {
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
161 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
166 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
167 opType != OPERATOR_GRAD &&
168 opType != OPERATOR_MAX,
169 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
177 template<
typename DeviceType = void,
178 typename outputValueType = double,
179 typename pointValueType =
double>
203 const EOperator operatorType = OPERATOR_VALUE )
const override {
204 #ifdef HAVE_INTREPID2_DEBUG
206 Intrepid2::getValues_HGRAD_Args(outputValues,
212 Impl::Basis_HGRAD_TET_C2_FEM::
213 getValues<DeviceType>(space,
222 #ifdef HAVE_INTREPID2_DEBUG
224 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
227 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
230 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
231 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
233 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
240 #ifdef HAVE_INTREPID2_DEBUG
242 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
245 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
248 Kokkos::deep_copy(dofCoeffs, 1.0);
254 return "Intrepid2_HGRAD_TET_C2_FEM";
265 BasisPtr<DeviceType,outputValueType,pointValueType>
267 if(subCellDim == 1) {
269 }
else if(subCellDim == 2) {
273 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
276 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
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...
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.
Basis_HGRAD_TET_C2_FEM()
Constructor.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
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...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual const char * getName() const override
Returns basis name.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TET cells.
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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< 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.
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.