Intrepid2
Intrepid2_HVOL_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 
15 #ifndef __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
16 #define __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
17 
18 #include "Intrepid2_Basis.hpp"
20 
21 namespace Intrepid2 {
22 
23  namespace Impl {
24 
29  public:
30  typedef struct Quadrilateral<4> cell_topology_type;
34  template<EOperator opType>
35  struct Serial {
36  template<typename outputValueViewType,
37  typename inputPointViewType,
38  typename workViewType,
39  typename vinvViewType>
40  KOKKOS_INLINE_FUNCTION
41  static void
42  getValues( outputValueViewType outputValues,
43  const inputPointViewType inputPoints,
44  workViewType work,
45  const vinvViewType vinv,
46  const ordinal_type operatorDn = 0 );
47 
48  KOKKOS_INLINE_FUNCTION
49  static ordinal_type
50  getWorkSizePerPoint(ordinal_type order) {return 3*getPnCardinality<1>(order); }
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( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
59  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
60  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
61  const EOperator operatorType);
62 
66  template<typename outputValueViewType,
67  typename inputPointViewType,
68  typename vinvViewType,
69  typename workViewType,
70  EOperator opType,
71  ordinal_type numPtsEval>
72  struct Functor {
73  outputValueViewType _outputValues;
74  const inputPointViewType _inputPoints;
75  const vinvViewType _vinv;
76  workViewType _work;
77  const ordinal_type _opDn;
78 
79  KOKKOS_INLINE_FUNCTION
80  Functor( outputValueViewType outputValues_,
81  inputPointViewType inputPoints_,
82  vinvViewType vinv_,
83  workViewType work_,
84  const ordinal_type opDn_ = 0 )
85  : _outputValues(outputValues_), _inputPoints(inputPoints_),
86  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
87 
88  KOKKOS_INLINE_FUNCTION
89  void operator()(const size_type iter) const {
90  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
91  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
92 
93  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
94  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
95 
96  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
97 
98  auto vcprop = Kokkos::common_view_alloc_prop(_work);
99  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
100 
101  switch (opType) {
102  case OPERATOR_VALUE : {
103  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
104  Serial<opType>::getValues( output, input, work, _vinv );
105  break;
106  }
107  case OPERATOR_Dn : {
108  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
109  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
110  break;
111  }
112  default: {
113  INTREPID2_TEST_FOR_ABORT( true,
114  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Functor) operator is not supported");
115 
116  }
117  }
118  }
119  };
120  };
121  }
122 
129  template<typename DeviceType = void,
130  typename outputValueType = double,
131  typename pointValueType = double>
133  : public Basis<DeviceType,outputValueType,pointValueType> {
134  public:
136 
137  using typename BasisBase::OrdinalTypeArray1DHost;
138  using typename BasisBase::OrdinalTypeArray2DHost;
139  using typename BasisBase::OrdinalTypeArray3DHost;
140 
141  using typename BasisBase::OutputViewType;
142  using typename BasisBase::PointViewType ;
143  using typename BasisBase::ScalarViewType;
144 
147  Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order,
148  const EPointType pointType = POINTTYPE_EQUISPACED);
149 
150  using BasisBase::getValues;
151 
152  virtual
153  void
154  getValues( OutputViewType outputValues,
155  const PointViewType inputPoints,
156  const EOperator operatorType = OPERATOR_VALUE ) const override {
157 #ifdef HAVE_INTREPID2_DEBUG
158  Intrepid2::getValues_HVOL_Args(outputValues,
159  inputPoints,
160  operatorType,
161  this->getBaseCellTopology(),
162  this->getCardinality() );
163 #endif
164  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
165  Impl::Basis_HVOL_QUAD_Cn_FEM::
166  getValues<DeviceType,numPtsPerEval>( outputValues,
167  inputPoints,
168  this->vinv_,
169  operatorType );
170  }
171 
172  virtual void
173  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
174  ordinal_type& perThreadSpaceSize,
175  const PointViewType inputPoints,
176  const EOperator operatorType = OPERATOR_VALUE) const override;
177 
178  KOKKOS_INLINE_FUNCTION
179  virtual void
180  getValues(
181  OutputViewType outputValues,
182  const PointViewType inputPoints,
183  const EOperator operatorType,
184  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
185  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
186  const ordinal_type subcellDim = -1,
187  const ordinal_type subcellOrdinal = -1) const override;
188 
189 
190  virtual
191  void
192  getDofCoords( ScalarViewType dofCoords ) const override {
193 #ifdef HAVE_INTREPID2_DEBUG
194  // Verify rank of output array.
195  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
196  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
197  // Verify 0th dimension of output array.
198  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
199  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
200  // Verify 1st dimension of output array.
201  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
202  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
203 #endif
204  Kokkos::deep_copy(dofCoords, this->dofCoords_);
205  }
206 
207  virtual
208  void
209  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
210 #ifdef HAVE_INTREPID2_DEBUG
211  // Verify rank of output array.
212  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
213  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
214  // Verify 0th dimension of output array.
215  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
216  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
217 #endif
218  Kokkos::deep_copy(dofCoeffs, 1.0);
219  }
220 
221  virtual
222  const char*
223  getName() const override {
224  return "Intrepid2_HVOL_QUAD_Cn_FEM";
225  }
226 
227  virtual
228  bool
229  requireOrientation() const override {
230  return false;
231  }
232 
233  virtual HostBasisPtr<outputValueType,pointValueType>
234  getHostBasis() const override{
236  }
237 
238  private:
239 
241  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
242  EPointType pointType_;
243  };
244 
245 }// namespace Intrepid2
246 
248 
249 #endif
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
small utility functions
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Definition file for FEM basis functions of degree n for H(vol) functions on QUAD. ...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, 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...
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.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements ...
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)
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
virtual const char * getName() const override
Returns basis name.
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.