48 #ifndef __INTREPID2_HVOL_TRI_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_TRI_CN_FEM_HPP__
54 #include "Teuchos_LAPACK.hpp"
79 typedef struct Triangle<3> cell_topology_type;
83 template<EOperator opType>
85 template<
typename outputValueViewType,
86 typename inputPointViewType,
87 typename workViewType,
88 typename vinvViewType>
89 KOKKOS_INLINE_FUNCTION
91 getValues( outputValueViewType outputValues,
92 const inputPointViewType inputPoints,
94 const vinvViewType vinv );
97 KOKKOS_INLINE_FUNCTION
99 getWorkSizePerPoint(ordinal_type order) {
100 auto cardinality = getPnCardinality<2>(order);
105 return 5*cardinality;
107 return getDkCardinality<opType,2>()*cardinality;
112 template<
typename DeviceType, ordinal_type numPtsPerEval,
113 typename outputValueValueType,
class ...outputValueProperties,
114 typename inputPointValueType,
class ...inputPointProperties,
115 typename vinvValueType,
class ...vinvProperties>
117 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
118 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
119 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
120 const EOperator operatorType);
125 template<
typename outputValueViewType,
126 typename inputPointViewType,
127 typename vinvViewType,
128 typename workViewType,
130 ordinal_type numPtsEval>
132 outputValueViewType _outputValues;
133 const inputPointViewType _inputPoints;
134 const vinvViewType _vinv;
137 KOKKOS_INLINE_FUNCTION
138 Functor( outputValueViewType outputValues_,
139 inputPointViewType inputPoints_,
142 : _outputValues(outputValues_), _inputPoints(inputPoints_),
143 _vinv(vinv_), _work(work_) {}
145 KOKKOS_INLINE_FUNCTION
146 void operator()(
const size_type iter)
const {
150 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
151 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
153 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
155 auto vcprop = Kokkos::common_view_alloc_prop(_work);
156 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
159 case OPERATOR_VALUE : {
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
166 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
171 INTREPID2_TEST_FOR_ABORT(
true,
172 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::Functor) operator is not supported");
181 template<
typename DeviceType = void,
182 typename outputValueType = double,
183 typename pointValueType =
double>
185 :
public Basis<DeviceType,outputValueType,pointValueType> {
196 const EPointType pointType = POINTTYPE_EQUISPACED);
209 getValues( OutputViewType outputValues,
210 const PointViewType inputPoints,
211 const EOperator operatorType = OPERATOR_VALUE)
const override {
212 #ifdef HAVE_INTREPID2_DEBUG
213 Intrepid2::getValues_HVOL_Args(outputValues,
220 Impl::Basis_HVOL_TRI_Cn_FEM::
221 getValues<DeviceType,numPtsPerEval>( outputValues,
229 getDofCoords( ScalarViewType dofCoords )
const override {
230 #ifdef HAVE_INTREPID2_DEBUG
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
235 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
241 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
246 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
247 #ifdef HAVE_INTREPID2_DEBUG
249 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
252 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
253 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
255 Kokkos::deep_copy(dofCoeffs, 1.0);
259 getVandermondeInverse( ScalarViewType vinv )
const {
261 Kokkos::deep_copy(vinv, this->
vinv_);
267 return "Intrepid2_HVOL_TRI_Cn_FEM";
276 virtual HostBasisPtr<outputValueType,pointValueType>
285 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
286 EPointType pointType_;
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
ordinal_type getCardinality() const
Returns cardinality of the basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
virtual const char * getName() const override
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Definition file for FEM basis functions of degree n for H(vol) functions on TRI.
Header file for the abstract base class Intrepid2::Basis.
virtual bool requireOrientation() const override
True if orientation is required.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.