49 #ifndef __INTREPID2_HGRAD_LINE_CN_FEM_HPP__ 50 #define __INTREPID2_HGRAD_LINE_CN_FEM_HPP__ 56 #include "Teuchos_LAPACK.hpp" 87 template<EOperator opType>
89 template<
typename outputValueViewType,
90 typename inputPointViewType,
91 typename workViewType,
92 typename vinvViewType>
93 KOKKOS_INLINE_FUNCTION
95 getValues( outputValueViewType outputValues,
96 const inputPointViewType inputPoints,
98 const vinvViewType vinv,
99 const ordinal_type operatorDn = 0 );
102 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
103 typename outputValueValueType,
class ...outputValueProperties,
104 typename inputPointValueType,
class ...inputPointProperties,
105 typename vinvValueType,
class ...vinvProperties>
107 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
108 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
109 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
115 template<
typename outputValueViewType,
116 typename inputPointViewType,
117 typename vinvViewType,
118 typename workViewType,
120 ordinal_type numPtsEval>
122 outputValueViewType _outputValues;
123 const inputPointViewType _inputPoints;
124 const vinvViewType _vinv;
126 const ordinal_type _opDn;
128 KOKKOS_INLINE_FUNCTION
129 Functor( outputValueViewType outputValues_,
130 inputPointViewType inputPoints_,
133 const ordinal_type opDn_ = 0 )
134 : _outputValues(outputValues_), _inputPoints(inputPoints_),
135 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
137 KOKKOS_INLINE_FUNCTION
138 void operator()(
const size_type iter)
const {
142 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
143 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
145 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
147 auto vcprop = Kokkos::common_view_alloc_prop(_work);
148 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
151 case OPERATOR_VALUE : {
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
157 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
162 INTREPID2_TEST_FOR_ABORT(
true,
163 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Functor) operator is not supported");
172 template<
typename ExecSpaceType = void,
173 typename outputValueType = double,
174 typename pointValueType =
double>
176 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
190 Kokkos::DynRankView<typename scalarViewType::value_type,ExecSpaceType>
vinv_;
196 const EPointType pointType = POINTTYPE_EQUISPACED);
204 const EOperator operatorType = OPERATOR_VALUE )
const {
205 #ifdef HAVE_INTREPID2_DEBUG 212 constexpr ordinal_type numPtsPerEval = 1;
213 Impl::Basis_HGRAD_LINE_Cn_FEM::
214 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
223 #ifdef HAVE_INTREPID2_DEBUG 225 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
228 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
231 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
234 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
240 #ifdef HAVE_INTREPID2_DEBUG 242 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
245 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
248 Kokkos::deep_copy(dofCoeffs, 1.0);
254 return "Intrepid2_HGRAD_LINE_Cn_FEM";
258 getVandermondeInverse( scalarViewType vinv )
const {
260 Kokkos::deep_copy(vinv, this->vinv_);
263 Kokkos::DynRankView<typename scalarViewType::const_value_type,ExecSpaceType>
264 getVandermondeInverse()
const {
269 getWorkSizePerPoint(
const EOperator operatorType)
const {
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
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< typename scalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
ordinal_type getCardinality() const
Returns cardinality of the basis.
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.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
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 ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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...
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
virtual const char * getName() const
Returns basis name.
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...
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
EPointType
Enumeration of types of point distributions in Intrepid.
Definition file for FEM basis functions of degree n for H(grad) functions on LINE.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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.