Intrepid2
Intrepid2_HDIV_TRI_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_HDIV_TRI_IN_FEM_HPP__
17 #define __INTREPID2_HDIV_TRI_IN_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
22 
23 #include "Intrepid2_PointTools.hpp"
24 #include "Teuchos_LAPACK.hpp"
25 
26 namespace Intrepid2 {
27 
56 #define CardinalityHDivTri(order) (order*(order+2))
57 
58 namespace Impl {
59 
64 public:
65  typedef struct Triangle<3> cell_topology_type;
66 
70  template<EOperator opType>
71  struct Serial {
72  template<typename outputValueViewType,
73  typename inputPointViewType,
74  typename workViewType,
75  typename vinvViewType>
76  KOKKOS_INLINE_FUNCTION
77  static void
78  getValues( outputValueViewType outputValues,
79  const inputPointViewType inputPoints,
80  workViewType work,
81  const vinvViewType vinv );
82 
83  KOKKOS_INLINE_FUNCTION
84  static ordinal_type
85  getWorkSizePerPoint(ordinal_type order) {
86  auto cardinality = CardinalityHDivTri(order);
87  switch (opType) {
88  case OPERATOR_GRAD:
89  case OPERATOR_DIV:
90  case OPERATOR_D1:
91  return 5*cardinality;
92  default:
93  return getDkCardinality<opType,2>()*cardinality;
94  }
95  }
96  };
97 
98  template<typename DeviceType, ordinal_type numPtsPerEval,
99  typename outputValueValueType, class ...outputValueProperties,
100  typename inputPointValueType, class ...inputPointProperties,
101  typename vinvValueType, class ...vinvProperties>
102  static void
103  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
104  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
105  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
106  const EOperator operatorType);
107 
111  template<typename outputValueViewType,
112  typename inputPointViewType,
113  typename vinvViewType,
114  typename workViewType,
115  EOperator opType,
116  ordinal_type numPtsEval>
117  struct Functor {
118  outputValueViewType _outputValues;
119  const inputPointViewType _inputPoints;
120  const vinvViewType _coeffs;
121  workViewType _work;
122 
123  KOKKOS_INLINE_FUNCTION
124  Functor( outputValueViewType outputValues_,
125  inputPointViewType inputPoints_,
126  vinvViewType coeffs_,
127  workViewType work_)
128  : _outputValues(outputValues_), _inputPoints(inputPoints_),
129  _coeffs(coeffs_), _work(work_) {}
130 
131  KOKKOS_INLINE_FUNCTION
132  void operator()(const size_type iter) const {
133  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
134  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
135 
136  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
137  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
138 
139  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
140 
141  auto vcprop = Kokkos::common_view_alloc_prop(_work);
142  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
143 
144 
145  switch (opType) {
146  case OPERATOR_VALUE : {
147  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
148  Serial<opType>::getValues( output, input, work, _coeffs );
149  break;
150  }
151  case OPERATOR_DIV: {
152  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
153  Serial<opType>::getValues( output, input, work, _coeffs );
154  break;
155  }
156  default: {
157  INTREPID2_TEST_FOR_ABORT( true,
158  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
159 
160  }
161  }
162  }
163  };
164 };
165 }
166 
167 template<typename DeviceType = void,
168  typename outputValueType = double,
169  typename pointValueType = double>
171  : public Basis<DeviceType,outputValueType,pointValueType> {
172  public:
174  using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
175  using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
176  using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
177 
180  Basis_HDIV_TRI_In_FEM(const ordinal_type order,
181  const EPointType pointType = POINTTYPE_EQUISPACED);
182 
184 
185  using OutputViewType = typename BasisBase::OutputViewType;
186  using PointViewType = typename BasisBase::PointViewType;
187  using ScalarViewType = typename BasisBase::ScalarViewType;
188  using scalarType = typename BasisBase::scalarType;
189  using BasisBase::getValues;
190 
191  virtual
192  void
193  getValues( OutputViewType outputValues,
194  const PointViewType inputPoints,
195  const EOperator operatorType = OPERATOR_VALUE) const override {
196 #ifdef HAVE_INTREPID2_DEBUG
197  Intrepid2::getValues_HDIV_Args(outputValues,
198  inputPoints,
199  operatorType,
200  this->getBaseCellTopology(),
201  this->getCardinality() );
202 #endif
203  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
204  Impl::Basis_HDIV_TRI_In_FEM::
205  getValues<DeviceType,numPtsPerEval>( outputValues,
206  inputPoints,
207  this->coeffs_,
208  operatorType);
209  }
210 
211  virtual void
212  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
213  ordinal_type& perThreadSpaceSize,
214  const PointViewType inputPointsconst,
215  const EOperator operatorType = OPERATOR_VALUE) const override;
216 
217  KOKKOS_INLINE_FUNCTION
218  virtual void
219  getValues(
220  OutputViewType outputValues,
221  const PointViewType inputPoints,
222  const EOperator operatorType,
223  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
224  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
225  const ordinal_type subcellDim = -1,
226  const ordinal_type subcellOrdinal = -1) const override;
227 
228  virtual
229  void
230  getDofCoords( ScalarViewType dofCoords ) const override {
231 #ifdef HAVE_INTREPID2_DEBUG
232  // Verify rank of output array.
233  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
234  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
235  // Verify 0th dimension of output array.
236  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
237  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
238  // Verify 1st dimension of output array.
239  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
240  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
241 #endif
242  Kokkos::deep_copy(dofCoords, this->dofCoords_);
243  }
244 
245  virtual
246  void
247  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
248 #ifdef HAVE_INTREPID2_DEBUG
249  // Verify rank of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
252  // Verify 0th dimension of output array.
253  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
254  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
255  // Verify 1st dimension of output array.
256  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
257  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
258 #endif
259  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
260  }
261 
262  void
263  getExpansionCoeffs( ScalarViewType coeffs ) const {
264  // has to be same rank and dimensions
265  Kokkos::deep_copy(coeffs, this->coeffs_);
266  }
267 
268  virtual
269  const char*
270  getName() const override {
271  return "Intrepid2_HDIV_TRI_In_FEM";
272  }
273 
274  virtual
275  bool
276  requireOrientation() const override {
277  return true;
278  }
279 
289  BasisPtr<DeviceType,outputValueType,pointValueType>
290  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
291  if(subCellDim == 1) {
292  return Teuchos::rcp(new
294  (this->basisDegree_-1, pointType_));
295  }
296  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
297  }
298 
299  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
300  getHostBasis() const override{
302  }
303  private:
304 
307  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
308 
310  EPointType pointType_;
311 
312 };
313 
314 }// namespace Intrepid2
315 
317 
318 #endif
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
small utility functions
ordinal_type getCardinality() const
Returns cardinality of the basis.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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...
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
EPointType pointType_
type of lattice used for creating the DoF coordinates
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(div) functions on TRI 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.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
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< 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.
virtual const char * getName() const override
Returns basis name.
Header file for the abstract base class Intrepid2::Basis.
Basis_HDIV_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.