49 #ifndef __INTREPID2_HCURL_WEDGE_I1_FEM_HPP__ 50 #define __INTREPID2_HCURL_WEDGE_I1_FEM_HPP__ 123 template<EOperator opType>
125 template<
typename OutputViewType,
126 typename inputViewType>
127 KOKKOS_INLINE_FUNCTION
129 getValues( OutputViewType output,
130 const inputViewType input );
134 template<
typename DeviceType,
135 typename outputValueValueType,
class ...outputValueProperties,
136 typename inputPointValueType,
class ...inputPointProperties>
138 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
139 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
145 template<
typename outputValueViewType,
146 typename inputPointViewType,
149 outputValueViewType _outputValues;
150 const inputPointViewType _inputPoints;
152 KOKKOS_INLINE_FUNCTION
153 Functor( outputValueViewType outputValues_,
154 inputPointViewType inputPoints_ )
155 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
157 KOKKOS_INLINE_FUNCTION
158 void operator()(
const ordinal_type pt)
const {
160 case OPERATOR_VALUE : {
161 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
162 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
166 case OPERATOR_CURL : {
167 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
168 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
173 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
174 opType != OPERATOR_CURL,
175 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
183 template<
typename DeviceType = void,
184 typename outputValueType = double,
185 typename pointValueType =
double>
204 getValues( OutputViewType outputValues,
205 const PointViewType inputPoints,
206 const EOperator operatorType = OPERATOR_VALUE )
const override {
207 #ifdef HAVE_INTREPID2_DEBUG 215 Impl::Basis_HCURL_WEDGE_I1_FEM::
216 getValues<DeviceType>( outputValues,
223 getDofCoords( ScalarViewType dofCoords )
const override {
224 #ifdef HAVE_INTREPID2_DEBUG 226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
229 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
basisCellTopology_.getDimension(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
235 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
240 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
241 #ifdef HAVE_INTREPID2_DEBUG 243 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
246 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
249 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
252 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
258 return "Intrepid2_HCURL_WEDGE_I1_FEM";
280 return Teuchos::rcp(
new 282 else if(subCellDim == 2) {
283 if((subCellOrd >=0) && (subCellOrd<3))
284 return Teuchos::rcp(
new 286 if((subCellOrd==3)||(subCellOrd==4))
287 return Teuchos::rcp(
new 291 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Wedge cell.
See Intrepid2::Basis_HCURL_WEDGE_I1_FEM.
void getValues_HCURL_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 HCURL-conforming FEM basis...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual bool requireOrientation() const override
True if orientation is required.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Quadrilateral cell...
See Intrepid2::Basis_HCURL_WEDGE_I1_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual const char * getName() const override
Returns basis name.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
See Intrepid2::Basis_HCURL_WEDGE_I1_FEM.
Basis_HCURL_WEDGE_I1_FEM()
Constructor.
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...
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
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(curl) functions on WEDGE cells...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
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...