16#ifndef __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
17#define __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
119 template<
bool serendipity>
122 typedef struct Hexahedron<serendipity ? 20 : 27> cell_topology_type;
126 template<EOperator opType>
128 template<
typename OutputViewType,
129 typename inputViewType>
130 KOKKOS_INLINE_FUNCTION
132 getValues( OutputViewType output,
133 const inputViewType input );
137 template<
typename DeviceType,
138 typename outputValueValueType,
class ...outputValueProperties,
139 typename inputPointValueType,
class ...inputPointProperties>
141 getValues(
const typename DeviceType::execution_space& space,
142 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
143 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
144 const EOperator operatorType);
149 template<
typename outputValueViewType,
150 typename inputPointViewType,
153 outputValueViewType _outputValues;
154 const inputPointViewType _inputPoints;
156 KOKKOS_INLINE_FUNCTION
157 Functor( outputValueViewType outputValues_,
158 inputPointViewType inputPoints_ )
159 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
161 KOKKOS_INLINE_FUNCTION
162 void operator()(
const ordinal_type pt)
const {
164 case OPERATOR_VALUE : {
165 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
166 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
167 Serial<opType>::getValues( output, input );
175 case OPERATOR_MAX : {
176 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
177 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
178 Serial<opType>::getValues( output, input );
182 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
183 opType != OPERATOR_GRAD &&
184 opType != OPERATOR_D1 &&
185 opType != OPERATOR_D2 &&
186 opType != OPERATOR_D3 &&
187 opType != OPERATOR_D4 &&
188 opType != OPERATOR_MAX,
189 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::Serial::getValues) operator is not supported");
197 template<
bool serendipity,
199 typename outputValueType,
200 typename pointValueType>
203 using BasisBase = Basis<DeviceType,outputValueType,pointValueType>;
225 const EOperator operatorType = OPERATOR_VALUE )
const override {
226#ifdef HAVE_INTREPID2_DEBUG
234 if constexpr (serendipity)
250 ordinal_type& perThreadSpaceSize,
252 const EOperator operatorType = OPERATOR_VALUE)
const override;
254 KOKKOS_INLINE_FUNCTION
259 const EOperator operatorType,
260 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
261 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
262 const ordinal_type subcellDim = -1,
263 const ordinal_type subcellOrdinal = -1)
const override;
268#ifdef HAVE_INTREPID2_DEBUG
270 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
271 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
273 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
274 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
276 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
277 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
279 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
285#ifdef HAVE_INTREPID2_DEBUG
287 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
288 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
290 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
291 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
293 Kokkos::deep_copy(dofCoeffs, 1.0);
299 return serendipity ?
"Intrepid2_HGRAD_HEX_I2_FEM" :
"Intrepid2_HGRAD_HEX_C2_FEM";
312 if(subCellDim == 1) {
313 return Teuchos::rcp(
new
315 }
else if(subCellDim == 2) {
316 return Teuchos::rcp(
new
319 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
328 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
329 using Basis_HGRAD_HEX_C2_FEM = Basis_HGRAD_HEX_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
331 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
332 using Basis_HGRAD_HEX_I2_FEM = Basis_HGRAD_HEX_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
Header file for the abstract base class Intrepid2::Basis.
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....
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Definition file for FEM basis functions of degree 2 for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
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 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...
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...
Basis_HGRAD_HEX_DEG2_FEM()
Constructor.
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 const char * getName() const override
Returns basis name.
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 Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
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
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
typename DeviceType::execution_space ExecutionSpace
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.