48 #ifndef __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__ 49 #define __INTREPID2_HDIV_WEDGE_I1_FEM_HPP__ 97 template<EOperator opType>
99 template<
typename OutputViewType,
100 typename inputViewType>
101 KOKKOS_INLINE_FUNCTION
103 getValues( OutputViewType output,
104 const inputViewType input );
108 template<
typename DeviceType,
109 typename outputValueValueType,
class ...outputValueProperties,
110 typename inputPointValueType,
class ...inputPointProperties>
112 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
113 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
119 template<
typename outputValueViewType,
120 typename inputPointViewType,
123 outputValueViewType _outputValues;
124 const inputPointViewType _inputPoints;
126 KOKKOS_INLINE_FUNCTION
127 Functor( outputValueViewType outputValues_,
128 inputPointViewType inputPoints_ )
129 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
131 KOKKOS_INLINE_FUNCTION
132 void operator()(
const ordinal_type pt)
const {
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
136 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
140 case OPERATOR_DIV : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
142 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
147 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
148 opType != OPERATOR_DIV ,
149 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
157 template<
typename DeviceType = void,
158 typename outputValueType = double,
159 typename pointValueType =
double>
178 getValues( OutputViewType outputValues,
179 const PointViewType inputPoints,
180 const EOperator operatorType = OPERATOR_VALUE )
const override {
181 #ifdef HAVE_INTREPID2_DEBUG 189 Impl::Basis_HDIV_WEDGE_I1_FEM::
190 getValues<DeviceType>( outputValues,
197 getDofCoords( ScalarViewType dofCoords )
const override {
198 #ifdef HAVE_INTREPID2_DEBUG 200 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
201 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
203 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
204 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
206 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
basisCellTopology_.getDimension(), std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
209 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
215 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
216 #ifdef HAVE_INTREPID2_DEBUG 218 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
221 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
224 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HDIV_WEDGE_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
227 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
233 return "Intrepid2_HDIV_WEDGE_I1_FEM";
254 if(subCellDim == 2) {
255 if((subCellOrd >=0) && (subCellOrd<3))
256 return Teuchos::rcp(
new 258 if((subCellOrd==3)||(subCellOrd==4))
259 return Teuchos::rcp(
new 262 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
See Intrepid2::Basis_HDIV_WEDGE_I1_FEM.
See Intrepid2::Basis_HDIV_WEDGE_I1_FEM.
Basis_HDIV_WEDGE_I1_FEM()
Constructor.
void getValues_HDIV_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 HDIV-conforming FEM basis...
See Intrepid2::Basis_HDIV_WEDGE_I1_FEM.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
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...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
Definition file for FEM basis functions of degree 1 for H(div) functions on WEDGE cells...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
virtual bool requireOrientation() const override
True if orientation is required.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge 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.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...