49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__ 50 #define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__ 57 #include "Teuchos_LAPACK.hpp" 96 template<EOperator opType>
98 template<
typename outputValueViewType,
99 typename inputPointViewType,
100 typename workViewType,
101 typename vinvViewType>
102 KOKKOS_INLINE_FUNCTION
104 getValues( outputValueViewType outputValues,
105 const inputPointViewType inputPoints,
107 const vinvViewType vinv );
110 template<
typename DeviceType, ordinal_type numPtsPerEval,
111 typename outputValueValueType,
class ...outputValueProperties,
112 typename inputPointValueType,
class ...inputPointProperties,
113 typename vinvValueType,
class ...vinvProperties>
115 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
116 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
117 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
123 template<
typename outputValueViewType,
124 typename inputPointViewType,
125 typename vinvViewType,
126 typename workViewType,
128 ordinal_type numPtsEval>
130 outputValueViewType _outputValues;
131 const inputPointViewType _inputPoints;
132 const vinvViewType _vinv;
135 KOKKOS_INLINE_FUNCTION
136 Functor( outputValueViewType outputValues_,
137 inputPointViewType inputPoints_,
140 : _outputValues(outputValues_), _inputPoints(inputPoints_),
141 _vinv(vinv_), _work(work_) {}
143 KOKKOS_INLINE_FUNCTION
144 void operator()(
const size_type iter)
const {
148 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
149 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
151 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
153 auto vcprop = Kokkos::common_view_alloc_prop(_work);
154 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
157 case OPERATOR_VALUE : {
158 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
165 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
170 INTREPID2_TEST_FOR_ABORT(
true,
171 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
180 template<
typename DeviceType = void,
181 typename outputValueType = double,
182 typename pointValueType =
double>
184 :
public Basis<DeviceType,outputValueType,pointValueType> {
200 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
209 const EPointType pointType = POINTTYPE_EQUISPACED);
215 getValues( OutputViewType outputValues,
216 const PointViewType inputPoints,
217 const EOperator operatorType = OPERATOR_VALUE)
const override {
218 #ifdef HAVE_INTREPID2_DEBUG 226 Impl::Basis_HGRAD_TRI_Cn_FEM::
227 getValues<DeviceType,numPtsPerEval>( outputValues,
235 getDofCoords( ScalarViewType dofCoords )
const override {
236 #ifdef HAVE_INTREPID2_DEBUG 238 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
241 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
247 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
252 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
253 #ifdef HAVE_INTREPID2_DEBUG 255 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
258 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
259 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
261 Kokkos::deep_copy(dofCoeffs, 1.0);
268 return "Intrepid2_HGRAD_TRI_Cn_FEM";
278 getVandermondeInverse( ScalarViewType vinv )
const {
280 Kokkos::deep_copy(vinv, this->vinv_);
283 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
284 getVandermondeInverse()
const {
289 getWorkSizePerPoint(
const EOperator operatorType)
const {
290 auto cardinality = getPnCardinality<2>(this->
basisDegree_);
291 switch (operatorType) {
295 return 5*cardinality;
309 BasisPtr<DeviceType,outputValueType,pointValueType>
311 if(subCellDim == 1) {
312 return Teuchos::rcp(
new 316 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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.
virtual const char * getName() const override
Returns basis name.
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 (...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns multiplicities of dx, dy, and dz based on the enumeration of the partial derivative, its order and the space dimension. Inverse of the getDkEnumeration() method.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EPointType pointType_
type of lattice used for creating the DoF coordinates
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
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...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
EPointType
Enumeration of types of point distributions in Intrepid.
Triangle topology, 3 nodes.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Definition file for FEM basis functions of degree n for H(grad) functions on TRI cells.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Header file for the abstract base class Intrepid2::Basis.