49 #ifndef __INTREPID2_HDIV_TRI_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TRI_IN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
88 #define CardinalityHDivTri(order) (order*(order+2))
97 typedef struct Triangle<3> cell_topology_type;
102 template<EOperator opType>
104 template<
typename outputValueViewType,
105 typename inputPointViewType,
106 typename workViewType,
107 typename vinvViewType>
108 KOKKOS_INLINE_FUNCTION
110 getValues( outputValueViewType outputValues,
111 const inputPointViewType inputPoints,
113 const vinvViewType vinv );
115 KOKKOS_INLINE_FUNCTION
117 getWorkSizePerPoint(ordinal_type order) {
118 auto cardinality = CardinalityHDivTri(order);
123 return 5*cardinality;
125 return getDkCardinality<opType,2>()*cardinality;
130 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
131 typename outputValueValueType,
class ...outputValueProperties,
132 typename inputPointValueType,
class ...inputPointProperties,
133 typename vinvValueType,
class ...vinvProperties>
135 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
136 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
137 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
138 const EOperator operatorType);
143 template<
typename outputValueViewType,
144 typename inputPointViewType,
145 typename vinvViewType,
146 typename workViewType,
148 ordinal_type numPtsEval>
150 outputValueViewType _outputValues;
151 const inputPointViewType _inputPoints;
152 const vinvViewType _coeffs;
155 KOKKOS_INLINE_FUNCTION
156 Functor( outputValueViewType outputValues_,
157 inputPointViewType inputPoints_,
158 vinvViewType coeffs_,
160 : _outputValues(outputValues_), _inputPoints(inputPoints_),
161 _coeffs(coeffs_), _work(work_) {}
163 KOKKOS_INLINE_FUNCTION
164 void operator()(
const size_type iter)
const {
168 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
169 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
171 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
173 auto vcprop = Kokkos::common_view_alloc_prop(_work);
174 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
178 case OPERATOR_VALUE : {
179 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
184 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
189 INTREPID2_TEST_FOR_ABORT(
true,
190 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
199 template<
typename ExecSpaceType = void,
200 typename outputValueType = double,
201 typename pointValueType =
double>
203 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
212 const EPointType pointType = POINTTYPE_EQUISPACED);
226 const EOperator operatorType = OPERATOR_VALUE)
const {
227 #ifdef HAVE_INTREPID2_DEBUG
228 Intrepid2::getValues_HDIV_Args(outputValues,
235 Impl::Basis_HDIV_TRI_In_FEM::
236 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
245 #ifdef HAVE_INTREPID2_DEBUG
247 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
250 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
253 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
254 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
256 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
262 #ifdef HAVE_INTREPID2_DEBUG
264 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
265 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
267 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
268 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
270 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
271 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
273 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
277 getExpansionCoeffs( scalarViewType coeffs )
const {
279 Kokkos::deep_copy(coeffs, this->
coeffs_);
285 return "Intrepid2_HDIV_TRI_In_FEM";
298 Kokkos::DynRankView<scalarType,ExecSpaceType>
coeffs_;
virtual bool requireOrientation() const
True if orientation is required.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
virtual const char * getName() const
Returns basis name.
Basis_HDIV_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
Definition file for FEM basis functions of degree n for H(div) functions on TRI cells.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
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 getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.