49#ifndef __INTREPID2_HDIV_QUAD_IN_FEM_HPP__
50#define __INTREPID2_HDIV_QUAD_IN_FEM_HPP__
69 template<EOperator opType>
71 template<
typename outputValueViewType,
72 typename inputPointViewType,
73 typename workViewType,
74 typename vinvViewType>
75 KOKKOS_INLINE_FUNCTION
77 getValues( outputValueViewType outputValues,
78 const inputPointViewType inputPoints,
80 const vinvViewType vinvLine,
81 const vinvViewType vinvBubble );
83 KOKKOS_INLINE_FUNCTION
85 getWorkSizePerPoint(ordinal_type order) {
86 return 2*getPnCardinality<1>(order)+getPnCardinality<1>(order-1);
90 template<
typename DeviceType, ordinal_type numPtsPerEval,
91 typename outputValueValueType,
class ...outputValueProperties,
92 typename inputPointValueType,
class ...inputPointProperties,
93 typename vinvValueType,
class ...vinvProperties>
95 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
96 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
97 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
98 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
99 const EOperator operatorType );
104 template<
typename outputValueViewType,
105 typename inputPointViewType,
106 typename vinvViewType,
107 typename workViewType,
109 ordinal_type numPtsEval>
111 outputValueViewType _outputValues;
112 const inputPointViewType _inputPoints;
113 const vinvViewType _vinvLine;
114 const vinvViewType _vinvBubble;
117 KOKKOS_INLINE_FUNCTION
118 Functor( outputValueViewType outputValues_,
119 inputPointViewType inputPoints_,
120 vinvViewType vinvLine_,
121 vinvViewType vinvBubble_,
123 : _outputValues(outputValues_), _inputPoints(inputPoints_),
124 _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
126 KOKKOS_INLINE_FUNCTION
127 void operator()(
const size_type iter)
const {
131 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
132 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
134 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
136 auto vcprop = Kokkos::common_view_alloc_prop(_work);
137 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
140 case OPERATOR_VALUE : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
145 case OPERATOR_DIV : {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
151 INTREPID2_TEST_FOR_ABORT(
true,
152 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::Functor) operator is not supported");
165 template<
typename DeviceType = void,
166 typename outputValueType = double,
167 typename pointValueType =
double>
169 :
public Basis<DeviceType,outputValueType,pointValueType> {
178 const EPointType pointType = POINTTYPE_EQUISPACED);
189 const PointViewType inputPoints,
190 const EOperator operatorType = OPERATOR_VALUE )
const override {
191#ifdef HAVE_INTREPID2_DEBUG
199 Impl::Basis_HDIV_QUAD_In_FEM::
200 getValues<DeviceType,numPtsPerEval>( outputValues,
209#ifdef HAVE_INTREPID2_DEBUG
211 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
214 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
217 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
220 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
226#ifdef HAVE_INTREPID2_DEBUG
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
231 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
234 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
237 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
243 return "Intrepid2_HDIV_QUAD_In_FEM";
264 return Teuchos::rcp(
new
268 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
278 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinvLine_, vinvBubble_;
279 EPointType pointType_;
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 n for H(div) functions on QUAD cells.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Implementation of the default H(div)-compatible FEM basis on Quadrilateral 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 void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
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...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,...
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 basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HDIV_QUAD_In_FEM.
See Intrepid2::Basis_HDIV_QUAD_In_FEM.
Quadrilateral topology, 4 nodes.