Intrepid2
Intrepid2_HCURL_QUAD_In_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_HCURL_QUAD_IN_FEM_HPP__
17 #define __INTREPID2_HCURL_QUAD_IN_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
22 
23 namespace Intrepid2 {
24 
25  namespace Impl {
26 
31  public:
32  typedef struct Quadrilateral<4> cell_topology_type;
33 
37  template<EOperator opType>
38  struct Serial {
39  template<typename outputValueViewType,
40  typename inputPointViewType,
41  typename workViewType,
42  typename vinvViewType>
43  KOKKOS_INLINE_FUNCTION
44  static void
45  getValues( outputValueViewType outputValues,
46  const inputPointViewType inputPoints,
47  workViewType work,
48  const vinvViewType vinvLine,
49  const vinvViewType vinvBubble );
50 
51  KOKKOS_INLINE_FUNCTION
52  static ordinal_type
53  getWorkSizePerPoint(ordinal_type order) {
54  return 2*getPnCardinality<1>(order)+getPnCardinality<1>(order-1);
55  }
56  };
57 
58  template<typename DeviceType, ordinal_type numPtsPerEval,
59  typename outputValueValueType, class ...outputValueProperties,
60  typename inputPointValueType, class ...inputPointProperties,
61  typename vinvValueType, class ...vinvProperties>
62  static void
63  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
64  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
65  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
66  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
67  const EOperator operatorType);
68 
72  template<typename outputValueViewType,
73  typename inputPointViewType,
74  typename vinvViewType,
75  typename workViewType,
76  EOperator opType,
77  ordinal_type numPtsEval>
78  struct Functor {
79  outputValueViewType _outputValues;
80  const inputPointViewType _inputPoints;
81  const vinvViewType _vinvLine;
82  const vinvViewType _vinvBubble;
83  workViewType _work;
84 
85  KOKKOS_INLINE_FUNCTION
86  Functor( outputValueViewType outputValues_,
87  inputPointViewType inputPoints_,
88  vinvViewType vinvLine_,
89  vinvViewType vinvBubble_,
90  workViewType work_)
91  : _outputValues(outputValues_), _inputPoints(inputPoints_),
92  _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
93 
94  KOKKOS_INLINE_FUNCTION
95  void operator()(const size_type iter) const {
96  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
97  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
98 
99  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
100  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
101 
102  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
103 
104  auto vcprop = Kokkos::common_view_alloc_prop(_work);
105  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
106 
107  switch (opType) {
108  case OPERATOR_VALUE : {
109  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
110  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
111  break;
112  }
113  case OPERATOR_CURL : {
114  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
115  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
116  break;
117  }
118  default: {
119  INTREPID2_TEST_FOR_ABORT( true,
120  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::Functor) operator is not supported");
121 
122  }
123  }
124  }
125  };
126  };
127  }
128 
133  template<typename DeviceType = void,
134  typename outputValueType = double,
135  typename pointValueType = double>
137  : public Basis<DeviceType,outputValueType,pointValueType> {
138  public:
140  using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
141  using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
142  using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
143 
146  Basis_HCURL_QUAD_In_FEM(const ordinal_type order,
147  const EPointType pointType = POINTTYPE_EQUISPACED);
148 
149  using OutputViewType = typename BasisBase::OutputViewType;
150  using PointViewType = typename BasisBase::PointViewType;
151  using ScalarViewType = typename BasisBase::ScalarViewType;
152 
153  using BasisBase::getValues;
154 
155  virtual
156  void
157  getValues( OutputViewType outputValues,
158  const PointViewType inputPoints,
159  const EOperator operatorType = OPERATOR_VALUE ) const override {
160 #ifdef HAVE_INTREPID2_DEBUG
161  Intrepid2::getValues_HCURL_Args(outputValues,
162  inputPoints,
163  operatorType,
164  this->getBaseCellTopology(),
165  this->getCardinality() );
166 #endif
167  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
168  Impl::Basis_HCURL_QUAD_In_FEM::
169  getValues<DeviceType,numPtsPerEval>( outputValues,
170  inputPoints,
171  this->vinvLine_,
172  this->vinvBubble_,
173  operatorType );
174  }
175 
176  virtual void
177  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
178  ordinal_type& perThreadSpaceSize,
179  const PointViewType inputPointsconst,
180  const EOperator operatorType = OPERATOR_VALUE) const override;
181 
182  KOKKOS_INLINE_FUNCTION
183  virtual void
184  getValues(
185  OutputViewType outputValues,
186  const PointViewType inputPoints,
187  const EOperator operatorType,
188  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
189  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
190  const ordinal_type subcellDim = -1,
191  const ordinal_type subcellOrdinal = -1) const override;
192 
193  virtual
194  void
195  getDofCoords( ScalarViewType dofCoords ) const override {
196 #ifdef HAVE_INTREPID2_DEBUG
197  // Verify rank of output array.
198  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
199  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
200  // Verify 0th dimension of output array.
201  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
202  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
203  // Verify 1st dimension of output array.
204  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
205  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
206 #endif
207  Kokkos::deep_copy(dofCoords, this->dofCoords_);
208  }
209 
210  virtual
211  void
212  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
213 #ifdef HAVE_INTREPID2_DEBUG
214  // Verify rank of output array.
215  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
216  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
217  // Verify 0th dimension of output array.
218  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
219  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
220  // Verify 1st dimension of output array.
221  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
222  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
223 #endif
224  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
225  }
226 
227  virtual
228  const char*
229  getName() const override {
230  return "Intrepid2_HCURL_QUAD_In_FEM";
231  }
232 
233  virtual
234  bool
235  requireOrientation() const override {
236  return true;
237  }
238 
248  BasisPtr<DeviceType,outputValueType,pointValueType>
249  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
250  if(subCellDim == 1)
251  return Teuchos::rcp( new
253  ( this->basisDegree_ - 1, POINTTYPE_GAUSS ) );
254 
255  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
256  }
257 
258  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
259  getHostBasis() const override{
261  }
262 
263  private:
264 
266  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinvLine_, vinvBubble_;
267  EPointType pointType_;
268  };
269 
270 }// namespace Intrepid2
271 
272 
274 
275 #endif
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
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.
Implementation of the default H(curl)-compatible FEM basis on Quadrilateral cell. ...
Basis_HCURL_QUAD_In_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 (...
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual bool requireOrientation() const override
True if orientation is required.
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_HCURL_QUAD_In_FEM.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(curl) functions on QUAD cells...
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.
virtual const char * getName() const override
Returns basis name.
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...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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.