49#ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
135 template<
bool serendipity>
138 typedef struct Wedge<serendipity ? 15 : 18> cell_topology_type;
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() );
181 Serial<opType>::getValues( output, input );
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() );
191 Serial<opType>::getValues( output, input );
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>
229 const PointViewType inputPoints,
230 const EOperator operatorType = OPERATOR_VALUE )
const override {
231#ifdef HAVE_INTREPID2_DEBUG
239 if constexpr (serendipity)
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_);
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 if constexpr (serendipity)
286 return "Intrepid2_HGRAD_WEDGE_I2_Serendipity_FEM";
288 return "Intrepid2_HGRAD_WEDGE_C2_FEM";
301 if(subCellDim == 1) {
302 return Teuchos::rcp(
new
304 }
else if(subCellDim == 2) {
310 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
319 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
320 using Basis_HGRAD_WEDGE_C2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
322 template<
typename DeviceType =
void,
typename outputValueType =
double,
typename po
intValueType =
double>
323 using Basis_HGRAD_WEDGE_I2_Serendipity_FEM = Basis_HGRAD_WEDGE_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.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Definition file for FEM basis functions of degree 2 for H(grad) functions on WEDGE cells.
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.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
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 void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Basis_HGRAD_WEDGE_DEG2_FEM()
Constructor.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM.