49 #ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__ 50 #define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__ 131 template<EOperator opType>
133 template<
typename OutputViewType,
134 typename inputViewType>
135 KOKKOS_INLINE_FUNCTION
137 getValues( OutputViewType output,
138 const inputViewType input );
142 template<
typename DeviceType,
143 typename outputValueValueType,
class ...outputValueProperties,
144 typename inputPointValueType,
class ...inputPointProperties>
146 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
147 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
153 template<
typename outputValueViewType,
154 typename inputPointViewType,
157 outputValueViewType _outputValues;
158 const inputPointViewType _inputPoints;
160 KOKKOS_INLINE_FUNCTION
161 Functor( outputValueViewType outputValues_,
162 inputPointViewType inputPoints_ )
163 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
164 KOKKOS_INLINE_FUNCTION
165 void operator()(
const ordinal_type pt)
const {
167 case OPERATOR_VALUE : {
168 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
169 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
177 case OPERATOR_MAX : {
178 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
179 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
184 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
185 opType != OPERATOR_GRAD &&
186 opType != OPERATOR_D2 &&
187 opType != OPERATOR_MAX,
188 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::Serial::getValues) operator is not supported");
196 template<
typename DeviceType = void,
197 typename outputValueType = double,
198 typename pointValueType =
double>
216 getValues( OutputViewType outputValues,
217 const PointViewType inputPoints,
218 const EOperator operatorType = OPERATOR_VALUE )
const override {
219 #ifdef HAVE_INTREPID2_DEBUG 227 Impl::Basis_HGRAD_WEDGE_C2_FEM::
228 getValues<DeviceType>( outputValues,
235 getDofCoords( ScalarViewType dofCoords )
const override {
236 #ifdef HAVE_INTREPID2_DEBUG 238 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
241 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
247 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
252 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
253 #ifdef HAVE_INTREPID2_DEBUG 255 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
258 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
259 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
261 Kokkos::deep_copy(dofCoeffs, 1.0);
267 return "Intrepid2_HGRAD_WEDGE_C2_FEM";
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
See Intrepid2::Basis_HGRAD_WEDGE_C2_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_C2_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_C2_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Basis_HGRAD_WEDGE_C2_FEM()
Constructor.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
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...
virtual const char * getName() const override
Returns basis name.
Definition file for FEM basis functions of degree 2 for H(grad) functions on WEDGE cells...
Wedge topology, 18 nodes.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
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.