16#ifndef __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
23#include "Teuchos_LAPACK.hpp"
54 template<EOperator opType>
56 template<
typename outputValueViewType,
57 typename inputPointViewType,
58 typename workViewType,
59 typename vinvViewType>
60 KOKKOS_INLINE_FUNCTION
62 getValues( outputValueViewType outputValues,
63 const inputPointViewType inputPoints,
65 const vinvViewType vinv,
66 const ordinal_type operatorDn = 0 );
69 template<
typename DeviceType, ordinal_type numPtsPerEval,
70 typename outputValueValueType,
class ...outputValueProperties,
71 typename inputPointValueType,
class ...inputPointProperties,
72 typename vinvValueType,
class ...vinvProperties>
74 getValues(
const typename DeviceType::execution_space& space,
75 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
76 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
77 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
78 const EOperator operatorType );
83 template<
typename outputValueViewType,
84 typename inputPointViewType,
85 typename vinvViewType,
86 typename workViewType,
88 ordinal_type numPtsEval>
90 outputValueViewType _outputValues;
91 const inputPointViewType _inputPoints;
92 const vinvViewType _vinv;
94 const ordinal_type _opDn;
96 KOKKOS_INLINE_FUNCTION
97 Functor( outputValueViewType outputValues_,
98 inputPointViewType inputPoints_,
101 const ordinal_type opDn_ = 0 )
102 : _outputValues(outputValues_), _inputPoints(inputPoints_),
103 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
105 KOKKOS_INLINE_FUNCTION
106 void operator()(
const size_type iter)
const {
110 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
111 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
113 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
115 auto vcprop = Kokkos::common_view_alloc_prop(_work);
116 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
119 case OPERATOR_VALUE : {
120 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
125 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
130 INTREPID2_TEST_FOR_ABORT(
true,
131 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Functor) operator is not supported");
140 template<
typename DeviceType = void,
141 typename outputValueType = double,
142 typename pointValueType =
double>
144 :
public Basis<DeviceType,outputValueType,pointValueType> {
163 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
164 EPointType pointType_;
169 const EPointType pointType = POINTTYPE_EQUISPACED);
178 const EOperator operatorType = OPERATOR_VALUE )
const override {
179#ifdef HAVE_INTREPID2_DEBUG
186 constexpr ordinal_type numPtsPerEval = 1;
187 Impl::Basis_HGRAD_LINE_Cn_FEM::
188 getValues<DeviceType,numPtsPerEval>(space,
197 ordinal_type& perThreadSpaceSize,
199 const EOperator operatorType = OPERATOR_VALUE)
const override;
201 KOKKOS_INLINE_FUNCTION
206 const EOperator operatorType,
207 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
208 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
209 const ordinal_type subcellDim = -1,
210 const ordinal_type subcellOrdinal = -1)
const override;
215#ifdef HAVE_INTREPID2_DEBUG
217 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
232#ifdef HAVE_INTREPID2_DEBUG
234 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
237 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
240 Kokkos::deep_copy(dofCoeffs, 1.0);
246 return "Intrepid2_HGRAD_LINE_Cn_FEM";
252 Kokkos::deep_copy(vinv, this->vinv_);
255 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
256 getVandermondeInverse()
const {
261 getWorkSizePerPoint(
const EOperator operatorType)
const {
269 virtual HostBasisPtr<outputValueType,pointValueType>
Header file for the abstract base class Intrepid2::Basis.
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....
Definition file for FEM basis functions of degree n for H(grad) functions on LINE.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
virtual const char * getName() const override
Returns basis name.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the 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< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, 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...
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.
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_LINE_Cn_FEM.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.