16#ifndef __INTREPID2_HDIV_TRI_IN_FEM_HPP__
17#define __INTREPID2_HDIV_TRI_IN_FEM_HPP__
24#include "Teuchos_LAPACK.hpp"
56#define CardinalityHDivTri(order) (order*(order+2))
70 template<EOperator opType>
72 template<
typename outputValueViewType,
73 typename inputPointViewType,
74 typename workViewType,
75 typename vinvViewType>
76 KOKKOS_INLINE_FUNCTION
78 getValues( outputValueViewType outputValues,
79 const inputPointViewType inputPoints,
81 const vinvViewType vinv );
83 KOKKOS_INLINE_FUNCTION
85 getWorkSizePerPoint(ordinal_type order) {
86 auto cardinality = CardinalityHDivTri(order);
93 return getDkCardinality<opType,2>()*cardinality;
98 template<
typename DeviceType, ordinal_type numPtsPerEval,
99 typename outputValueValueType,
class ...outputValueProperties,
100 typename inputPointValueType,
class ...inputPointProperties,
101 typename vinvValueType,
class ...vinvProperties>
103 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
104 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
105 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
106 const EOperator operatorType);
111 template<
typename outputValueViewType,
112 typename inputPointViewType,
113 typename vinvViewType,
114 typename workViewType,
116 ordinal_type numPtsEval>
118 outputValueViewType _outputValues;
119 const inputPointViewType _inputPoints;
120 const vinvViewType _coeffs;
123 KOKKOS_INLINE_FUNCTION
124 Functor( outputValueViewType outputValues_,
125 inputPointViewType inputPoints_,
126 vinvViewType coeffs_,
128 : _outputValues(outputValues_), _inputPoints(inputPoints_),
129 _coeffs(coeffs_), _work(work_) {}
131 KOKKOS_INLINE_FUNCTION
132 void operator()(
const size_type iter)
const {
136 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
137 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
139 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
141 auto vcprop = Kokkos::common_view_alloc_prop(_work);
142 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
146 case OPERATOR_VALUE : {
147 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
157 INTREPID2_TEST_FOR_ABORT(
true,
158 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
167template<
typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType =
double>
171 :
public Basis<DeviceType,outputValueType,pointValueType> {
181 const EPointType pointType = POINTTYPE_EQUISPACED);
193 getValues( OutputViewType outputValues,
194 const PointViewType inputPoints,
195 const EOperator operatorType = OPERATOR_VALUE)
const override {
196#ifdef HAVE_INTREPID2_DEBUG
204 Impl::Basis_HDIV_TRI_In_FEM::
205 getValues<DeviceType,numPtsPerEval>( outputValues,
212 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
213 ordinal_type& perThreadSpaceSize,
214 const PointViewType inputPointsconst,
215 const EOperator operatorType = OPERATOR_VALUE)
const override;
217 KOKKOS_INLINE_FUNCTION
220 OutputViewType outputValues,
221 const PointViewType inputPoints,
222 const EOperator operatorType,
223 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
224 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
225 const ordinal_type subcellDim = -1,
226 const ordinal_type subcellOrdinal = -1)
const override;
231#ifdef HAVE_INTREPID2_DEBUG
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
236 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
239 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
242 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
248#ifdef HAVE_INTREPID2_DEBUG
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
253 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
254 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
256 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
257 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
259 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
263 getExpansionCoeffs( ScalarViewType coeffs )
const {
265 Kokkos::deep_copy(coeffs, this->
coeffs_);
271 return "Intrepid2_HDIV_TRI_In_FEM";
291 if(subCellDim == 1) {
292 return Teuchos::rcp(
new
296 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
307 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
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 TRI cells.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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 getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual const char * getName() const override
Returns basis name.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual bool requireOrientation() const override
True if orientation is required.
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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.
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 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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
Triangle topology, 3 nodes.