49 #ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
50 #define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
135 template<
bool serendipity>
142 template<EOperator opType>
144 template<
typename OutputViewType,
145 typename inputViewType>
146 KOKKOS_INLINE_FUNCTION
148 getValues( OutputViewType output,
149 const inputViewType input );
153 template<
typename DeviceType,
154 typename outputValueValueType,
class ...outputValueProperties,
155 typename inputPointValueType,
class ...inputPointProperties>
157 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
158 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
159 const EOperator operatorType);
164 template<
typename outputValueViewType,
165 typename inputPointViewType,
168 outputValueViewType _outputValues;
169 const inputPointViewType _inputPoints;
171 KOKKOS_INLINE_FUNCTION
172 Functor( outputValueViewType outputValues_,
173 inputPointViewType inputPoints_ )
174 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
175 KOKKOS_INLINE_FUNCTION
176 void operator()(
const ordinal_type pt)
const {
178 case OPERATOR_VALUE : {
179 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
180 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
188 case OPERATOR_MAX : {
189 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
190 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
195 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
196 opType != OPERATOR_GRAD &&
197 opType != OPERATOR_D2 &&
198 opType != OPERATOR_MAX,
199 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::Serial::getValues) operator is not supported");
207 template<
bool serendipity,
208 typename DeviceType = void,
209 typename outputValueType = double,
210 typename pointValueType =
double>
228 getValues( OutputViewType outputValues,
229 const PointViewType inputPoints,
230 const EOperator operatorType = OPERATOR_VALUE )
const override {
231 #ifdef HAVE_INTREPID2_DEBUG
233 Intrepid2::getValues_HGRAD_Args(outputValues,
239 if constexpr (serendipity)
241 getValues<DeviceType>( outputValues,
246 getValues<DeviceType>( outputValues,
253 getDofCoords( ScalarViewType dofCoords )
const override {
254 #ifdef HAVE_INTREPID2_DEBUG
256 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
257 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
259 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
260 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
262 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
265 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
270 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
271 #ifdef HAVE_INTREPID2_DEBUG
273 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
274 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
276 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
277 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
279 Kokkos::deep_copy(dofCoeffs, 1.0);
285 return serendipity ?
"Intrepid2_HGRAD_WEDGE_I2_FEM" :
"Intrepid2_HGRAD_WEDGE_C2_FEM";
296 BasisPtr<DeviceType,outputValueType,pointValueType>
298 if(subCellDim == 1) {
299 return Teuchos::rcp(
new
301 }
else if(subCellDim == 2) {
307 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
310 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
316 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
317 using Basis_HGRAD_WEDGE_C2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
319 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
320 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.