49 #ifndef __INTREPID2_HDIV_TET_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TET_IN_FEM_HPP__
57 #include "Teuchos_LAPACK.hpp"
87 #define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
100 template<EOperator opType>
102 template<
typename outputValueViewType,
103 typename inputPointViewType,
104 typename workViewType,
105 typename vinvViewType>
106 KOKKOS_INLINE_FUNCTION
108 getValues( outputValueViewType outputValues,
109 const inputPointViewType inputPoints,
111 const vinvViewType vinv );
114 KOKKOS_INLINE_FUNCTION
116 getWorkSizePerPoint(ordinal_type order) {
117 auto cardinality = CardinalityHDivTet(order);
121 return 7*cardinality;
123 return getDkCardinality<opType,3>()*cardinality;
128 template<
typename DeviceType, ordinal_type numPtsPerEval,
129 typename outputValueValueType,
class ...outputValueProperties,
130 typename inputPointValueType,
class ...inputPointProperties,
131 typename vinvValueType,
class ...vinvProperties>
133 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
134 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
135 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
136 const EOperator operatorType);
141 template<
typename outputValueViewType,
142 typename inputPointViewType,
143 typename vinvViewType,
144 typename workViewType,
146 ordinal_type numPtsEval>
148 outputValueViewType _outputValues;
149 const inputPointViewType _inputPoints;
150 const vinvViewType _coeffs;
153 KOKKOS_INLINE_FUNCTION
154 Functor( outputValueViewType outputValues_,
155 inputPointViewType inputPoints_,
156 vinvViewType coeffs_,
158 : _outputValues(outputValues_), _inputPoints(inputPoints_),
159 _coeffs(coeffs_), _work(work_) {}
161 KOKKOS_INLINE_FUNCTION
162 void operator()(
const size_type iter)
const {
166 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
167 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
169 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
171 auto vcprop = Kokkos::common_view_alloc_prop(_work);
172 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
175 case OPERATOR_VALUE : {
176 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
181 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
186 INTREPID2_TEST_FOR_ABORT(
true,
187 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
196 template<
typename DeviceType = void,
197 typename outputValueType = double,
198 typename pointValueType =
double>
200 :
public Basis<DeviceType,outputValueType,pointValueType> {
209 const EPointType pointType = POINTTYPE_EQUISPACED);
222 getValues( OutputViewType outputValues,
223 const PointViewType inputPoints,
224 const EOperator operatorType = OPERATOR_VALUE)
const override {
225 #ifdef HAVE_INTREPID2_DEBUG
226 Intrepid2::getValues_HDIV_Args(outputValues,
233 Impl::Basis_HDIV_TET_In_FEM::
234 getValues<DeviceType,numPtsPerEval>( outputValues,
242 getDofCoords( ScalarViewType dofCoords )
const override {
243 #ifdef HAVE_INTREPID2_DEBUG
245 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
248 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
249 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
251 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
254 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
259 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
260 #ifdef HAVE_INTREPID2_DEBUG
262 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
265 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
266 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
268 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
269 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
271 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
275 getExpansionCoeffs( ScalarViewType coeffs )
const {
277 Kokkos::deep_copy(coeffs, this->
coeffs_);
283 return "Intrepid2_HDIV_TET_In_FEM";
301 BasisPtr<DeviceType,outputValueType,pointValueType>
304 if(subCellDim == 2) {
305 return Teuchos::rcp(
new
309 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
312 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
320 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
ordinal_type getCardinality() const
Returns cardinality of the basis.
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, 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 bool requireOrientation() const override
True if orientation is required.
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...
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.
See Intrepid2::Basis_HDIV_TET_In_FEM.
Header file for the abstract base class Intrepid2::Basis.