49 #ifndef __INTREPID2_HGRAD_LINE_CN_FEM_JACOBI_HPP__
50 #define __INTREPID2_HGRAD_LINE_CN_FEM_JACOBI_HPP__
142 template<EOperator opType>
144 template<
typename OutputViewType,
145 typename inputViewType>
146 KOKKOS_INLINE_FUNCTION
148 getValues( OutputViewType output,
149 const inputViewType input,
150 const ordinal_type order,
153 const ordinal_type opDn = 0 );
156 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
157 typename outputValueValueType,
class ...outputValueProperties,
158 typename inputPointValueType,
class ...inputPointProperties>
160 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
161 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
162 const ordinal_type order,
165 const EOperator operatorType );
170 template<
typename outputValueViewType,
171 typename inputPointViewType,
173 ordinal_type numPtsEval>
175 outputValueViewType _outputValues;
176 const inputPointViewType _inputPoints;
177 const ordinal_type _order;
178 const double _alpha, _beta;
179 const ordinal_type _opDn;
181 KOKKOS_INLINE_FUNCTION
182 Functor( outputValueViewType outputValues_,
183 inputPointViewType inputPoints_,
184 const ordinal_type order_,
187 const ordinal_type opDn_ = 0 )
188 : _outputValues(outputValues_), _inputPoints(inputPoints_),
189 _order(order_), _alpha(alpha_), _beta(beta_), _opDn(opDn_) {}
191 KOKKOS_INLINE_FUNCTION
192 void operator()(
const size_type iter)
const {
196 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
197 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
199 if (input.extent(0)) {
201 case OPERATOR_VALUE : {
202 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
206 case OPERATOR_GRAD : {
207 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
212 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
217 INTREPID2_TEST_FOR_ABORT(
true,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI::Functor) operator is not supported");
228 template<
typename ExecSpaceType = void,
229 typename outputValueType = double,
230 typename pointValueType =
double>
232 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
234 typedef double value_type;
242 const double alpha = 0,
243 const double beta = 0 );
255 const EOperator operatorType = OPERATOR_VALUE )
const {
256 #ifdef HAVE_INTREPID2_DEBUG
257 Intrepid2::getValues_HGRAD_Args(outputValues,
263 constexpr ordinal_type numPtsPerEval = 1;
264 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
265 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
274 double alpha_, beta_;
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI.
Basis_HGRAD_LINE_Cn_FEM_JACOBI(const ordinal_type order, const double alpha=0, const double beta=0)
Constructor.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type getDegree() const
Returns the degree of the basis.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI.
Definition file for FEM orthogonal basis functions of degree n for H(grad) functions on LINE...
Header file for the abstract base class Intrepid2::Basis.