16 #ifndef __INTREPID2_HDIV_TET_IN_FEM_HPP__
17 #define __INTREPID2_HDIV_TET_IN_FEM_HPP__
24 #include "Teuchos_LAPACK.hpp"
54 #define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
67 template<EOperator opType>
69 template<
typename outputValueViewType,
70 typename inputPointViewType,
71 typename workViewType,
72 typename vinvViewType>
73 KOKKOS_INLINE_FUNCTION
75 getValues( outputValueViewType outputValues,
76 const inputPointViewType inputPoints,
78 const vinvViewType vinv );
81 KOKKOS_INLINE_FUNCTION
83 getWorkSizePerPoint(ordinal_type order) {
84 auto cardinality = CardinalityHDivTet(order);
90 return getDkCardinality<opType,3>()*cardinality;
95 template<
typename DeviceType, ordinal_type numPtsPerEval,
96 typename outputValueValueType,
class ...outputValueProperties,
97 typename inputPointValueType,
class ...inputPointProperties,
98 typename vinvValueType,
class ...vinvProperties>
100 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
101 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
102 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
103 const EOperator operatorType);
108 template<
typename outputValueViewType,
109 typename inputPointViewType,
110 typename vinvViewType,
111 typename workViewType,
113 ordinal_type numPtsEval>
115 outputValueViewType _outputValues;
116 const inputPointViewType _inputPoints;
117 const vinvViewType _coeffs;
120 KOKKOS_INLINE_FUNCTION
121 Functor( outputValueViewType outputValues_,
122 inputPointViewType inputPoints_,
123 vinvViewType coeffs_,
125 : _outputValues(outputValues_), _inputPoints(inputPoints_),
126 _coeffs(coeffs_), _work(work_) {}
128 KOKKOS_INLINE_FUNCTION
129 void operator()(
const size_type iter)
const {
133 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
134 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
136 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
138 auto vcprop = Kokkos::common_view_alloc_prop(_work);
139 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
142 case OPERATOR_VALUE : {
143 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
148 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
153 INTREPID2_TEST_FOR_ABORT(
true,
154 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
163 template<
typename DeviceType = void,
164 typename outputValueType = double,
165 typename pointValueType =
double>
167 :
public Basis<DeviceType,outputValueType,pointValueType> {
177 const EPointType pointType = POINTTYPE_EQUISPACED);
187 getValues( OutputViewType outputValues,
188 const PointViewType inputPoints,
189 const EOperator operatorType = OPERATOR_VALUE)
const override {
190 #ifdef HAVE_INTREPID2_DEBUG
191 Intrepid2::getValues_HDIV_Args(outputValues,
198 Impl::Basis_HDIV_TET_In_FEM::
199 getValues<DeviceType,numPtsPerEval>( outputValues,
206 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
207 ordinal_type& perThreadSpaceSize,
208 const PointViewType inputPointsconst,
209 const EOperator operatorType = OPERATOR_VALUE)
const override;
211 KOKKOS_INLINE_FUNCTION
214 OutputViewType outputValues,
215 const PointViewType inputPoints,
216 const EOperator operatorType,
217 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
218 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
219 const ordinal_type subcellDim = -1,
220 const ordinal_type subcellOrdinal = -1)
const override;
224 getDofCoords( ScalarViewType dofCoords )
const override {
225 #ifdef HAVE_INTREPID2_DEBUG
227 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
230 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
231 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
236 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
241 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
242 #ifdef HAVE_INTREPID2_DEBUG
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
247 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
253 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
257 getExpansionCoeffs( ScalarViewType coeffs )
const {
259 Kokkos::deep_copy(coeffs, this->
coeffs_);
265 return "Intrepid2_HDIV_TET_In_FEM";
283 BasisPtr<DeviceType,outputValueType,pointValueType>
286 if(subCellDim == 2) {
287 return Teuchos::rcp(
new
291 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
294 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
301 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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...
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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.
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HDIV_TET_In_FEM.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
See Intrepid2::Basis_HDIV_TET_In_FEM.
virtual const char * getName() const override
Returns basis name.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
See Intrepid2::Basis_HDIV_TET_In_FEM.
Header file for the abstract base class Intrepid2::Basis.