Intrepid2
Intrepid2_HGRAD_QUAD_Cn_FEM.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
16 #ifndef __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
17 #define __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
24  namespace Impl {
25 
30  public:
31  typedef struct Quadrilateral<4> cell_topology_type;
32 
38  template<EOperator opType>
39  struct Serial {
40  template<typename outputValueViewType,
41  typename inputPointViewType,
42  typename workViewType,
43  typename vinvViewType>
44  KOKKOS_INLINE_FUNCTION
45  static void
46  getValues( outputValueViewType outputValues,
47  const inputPointViewType inputPoints,
48  workViewType work,
49  const vinvViewType vinv,
50  const ordinal_type operatorDn = 0 );
51  };
52 
53  template<typename DeviceType, ordinal_type numPtsPerEval,
54  typename outputValueValueType, class ...outputValueProperties,
55  typename inputPointValueType, class ...inputPointProperties,
56  typename vinvValueType, class ...vinvProperties>
57  static void
58  getValues( const typename DeviceType::execution_space& space,
59  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
60  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
61  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
62  const EOperator operatorType);
63 
64 
68  template<typename outputValueViewType,
69  typename inputPointViewType,
70  typename vinvViewType,
71  typename workViewType,
72  EOperator opType,
73  ordinal_type numPtsEval>
74  struct Functor {
75  outputValueViewType _outputValues;
76  const inputPointViewType _inputPoints;
77  const vinvViewType _vinv;
78  workViewType _work;
79  const ordinal_type _opDn;
80 
81  KOKKOS_INLINE_FUNCTION
82  Functor( outputValueViewType outputValues_,
83  inputPointViewType inputPoints_,
84  vinvViewType vinv_,
85  workViewType work_,
86  const ordinal_type opDn_ = 0 )
87  : _outputValues(outputValues_), _inputPoints(inputPoints_),
88  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
89 
90  KOKKOS_INLINE_FUNCTION
91  void operator()(const size_type iter) const {
92  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
93  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
94 
95  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
96  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
97 
98  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
99 
100  auto vcprop = Kokkos::common_view_alloc_prop(_work);
101  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
102 
103  switch (opType) {
104  case OPERATOR_VALUE : {
105  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
106  Serial<opType>::getValues( output, input, work, _vinv );
107  break;
108  }
109  case OPERATOR_CURL :
110  case OPERATOR_Dn : {
111  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
112  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
113  break;
114  }
115  default: {
116  INTREPID2_TEST_FOR_ABORT( true,
117  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
118 
119  }
120  }
121  }
122  };
123  };
124  }
125 
131  template<typename DeviceType = void,
132  typename outputValueType = double,
133  typename pointValueType = double>
135  : public Basis<DeviceType,outputValueType,pointValueType> {
136  public:
138  using typename BasisBase::ExecutionSpace;
139  using typename BasisBase::OrdinalTypeArray1DHost;
140  using typename BasisBase::OrdinalTypeArray2DHost;
141  using typename BasisBase::OrdinalTypeArray3DHost;
142 
143  using typename BasisBase::OutputViewType;
144  using typename BasisBase::PointViewType;
145  using typename BasisBase::ScalarViewType;
146 
147  private:
149  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
150 
152  EPointType pointType_;
153 
154  public:
157  Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order,
158  const EPointType pointType = POINTTYPE_EQUISPACED);
159 
160  using BasisBase::getValues;
161 
162  virtual
163  void
164  getValues( const ExecutionSpace& space,
165  OutputViewType outputValues,
166  const PointViewType inputPoints,
167  const EOperator operatorType = OPERATOR_VALUE ) const override {
168 #ifdef HAVE_INTREPID2_DEBUG
169  Intrepid2::getValues_HGRAD_Args(outputValues,
170  inputPoints,
171  operatorType,
172  this->getBaseCellTopology(),
173  this->getCardinality() );
174 #endif
175  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
176  Impl::Basis_HGRAD_QUAD_Cn_FEM::
177  getValues<DeviceType,numPtsPerEval>(space,
178  outputValues,
179  inputPoints,
180  this->vinv_,
181  operatorType);
182  }
183 
184  virtual void
185  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
186  ordinal_type& perThreadSpaceSize,
187  const PointViewType inputPointsconst,
188  const EOperator operatorType = OPERATOR_VALUE) const override;
189 
190  KOKKOS_INLINE_FUNCTION
191  virtual void
192  getValues(
193  OutputViewType outputValues,
194  const PointViewType inputPoints,
195  const EOperator operatorType,
196  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
197  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
198  const ordinal_type subcellDim = -1,
199  const ordinal_type subcellOrdinal = -1) const override;
200 
201  virtual
202  void
203  getDofCoords( ScalarViewType dofCoords ) const override {
204 #ifdef HAVE_INTREPID2_DEBUG
205  // Verify rank of output array.
206  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
207  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
208  // Verify 0th dimension of output array.
209  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
210  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
211  // Verify 1st dimension of output array.
212  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
213  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
214 #endif
215  Kokkos::deep_copy(dofCoords, this->dofCoords_);
216  }
217 
218  virtual
219  void
220  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
221 #ifdef HAVE_INTREPID2_DEBUG
222  // Verify rank of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
225  // Verify 0th dimension of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
228 #endif
229  Kokkos::deep_copy(dofCoeffs, 1.0);
230  }
231 
232  virtual
233  const char*
234  getName() const override {
235  return "Intrepid2_HGRAD_QUAD_Cn_FEM";
236  }
237 
238  virtual
239  bool
240  requireOrientation() const override {
241  return (this->basisDegree_ > 2);
242  }
243 
244  Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
245  getVandermondeInverse() const {
246  return vinv_;
247  }
248 
249  ordinal_type
250  getWorkSizePerPoint(const EOperator operatorType) const {
251  return 3*getPnCardinality<1>(this->basisDegree_);
252  }
253 
262  BasisPtr<DeviceType,outputValueType,pointValueType>
263  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
264  if(subCellDim == 1) {
265  return Teuchos::rcp(new
267  (this->basisDegree_,pointType_));
268  }
269  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
270  }
271 
272  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
273  getHostBasis() const override{
275  }
276  };
277 
278 }// namespace Intrepid2
279 
281 
282 #endif
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
small utility functions
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...
ordinal_type getCardinality() const
Returns cardinality of the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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...
Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual const char * getName() const override
Returns basis name.
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 (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree n for H(grad) functions on QUAD cells...
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM work is a rank 1 view having the same value_type of inputPoint...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.