16 #ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
17 #define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
102 template<
bool serendipity>
109 template<EOperator opType>
111 template<
typename OutputViewType,
112 typename inputViewType>
113 KOKKOS_INLINE_FUNCTION
115 getValues( OutputViewType output,
116 const inputViewType input );
120 template<
typename DeviceType,
121 typename outputValueValueType,
class ...outputValueProperties,
122 typename inputPointValueType,
class ...inputPointProperties>
124 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
125 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
126 const EOperator operatorType);
131 template<
typename outputValueViewType,
132 typename inputPointViewType,
135 outputValueViewType _outputValues;
136 const inputPointViewType _inputPoints;
138 KOKKOS_INLINE_FUNCTION
139 Functor( outputValueViewType outputValues_,
140 inputPointViewType inputPoints_ )
141 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
142 KOKKOS_INLINE_FUNCTION
143 void operator()(
const ordinal_type pt)
const {
145 case OPERATOR_VALUE : {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
147 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
155 case OPERATOR_MAX : {
156 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
157 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
162 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
163 opType != OPERATOR_GRAD &&
164 opType != OPERATOR_D2 &&
165 opType != OPERATOR_MAX,
166 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::Serial::getValues) operator is not supported");
174 template<
bool serendipity,
175 typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType =
double>
195 getValues( OutputViewType outputValues,
196 const PointViewType inputPoints,
197 const EOperator operatorType = OPERATOR_VALUE )
const override {
198 #ifdef HAVE_INTREPID2_DEBUG
200 Intrepid2::getValues_HGRAD_Args(outputValues,
206 if constexpr (serendipity)
208 getValues<DeviceType>( outputValues,
213 getValues<DeviceType>( outputValues,
219 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
220 ordinal_type& perThreadSpaceSize,
221 const PointViewType inputPointsconst,
222 const EOperator operatorType = OPERATOR_VALUE)
const override;
224 KOKKOS_INLINE_FUNCTION
227 OutputViewType outputValues,
228 const PointViewType inputPoints,
229 const EOperator operatorType,
230 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
231 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
232 const ordinal_type subcellDim = -1,
233 const ordinal_type subcellOrdinal = -1)
const override;
237 getDofCoords( ScalarViewType dofCoords )
const override {
238 #ifdef HAVE_INTREPID2_DEBUG
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
249 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
254 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
255 #ifdef HAVE_INTREPID2_DEBUG
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
260 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
261 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
263 Kokkos::deep_copy(dofCoeffs, 1.0);
269 return serendipity ?
"Intrepid2_HGRAD_WEDGE_I2_FEM" :
"Intrepid2_HGRAD_WEDGE_C2_FEM";
280 BasisPtr<DeviceType,outputValueType,pointValueType>
282 if(subCellDim == 1) {
283 return Teuchos::rcp(
new
285 }
else if(subCellDim == 2) {
291 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
294 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
300 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
301 using Basis_HGRAD_WEDGE_C2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
303 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
304 using Basis_HGRAD_WEDGE_I2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
ordinal_type getCardinality() const
Returns cardinality of the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
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...
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.
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.
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 Wedge 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...
virtual const char * getName() const override
Returns basis name.
Basis_HGRAD_WEDGE_DEG2_FEM()
Constructor.
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.
Definition file for FEM basis functions of degree 2 for H(grad) functions on WEDGE cells...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.