15 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_HPP__
16 #define __INTREPID2_HVOL_LINE_CN_FEM_HPP__
22 #include "Teuchos_LAPACK.hpp"
48 typedef struct Line<2> cell_topology_type;
52 template<EOperator opType>
54 template<
typename outputValueViewType,
55 typename inputPointViewType,
56 typename workViewType,
57 typename vinvViewType>
58 KOKKOS_INLINE_FUNCTION
60 getValues( outputValueViewType outputValues,
61 const inputPointViewType inputPoints,
63 const vinvViewType vinv,
64 const ordinal_type operatorDn = 0 );
66 KOKKOS_INLINE_FUNCTION
68 getWorkSizePerPoint(ordinal_type order) {
return getPnCardinality<1>(order);}
71 template<
typename DeviceType, ordinal_type numPtsPerEval,
72 typename outputValueValueType,
class ...outputValueProperties,
73 typename inputPointValueType,
class ...inputPointProperties,
74 typename vinvValueType,
class ...vinvProperties>
76 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
77 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
78 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
79 const EOperator operatorType );
84 template<
typename outputValueViewType,
85 typename inputPointViewType,
86 typename vinvViewType,
87 typename workViewType,
89 ordinal_type numPtsEval>
91 outputValueViewType _outputValues;
92 const inputPointViewType _inputPoints;
93 const vinvViewType _vinv;
95 const ordinal_type _opDn;
97 KOKKOS_INLINE_FUNCTION
98 Functor( outputValueViewType outputValues_,
99 inputPointViewType inputPoints_,
102 const ordinal_type opDn_ = 0 )
103 : _outputValues(outputValues_), _inputPoints(inputPoints_),
104 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
106 KOKKOS_INLINE_FUNCTION
107 void operator()(
const size_type iter)
const {
111 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
112 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
114 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
116 auto vcprop = Kokkos::common_view_alloc_prop(_work);
117 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
120 case OPERATOR_VALUE : {
121 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
126 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
131 INTREPID2_TEST_FOR_ABORT(
true,
132 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Functor) operator is not supported");
141 template<
typename DeviceType = void,
142 typename outputValueType = double,
143 typename pointValueType =
double>
145 :
public Basis<DeviceType,outputValueType,pointValueType> {
162 const EPointType pointType = POINTTYPE_EQUISPACED);
170 const EOperator operatorType = OPERATOR_VALUE )
const override {
171 #ifdef HAVE_INTREPID2_DEBUG
172 Intrepid2::getValues_HVOL_Args(outputValues,
178 constexpr ordinal_type numPtsPerEval = 1;
179 Impl::Basis_HVOL_LINE_Cn_FEM::
180 getValues<DeviceType,numPtsPerEval>( outputValues,
188 ordinal_type& perThreadSpaceSize,
190 const EOperator operatorType = OPERATOR_VALUE)
const override;
192 KOKKOS_INLINE_FUNCTION
197 const EOperator operatorType,
198 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
199 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
200 const ordinal_type subcellDim = -1,
201 const ordinal_type subcellOrdinal = -1)
const override;
206 #ifdef HAVE_INTREPID2_DEBUG
208 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
211 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
217 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
223 #ifdef HAVE_INTREPID2_DEBUG
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
228 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
231 Kokkos::deep_copy(dofCoeffs, 1.0);
237 Kokkos::deep_copy(vinv, this->
vinv_);
243 return "Intrepid2_HVOL_LINE_Cn_FEM";
246 virtual HostBasisPtr<outputValueType,pointValueType>
255 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
256 EPointType pointType_;
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
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::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
virtual const char * getName() const override
Returns basis name.
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...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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...
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.