48 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_LINE_CN_FEM_HPP__
55 #include "Teuchos_LAPACK.hpp"
81 typedef struct Line<2> cell_topology_type;
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 ExecSpaceType, 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,
112 const EOperator operatorType );
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 ExecSpaceType = void,
175 typename outputValueType = double,
176 typename pointValueType =
double>
178 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
187 const EPointType pointType = POINTTYPE_EQUISPACED);
199 const EOperator operatorType = OPERATOR_VALUE )
const {
200 #ifdef HAVE_INTREPID2_DEBUG
201 Intrepid2::getValues_HVOL_Args(outputValues,
207 constexpr ordinal_type numPtsPerEval = 1;
208 Impl::Basis_HVOL_LINE_Cn_FEM::
209 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
218 #ifdef HAVE_INTREPID2_DEBUG
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
229 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
235 #ifdef HAVE_INTREPID2_DEBUG
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
240 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
243 Kokkos::deep_copy(dofCoeffs, 1.0);
247 getVandermondeInverse( scalarViewType vinv )
const {
249 Kokkos::deep_copy(vinv, this->
vinv_);
255 return "Intrepid2_HVOL_LINE_Cn_FEM";
262 Kokkos::DynRankView<typename scalarViewType::value_type,ExecSpaceType>
vinv_;
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
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.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< typename scalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.