49 #ifndef __INTREPID2_HDIV_TET_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TET_IN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
86 #define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
99 template<EOperator opType>
101 template<
typename outputValueViewType,
102 typename inputPointViewType,
103 typename workViewType,
104 typename vinvViewType>
105 KOKKOS_INLINE_FUNCTION
107 getValues( outputValueViewType outputValues,
108 const inputPointViewType inputPoints,
110 const vinvViewType vinv );
113 KOKKOS_INLINE_FUNCTION
115 getWorkSizePerPoint(ordinal_type order) {
116 auto cardinality = CardinalityHDivTet(order);
120 return 7*cardinality;
122 return getDkCardinality<opType,3>()*cardinality;
127 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
128 typename outputValueValueType,
class ...outputValueProperties,
129 typename inputPointValueType,
class ...inputPointProperties,
130 typename vinvValueType,
class ...vinvProperties>
132 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
133 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
134 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
135 const EOperator operatorType);
140 template<
typename outputValueViewType,
141 typename inputPointViewType,
142 typename vinvViewType,
143 typename workViewType,
145 ordinal_type numPtsEval>
147 outputValueViewType _outputValues;
148 const inputPointViewType _inputPoints;
149 const vinvViewType _coeffs;
152 KOKKOS_INLINE_FUNCTION
153 Functor( outputValueViewType outputValues_,
154 inputPointViewType inputPoints_,
155 vinvViewType coeffs_,
157 : _outputValues(outputValues_), _inputPoints(inputPoints_),
158 _coeffs(coeffs_), _work(work_) {}
160 KOKKOS_INLINE_FUNCTION
161 void operator()(
const size_type iter)
const {
165 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
166 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
168 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
170 auto vcprop = Kokkos::common_view_alloc_prop(_work);
171 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
174 case OPERATOR_VALUE : {
175 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
180 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
185 INTREPID2_TEST_FOR_ABORT(
true,
186 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
195 template<
typename ExecSpaceType = void,
196 typename outputValueType = double,
197 typename pointValueType =
double>
199 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
208 const EPointType pointType = POINTTYPE_EQUISPACED);
223 const EOperator operatorType = OPERATOR_VALUE)
const {
224 #ifdef HAVE_INTREPID2_DEBUG
225 Intrepid2::getValues_HDIV_Args(outputValues,
232 Impl::Basis_HDIV_TET_In_FEM::
233 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
242 #ifdef HAVE_INTREPID2_DEBUG
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
247 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
253 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
259 #ifdef HAVE_INTREPID2_DEBUG
261 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
262 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
264 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
265 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
267 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
268 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
270 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
274 getExpansionCoeffs( ScalarViewType coeffs )
const {
276 Kokkos::deep_copy(coeffs, this->
coeffs_);
282 return "Intrepid2_HDIV_TET_In_FEM";
295 Kokkos::DynRankView<scalarType,ExecSpaceType>
coeffs_;
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual const char * getName() const
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
See Intrepid2::Basis_HDIV_TET_In_FEM.
virtual bool requireOrientation() const
True if orientation is required.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
See Intrepid2::Basis_HDIV_TET_In_FEM.
See Intrepid2::Basis_HDIV_TET_In_FEM.
Header file for the abstract base class Intrepid2::Basis.