49#ifndef __INTREPID2_HDIV_QUAD_I1_FEM_HPP__
50#define __INTREPID2_HDIV_QUAD_I1_FEM_HPP__
116 template<EOperator opType>
118 template<
typename OutputViewType,
119 typename inputViewType>
120 KOKKOS_INLINE_FUNCTION
122 getValues( OutputViewType output,
123 const inputViewType input );
127 template<
typename DeviceType,
128 typename outputValueValueType,
class ...outputValueProperties,
129 typename inputPointValueType,
class ...inputPointProperties>
131 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
132 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
133 const EOperator operatorType);
138 template<
typename outputValueViewType,
139 typename inputPointViewType,
142 outputValueViewType _outputValues;
143 const inputPointViewType _inputPoints;
145 KOKKOS_INLINE_FUNCTION
146 Functor( outputValueViewType outputValues_,
147 inputPointViewType inputPoints_ )
148 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
150 KOKKOS_INLINE_FUNCTION
151 void operator()(
const ordinal_type pt)
const {
153 case OPERATOR_VALUE : {
154 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
155 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt);
161 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
166 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
167 opType != OPERATOR_DIV,
168 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::Serial::getValues) operator is not supported");
176 template<
typename DeviceType = void,
177 typename outputValueType = double,
178 typename pointValueType =
double>
199 const PointViewType inputPoints,
200 const EOperator operatorType = OPERATOR_VALUE )
const override {
201#ifdef HAVE_INTREPID2_DEBUG
209 Impl::Basis_HDIV_QUAD_I1_FEM::
210 getValues<DeviceType>( outputValues,
218#ifdef HAVE_INTREPID2_DEBUG
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
229 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
235#ifdef HAVE_INTREPID2_DEBUG
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
240 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
243 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
246 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
252 return "Intrepid2_HDIV_QUAD_I1_FEM";
273 return Teuchos::rcp(
new
276 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Definition file for FEM basis functions of degree 1 for H(div) functions on QUAD cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on Quadrilateral cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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 getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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 const char * getName() const override
Returns basis name.
virtual bool requireOrientation() const override
True if orientation is required.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_HDIV_QUAD_I1_FEM()
Constructor.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral,...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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 ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
Quadrilateral topology, 4 nodes.