49 #ifndef __INTREPID2_HDIV_HEX_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_HEX_IN_FEM_HPP__
68 template<EOperator opType>
70 template<
typename outputValueViewType,
71 typename inputPointViewType,
72 typename workViewType,
73 typename vinvViewType>
74 KOKKOS_INLINE_FUNCTION
76 getValues( outputValueViewType outputValues,
77 const inputPointViewType inputPoints,
79 const vinvViewType vinvLine,
80 const vinvViewType vinvBubble );
82 KOKKOS_INLINE_FUNCTION
84 getWorkSizePerPoint(ordinal_type order) {
85 return 2*getPnCardinality<1>(order)+2*getPnCardinality<1>(order-1);
89 template<
typename DeviceType, ordinal_type numPtsPerEval,
90 typename outputValueValueType,
class ...outputValueProperties,
91 typename inputPointValueType,
class ...inputPointProperties,
92 typename vinvValueType,
class ...vinvProperties>
94 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
95 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
96 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
97 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
98 const EOperator operatorType );
103 template<
typename outputValueViewType,
104 typename inputPointViewType,
105 typename vinvViewType,
106 typename workViewType,
108 ordinal_type numPtsEval>
110 outputValueViewType _outputValues;
111 const inputPointViewType _inputPoints;
112 const vinvViewType _vinvLine;
113 const vinvViewType _vinvBubble;
116 KOKKOS_INLINE_FUNCTION
117 Functor( outputValueViewType outputValues_,
118 inputPointViewType inputPoints_,
119 vinvViewType vinvLine_,
120 vinvViewType vinvBubble_,
122 : _outputValues(outputValues_), _inputPoints(inputPoints_),
123 _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
125 KOKKOS_INLINE_FUNCTION
126 void operator()(
const size_type iter)
const {
130 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
131 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
133 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
135 auto vcprop = Kokkos::common_view_alloc_prop(_work);
136 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
139 case OPERATOR_VALUE : {
140 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
144 case OPERATOR_DIV : {
145 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
150 INTREPID2_TEST_FOR_ABORT(
true,
151 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::Functor) operator is not supported." );
167 template<
typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType =
double>
171 :
public Basis<DeviceType,outputValueType,pointValueType> {
180 const EPointType pointType = POINTTYPE_EQUISPACED);
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE )
const override {
193 #ifdef HAVE_INTREPID2_DEBUG
200 constexpr ordinal_type numPtsPerEval = 1;
201 Impl::Basis_HDIV_HEX_In_FEM::
202 getValues<DeviceType,numPtsPerEval>( outputValues,
212 #ifdef HAVE_INTREPID2_DEBUG
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
217 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
223 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
230 #ifdef HAVE_INTREPID2_DEBUG
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
235 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
241 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
247 return "Intrepid2_HDIV_HEX_In_FEM";
268 if(subCellDim == 2) {
269 return Teuchos::rcp(
new
273 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
283 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinvLine_, vinvBubble_;
284 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 HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Implementation of the default H(div)-compatible FEM basis on Hexahedron cell.
virtual bool requireOrientation() const override
True if orientation is required.
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...
Basis_HDIV_HEX_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
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< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
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.
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_HEX_In_FEM.
See Intrepid2::Basis_HDIV_HEX_In_FEM.
See Intrepid2::Basis_HDIV_HEX_In_FEM.
Hexahedron topology, 8 nodes.