49 #ifndef __INTREPID2_HDIV_QUAD_I1_FEM_HPP__ 50 #define __INTREPID2_HDIV_QUAD_I1_FEM_HPP__ 115 template<EOperator opType>
117 template<
typename outputViewType,
118 typename inputViewType>
119 KOKKOS_INLINE_FUNCTION
121 getValues( outputViewType output,
122 const inputViewType input );
126 template<
typename ExecSpaceType,
127 typename outputValueValueType,
class ...outputValueProperties,
128 typename inputPointValueType,
class ...inputPointProperties>
130 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
131 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
137 template<
typename outputValueViewType,
138 typename inputPointViewType,
141 outputValueViewType _outputValues;
142 const inputPointViewType _inputPoints;
144 KOKKOS_INLINE_FUNCTION
145 Functor( outputValueViewType outputValues_,
146 inputPointViewType inputPoints_ )
147 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
149 KOKKOS_INLINE_FUNCTION
150 void operator()(
const ordinal_type pt)
const {
152 case OPERATOR_VALUE : {
153 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
154 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
159 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt);
160 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
165 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
166 opType != OPERATOR_DIV,
167 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::Serial::getValues) operator is not supported");
175 template<
typename ExecSpaceType = void,
176 typename outputValueType = double,
177 typename pointValueType =
double>
199 const EOperator operatorType = OPERATOR_VALUE )
const {
200 #ifdef HAVE_INTREPID2_DEBUG 208 Impl::Basis_HDIV_QUAD_I1_FEM::
209 getValues<ExecSpaceType>( outputValues,
217 #ifdef HAVE_INTREPID2_DEBUG 219 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
222 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
basisCellTopology_.getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
234 #ifdef HAVE_INTREPID2_DEBUG 236 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
239 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
245 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
251 return "Intrepid2_HDIV_QUAD_I1_FEM";
Implementation of the default H(div)-compatible FEM basis of degree 1 on Quadrilateral cell...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
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...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
ordinal_type getCardinality() const
Returns cardinality of the basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Quadrilateral topology, 4 nodes.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Definition file for FEM basis functions of degree 1 for H(div) functions on QUAD cells.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Basis_HDIV_QUAD_I1_FEM()
Constructor.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
virtual bool requireOrientation() const
True if orientation is required.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...