48 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_HPP__ 49 #define __INTREPID2_HVOL_LINE_CN_FEM_HPP__ 55 #include "Teuchos_LAPACK.hpp" 85 template<EOperator opType>
87 template<
typename outputValueViewType,
88 typename inputPointViewType,
89 typename workViewType,
90 typename vinvViewType>
91 KOKKOS_INLINE_FUNCTION
93 getValues( outputValueViewType outputValues,
94 const inputPointViewType inputPoints,
96 const vinvViewType vinv,
97 const ordinal_type operatorDn = 0 );
99 KOKKOS_INLINE_FUNCTION
101 getWorkSizePerPoint(ordinal_type order) {
return getPnCardinality<1>(order);}
104 template<
typename DeviceType, ordinal_type numPtsPerEval,
105 typename outputValueValueType,
class ...outputValueProperties,
106 typename inputPointValueType,
class ...inputPointProperties,
107 typename vinvValueType,
class ...vinvProperties>
109 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
110 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
111 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
117 template<
typename outputValueViewType,
118 typename inputPointViewType,
119 typename vinvViewType,
120 typename workViewType,
122 ordinal_type numPtsEval>
124 outputValueViewType _outputValues;
125 const inputPointViewType _inputPoints;
126 const vinvViewType _vinv;
128 const ordinal_type _opDn;
130 KOKKOS_INLINE_FUNCTION
131 Functor( outputValueViewType outputValues_,
132 inputPointViewType inputPoints_,
135 const ordinal_type opDn_ = 0 )
136 : _outputValues(outputValues_), _inputPoints(inputPoints_),
137 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
139 KOKKOS_INLINE_FUNCTION
140 void operator()(
const size_type iter)
const {
144 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
145 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
147 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
149 auto vcprop = Kokkos::common_view_alloc_prop(_work);
150 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
153 case OPERATOR_VALUE : {
154 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
159 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
164 INTREPID2_TEST_FOR_ABORT(
true,
165 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Functor) operator is not supported");
174 template<
typename DeviceType = void,
175 typename outputValueType = double,
176 typename pointValueType =
double>
178 :
public Basis<DeviceType,outputValueType,pointValueType> {
194 const EPointType pointType = POINTTYPE_EQUISPACED);
200 getValues( OutputViewType outputValues,
201 const PointViewType inputPoints,
202 const EOperator operatorType = OPERATOR_VALUE )
const override {
203 #ifdef HAVE_INTREPID2_DEBUG 210 constexpr ordinal_type numPtsPerEval = 1;
211 Impl::Basis_HVOL_LINE_Cn_FEM::
212 getValues<DeviceType,numPtsPerEval>( outputValues,
220 getDofCoords( ScalarViewType dofCoords )
const override {
221 #ifdef HAVE_INTREPID2_DEBUG 223 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
226 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
229 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
232 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
237 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
238 #ifdef HAVE_INTREPID2_DEBUG 240 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
246 Kokkos::deep_copy(dofCoeffs, 1.0);
250 getVandermondeInverse( ScalarViewType vinv )
const {
252 Kokkos::deep_copy(vinv, this->
vinv_);
258 return "Intrepid2_HVOL_LINE_Cn_FEM";
270 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type), allowing host access to input and output views, and ensuring host execution of basis evaluation.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
virtual const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
void getValues_HVOL_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 HVOL-conforming FEM basis...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
EPointType
Enumeration of types of point distributions in Intrepid.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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 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...
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.
Header file for the abstract base class Intrepid2::Basis.