16 #ifndef __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
17 #define __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
74 template<
bool serendipity>
81 template<EOperator opType>
83 template<
typename OutputViewType,
84 typename inputViewType>
85 KOKKOS_INLINE_FUNCTION
87 getValues( OutputViewType output,
88 const inputViewType input );
92 template<
typename DeviceType,
93 typename outputValueValueType,
class ...outputValueProperties,
94 typename inputPointValueType,
class ...inputPointProperties>
96 getValues(
const typename DeviceType::execution_space& space,
97 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
98 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
99 const EOperator operatorType);
104 template<
typename outputValueViewType,
105 typename inputPointViewType,
108 outputValueViewType _outputValues;
109 const inputPointViewType _inputPoints;
111 KOKKOS_INLINE_FUNCTION
112 Functor( outputValueViewType outputValues_,
113 inputPointViewType inputPoints_ )
114 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
116 KOKKOS_INLINE_FUNCTION
117 void operator()(
const ordinal_type pt)
const {
119 case OPERATOR_VALUE : {
120 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
121 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
131 case OPERATOR_MAX : {
132 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
133 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
138 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
139 opType != OPERATOR_GRAD &&
140 opType != OPERATOR_CURL &&
141 opType != OPERATOR_D1 &&
142 opType != OPERATOR_D2 &&
143 opType != OPERATOR_D3 &&
144 opType != OPERATOR_D4 &&
145 opType != OPERATOR_MAX,
146 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::Serial::getValues) operator is not supported");
154 template<
bool serendipity,
156 typename outputValueType,
157 typename pointValueType>
181 const EOperator operatorType = OPERATOR_VALUE )
const override {
182 #ifdef HAVE_INTREPID2_DEBUG
184 Intrepid2::getValues_HGRAD_Args(outputValues,
191 template getValues<DeviceType>(space,
199 ordinal_type& perThreadSpaceSize,
201 const EOperator operatorType = OPERATOR_VALUE)
const override;
203 KOKKOS_INLINE_FUNCTION
208 const EOperator operatorType,
209 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
210 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
211 const ordinal_type subcellDim = -1,
212 const ordinal_type subcellOrdinal = -1)
const override;
217 #ifdef HAVE_INTREPID2_DEBUG
219 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
222 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
234 #ifdef HAVE_INTREPID2_DEBUG
236 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
239 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242 Kokkos::deep_copy(dofCoeffs, 1.0);
248 return serendipity ?
"Intrepid2_HGRAD_QUAD_I2_FEM" :
"Intrepid2_HGRAD_QUAD_C2_FEM";
259 BasisPtr<DeviceType, outputValueType, pointValueType>
261 if(subCellDim == 1) {
262 return Teuchos::rcp(
new
265 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
268 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
274 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
275 using Basis_HGRAD_QUAD_C2_FEM = Basis_HGRAD_QUAD_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
277 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
278 using Basis_HGRAD_QUAD_I2_FEM = Basis_HGRAD_QUAD_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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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...
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
See Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
See Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM.
Basis_HGRAD_QUAD_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.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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...
Definition file for FEM basis functions of degree 2 for H(grad) functions on QUAD cells...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
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.
See Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.