49#ifndef __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
50#define __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
56#include "Teuchos_LAPACK.hpp"
83 typedef struct Line<2> cell_topology_type;
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 DeviceType, 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,
110 const EOperator operatorType );
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 {
139 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
140 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
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 );
153 Serial<opType>::getValues( output, input, work, _vinv );
157 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
158 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
162 INTREPID2_TEST_FOR_ABORT(
true,
163 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Functor) operator is not supported");
172 template<
typename DeviceType = void,
173 typename outputValueType = double,
174 typename pointValueType =
double>
176 :
public Basis<DeviceType,outputValueType,pointValueType> {
178 using BasisBase = Basis<DeviceType,outputValueType,pointValueType>;
194 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
195 EPointType pointType_;
200 const EPointType pointType = POINTTYPE_EQUISPACED);
206 getValues( OutputViewType outputValues,
207 const PointViewType inputPoints,
208 const EOperator operatorType = OPERATOR_VALUE )
const override {
209#ifdef HAVE_INTREPID2_DEBUG
216 constexpr ordinal_type numPtsPerEval = 1;
217 Impl::Basis_HGRAD_LINE_Cn_FEM::
218 getValues<DeviceType,numPtsPerEval>( outputValues,
227#ifdef HAVE_INTREPID2_DEBUG
229 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
232 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
235 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
238 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
244#ifdef HAVE_INTREPID2_DEBUG
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
249 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
252 Kokkos::deep_copy(dofCoeffs, 1.0);
258 return "Intrepid2_HGRAD_LINE_Cn_FEM";
264 Kokkos::deep_copy(vinv, this->vinv_);
267 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
268 getVandermondeInverse()
const {
273 getWorkSizePerPoint(
const EOperator operatorType)
const {
281 virtual HostBasisPtr<outputValueType,pointValueType>
283 auto hostBasis = Teuchos::rcp(
new HostBasis(this->
basisDegree_, pointType_));
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.
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 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, typename Kokkos::HostSpace::device_type > vinv_
Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
ordinal_type basisDegree_
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
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
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.