49#ifndef __INTREPID2_HDIV_TRI_IN_FEM_HPP__
50#define __INTREPID2_HDIV_TRI_IN_FEM_HPP__
57#include "Teuchos_LAPACK.hpp"
89#define CardinalityHDivTri(order) (order*(order+2))
103 template<EOperator opType>
105 template<
typename outputValueViewType,
106 typename inputPointViewType,
107 typename workViewType,
108 typename vinvViewType>
109 KOKKOS_INLINE_FUNCTION
111 getValues( outputValueViewType outputValues,
112 const inputPointViewType inputPoints,
114 const vinvViewType vinv );
116 KOKKOS_INLINE_FUNCTION
118 getWorkSizePerPoint(ordinal_type order) {
119 auto cardinality = CardinalityHDivTri(order);
124 return 5*cardinality;
126 return getDkCardinality<opType,2>()*cardinality;
131 template<
typename DeviceType, ordinal_type numPtsPerEval,
132 typename outputValueValueType,
class ...outputValueProperties,
133 typename inputPointValueType,
class ...inputPointProperties,
134 typename vinvValueType,
class ...vinvProperties>
136 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
137 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
138 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
139 const EOperator operatorType);
144 template<
typename outputValueViewType,
145 typename inputPointViewType,
146 typename vinvViewType,
147 typename workViewType,
149 ordinal_type numPtsEval>
151 outputValueViewType _outputValues;
152 const inputPointViewType _inputPoints;
153 const vinvViewType _coeffs;
156 KOKKOS_INLINE_FUNCTION
157 Functor( outputValueViewType outputValues_,
158 inputPointViewType inputPoints_,
159 vinvViewType coeffs_,
161 : _outputValues(outputValues_), _inputPoints(inputPoints_),
162 _coeffs(coeffs_), _work(work_) {}
164 KOKKOS_INLINE_FUNCTION
165 void operator()(
const size_type iter)
const {
169 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
170 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
172 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
174 auto vcprop = Kokkos::common_view_alloc_prop(_work);
175 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
179 case OPERATOR_VALUE : {
180 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
185 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
190 INTREPID2_TEST_FOR_ABORT(
true,
191 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
200template<
typename DeviceType = void,
201 typename outputValueType = double,
202 typename pointValueType =
double>
204 :
public Basis<DeviceType,outputValueType,pointValueType> {
213 const EPointType pointType = POINTTYPE_EQUISPACED);
228 const PointViewType inputPoints,
229 const EOperator operatorType = OPERATOR_VALUE)
const override {
230#ifdef HAVE_INTREPID2_DEBUG
238 Impl::Basis_HDIV_TRI_In_FEM::
239 getValues<DeviceType,numPtsPerEval>( outputValues,
248#ifdef HAVE_INTREPID2_DEBUG
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
253 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
254 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
256 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
257 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
259 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
265#ifdef HAVE_INTREPID2_DEBUG
267 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
268 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
270 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
271 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
273 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
274 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
276 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
280 getExpansionCoeffs( ScalarViewType coeffs )
const {
282 Kokkos::deep_copy(coeffs, this->
coeffs_);
288 return "Intrepid2_HDIV_TRI_In_FEM";
308 if(subCellDim == 1) {
309 return Teuchos::rcp(
new
313 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
324 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 ...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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.
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.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
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.