16#ifndef __INTREPID2_HDIV_QUAD_I1_FEM_HPP__
17#define __INTREPID2_HDIV_QUAD_I1_FEM_HPP__
83 template<EOperator opType>
85 template<
typename OutputViewType,
86 typename inputViewType>
87 KOKKOS_INLINE_FUNCTION
89 getValues( OutputViewType output,
90 const inputViewType input );
94 template<
typename DeviceType,
95 typename outputValueValueType,
class ...outputValueProperties,
96 typename inputPointValueType,
class ...inputPointProperties>
98 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
99 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
100 const EOperator operatorType);
105 template<
typename outputValueViewType,
106 typename inputPointViewType,
109 outputValueViewType _outputValues;
110 const inputPointViewType _inputPoints;
112 KOKKOS_INLINE_FUNCTION
113 Functor( outputValueViewType outputValues_,
114 inputPointViewType inputPoints_ )
115 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
117 KOKKOS_INLINE_FUNCTION
118 void operator()(
const ordinal_type pt)
const {
120 case OPERATOR_VALUE : {
121 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
122 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
127 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt);
128 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
133 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
134 opType != OPERATOR_DIV,
135 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::Serial::getValues) operator is not supported");
143 template<
typename DeviceType = void,
144 typename outputValueType = double,
145 typename pointValueType =
double>
166 getValues( OutputViewType outputValues,
167 const PointViewType inputPoints,
168 const EOperator operatorType = OPERATOR_VALUE )
const override {
169#ifdef HAVE_INTREPID2_DEBUG
177 Impl::Basis_HDIV_QUAD_I1_FEM::
178 getValues<DeviceType>( outputValues,
184 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
185 ordinal_type& perThreadSpaceSize,
186 const PointViewType inputPointsconst,
187 const EOperator operatorType = OPERATOR_VALUE)
const override;
189 KOKKOS_INLINE_FUNCTION
192 OutputViewType outputValues,
193 const PointViewType inputPoints,
194 const EOperator operatorType,
195 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
196 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
197 const ordinal_type subcellDim = -1,
198 const ordinal_type subcellOrdinal = -1)
const override;
203#ifdef HAVE_INTREPID2_DEBUG
205 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
208 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
211 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
214 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
220#ifdef HAVE_INTREPID2_DEBUG
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
225 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
231 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
237 return "Intrepid2_HDIV_QUAD_I1_FEM";
258 return Teuchos::rcp(
new
261 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...
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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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.
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.