49 #ifndef __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
152 template<
bool serendipity>
159 template<EOperator opType>
161 template<
typename OutputViewType,
162 typename inputViewType>
163 KOKKOS_INLINE_FUNCTION
165 getValues( OutputViewType output,
166 const inputViewType input );
170 template<
typename DeviceType,
171 typename outputValueValueType,
class ...outputValueProperties,
172 typename inputPointValueType,
class ...inputPointProperties>
174 getValues(
const typename DeviceType::execution_space& space,
175 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
176 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
177 const EOperator operatorType);
182 template<
typename outputValueViewType,
183 typename inputPointViewType,
186 outputValueViewType _outputValues;
187 const inputPointViewType _inputPoints;
189 KOKKOS_INLINE_FUNCTION
190 Functor( outputValueViewType outputValues_,
191 inputPointViewType inputPoints_ )
192 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
194 KOKKOS_INLINE_FUNCTION
195 void operator()(
const ordinal_type pt)
const {
197 case OPERATOR_VALUE : {
198 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
199 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
208 case OPERATOR_MAX : {
209 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
210 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
215 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
216 opType != OPERATOR_GRAD &&
217 opType != OPERATOR_D1 &&
218 opType != OPERATOR_D2 &&
219 opType != OPERATOR_D3 &&
220 opType != OPERATOR_D4 &&
221 opType != OPERATOR_MAX,
222 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::Serial::getValues) operator is not supported");
230 template<
bool serendipity,
232 typename outputValueType,
233 typename pointValueType>
258 const EOperator operatorType = OPERATOR_VALUE )
const override {
259 #ifdef HAVE_INTREPID2_DEBUG
261 Intrepid2::getValues_HGRAD_Args(outputValues,
267 if constexpr (serendipity)
269 getValues<DeviceType>(space,
275 getValues<DeviceType>(space,
284 #ifdef HAVE_INTREPID2_DEBUG
286 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
287 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
289 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
290 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
292 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
293 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
295 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
301 #ifdef HAVE_INTREPID2_DEBUG
303 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
304 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
306 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
307 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
309 Kokkos::deep_copy(dofCoeffs, 1.0);
315 return serendipity ?
"Intrepid2_HGRAD_HEX_I2_FEM" :
"Intrepid2_HGRAD_HEX_C2_FEM";
326 BasisPtr<DeviceType,outputValueType,pointValueType>
328 if(subCellDim == 1) {
329 return Teuchos::rcp(
new
331 }
else if(subCellDim == 2) {
332 return Teuchos::rcp(
new
335 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
338 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
344 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
345 using Basis_HGRAD_HEX_C2_FEM = Basis_HGRAD_HEX_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
347 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
348 using Basis_HGRAD_HEX_I2_FEM = Basis_HGRAD_HEX_DEG2_FEM<true, DeviceType, 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.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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...
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 Hexahedron 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...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.
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 HEX cells.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
virtual const char * getName() const override
Returns basis name.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Basis_HGRAD_HEX_DEG2_FEM()
Constructor.
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...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.