16#ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
24#include "Teuchos_LAPACK.hpp"
63 template<EOperator OpType>
65 template<
typename OutputValueViewType,
66 typename InputPointViewType,
67 typename WorkViewType,
68 typename VinvViewType>
69 KOKKOS_INLINE_FUNCTION
71 getValues( OutputValueViewType outputValues,
72 const InputPointViewType inputPoints,
74 const VinvViewType vinv,
75 const ordinal_type order);
78 template<
typename DeviceType, ordinal_type numPtsPerEval,
79 typename OutputValueValueType,
class ...OutputValueProperties,
80 typename InputPointValueType,
class ...InputPointProperties,
81 typename VinvValueType,
class ...VinvProperties>
83 getValues(
const typename DeviceType::execution_space& space,
84 Kokkos::DynRankView<OutputValueValueType,OutputValueProperties...> outputValues,
85 const Kokkos::DynRankView<InputPointValueType, InputPointProperties...> inputPoints,
86 const Kokkos::DynRankView<VinvValueType, VinvProperties...> vinv,
87 const ordinal_type order,
88 const EOperator operatorType);
93 template<
typename OutputValueViewType,
94 typename InputPointViewType,
95 typename VinvViewType,
96 typename WorkViewType,
98 ordinal_type numPtsEval>
100 OutputValueViewType _outputValues;
101 const InputPointViewType _inputPoints;
102 const VinvViewType _vinv;
104 const ordinal_type _order;
106 KOKKOS_INLINE_FUNCTION
107 Functor( OutputValueViewType outputValues_,
108 InputPointViewType inputPoints_,
112 : _outputValues(outputValues_), _inputPoints(inputPoints_),
113 _vinv(vinv_), _work(work_), _order(order_) {}
115 KOKKOS_INLINE_FUNCTION
116 void operator()(
const size_type iter)
const {
120 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
121 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
123 typename WorkViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
125 auto vcprop = Kokkos::common_view_alloc_prop(_work);
126 WorkViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
129 case OPERATOR_VALUE : {
130 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
137 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
142 INTREPID2_TEST_FOR_ABORT(
true,
143 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
152 template<
typename DeviceType = void,
153 typename outputValueType = double,
154 typename pointValueType =
double>
156 :
public Basis<DeviceType,outputValueType,pointValueType> {
175 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
184 const EPointType pointType = POINTTYPE_EQUISPACED);
193 const EOperator operatorType = OPERATOR_VALUE)
const override {
194#ifdef HAVE_INTREPID2_DEBUG
202 Impl::Basis_HGRAD_TRI_Cn_FEM::
203 getValues<DeviceType,numPtsPerEval>(space,
213 ordinal_type& perThreadSpaceSize,
215 const EOperator operatorType = OPERATOR_VALUE)
const override;
217 KOKKOS_INLINE_FUNCTION
222 const EOperator operatorType,
223 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
224 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
225 const ordinal_type subcellDim = -1,
226 const ordinal_type subcellOrdinal = -1)
const override;
233#ifdef HAVE_INTREPID2_DEBUG
235 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
238 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
241 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
244 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
250#ifdef HAVE_INTREPID2_DEBUG
252 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
253 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
255 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
258 Kokkos::deep_copy(dofCoeffs, 1.0);
265 return "Intrepid2_HGRAD_TRI_Cn_FEM";
277 Kokkos::deep_copy(vinv, this->vinv_);
280 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
281 getVandermondeInverse()
const {
286 getWorkSizePerPoint(
const EOperator operatorType)
const {
287 auto cardinality = getPnCardinality<2>(this->
basisDegree_);
288 switch (operatorType) {
292 return 5*cardinality;
306 BasisPtr<DeviceType,outputValueType,pointValueType>
308 if(subCellDim == 1) {
309 return Teuchos::rcp(
new
313 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(grad) functions on TRI cells.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual bool requireOrientation() const override
True if orientation is required.
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.
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.
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.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Triangle topology, 3 nodes.