49#ifndef __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
50#define __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
71 template<EOperator opType>
73 template<
typename outputValueViewType,
74 typename inputPointViewType,
75 typename workViewType,
76 typename vinvViewType>
77 KOKKOS_INLINE_FUNCTION
79 getValues( outputValueViewType outputValues,
80 const inputPointViewType inputPoints,
82 const vinvViewType vinv,
83 const ordinal_type operatorDn = 0 );
86 template<
typename DeviceType, ordinal_type numPtsPerEval,
87 typename outputValueValueType,
class ...outputValueProperties,
88 typename inputPointValueType,
class ...inputPointProperties,
89 typename vinvValueType,
class ...vinvProperties>
91 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
94 const EOperator operatorType);
100 template<
typename outputValueViewType,
101 typename inputPointViewType,
102 typename vinvViewType,
103 typename workViewType,
105 ordinal_type numPtsEval>
107 outputValueViewType _outputValues;
108 const inputPointViewType _inputPoints;
109 const vinvViewType _vinv;
111 const ordinal_type _opDn;
113 KOKKOS_INLINE_FUNCTION
114 Functor( outputValueViewType outputValues_,
115 inputPointViewType inputPoints_,
118 const ordinal_type opDn_ = 0 )
119 : _outputValues(outputValues_), _inputPoints(inputPoints_),
120 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
122 KOKKOS_INLINE_FUNCTION
123 void operator()(
const size_type iter)
const {
127 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
128 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
130 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
132 auto vcprop = Kokkos::common_view_alloc_prop(_work);
133 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
136 case OPERATOR_VALUE : {
137 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
143 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
148 INTREPID2_TEST_FOR_ABORT(
true,
149 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
163 template<
typename DeviceType = void,
164 typename outputValueType = double,
165 typename pointValueType =
double>
167 :
public Basis<DeviceType,outputValueType,pointValueType> {
179 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
188 const EPointType pointType = POINTTYPE_EQUISPACED);
195 const PointViewType inputPoints,
196 const EOperator operatorType = OPERATOR_VALUE )
const override {
197#ifdef HAVE_INTREPID2_DEBUG
205 Impl::Basis_HGRAD_QUAD_Cn_FEM::
206 getValues<DeviceType,numPtsPerEval>( outputValues,
215#ifdef HAVE_INTREPID2_DEBUG
217 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
232#ifdef HAVE_INTREPID2_DEBUG
234 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
237 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
240 Kokkos::deep_copy(dofCoeffs, 1.0);
246 return "Intrepid2_HGRAD_QUAD_Cn_FEM";
255 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
256 getVandermondeInverse()
const {
261 getWorkSizePerPoint(
const EOperator operatorType)
const {
273 BasisPtr<DeviceType,outputValueType,pointValueType>
275 if(subCellDim == 1) {
276 return Teuchos::rcp(
new
280 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(grad) functions on QUAD cells.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
BasisPtr< typename Kokkos::HostSpace::device_type, 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 const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
EPointType pointType_
type of lattice used for creating the DoF coordinates
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
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
View type for 2d host array.
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 ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM work is a rank 1 view having the same value_type of inputPoint...
Quadrilateral topology, 4 nodes.