Intrepid2
Intrepid2_HVOL_LINE_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_LINE_CN_FEM_HPP__
16 #define __INTREPID2_HVOL_LINE_CN_FEM_HPP__
17 
18 #include "Intrepid2_Basis.hpp"
20 
21 #include "Intrepid2_PointTools.hpp"
22 #include "Teuchos_LAPACK.hpp"
23 
24 namespace Intrepid2 {
25 
41  namespace Impl {
42 
47  public:
48  typedef struct Line<2> cell_topology_type;
52  template<EOperator opType>
53  struct Serial {
54  template<typename outputValueViewType,
55  typename inputPointViewType,
56  typename workViewType,
57  typename vinvViewType>
58  KOKKOS_INLINE_FUNCTION
59  static void
60  getValues( outputValueViewType outputValues,
61  const inputPointViewType inputPoints,
62  workViewType work,
63  const vinvViewType vinv,
64  const ordinal_type operatorDn = 0 );
65 
66  KOKKOS_INLINE_FUNCTION
67  static ordinal_type
68  getWorkSizePerPoint(ordinal_type order) {return getPnCardinality<1>(order);}
69  };
70 
71  template<typename DeviceType, ordinal_type numPtsPerEval,
72  typename outputValueValueType, class ...outputValueProperties,
73  typename inputPointValueType, class ...inputPointProperties,
74  typename vinvValueType, class ...vinvProperties>
75  static void
76  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
77  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
78  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
79  const EOperator operatorType );
80 
84  template<typename outputValueViewType,
85  typename inputPointViewType,
86  typename vinvViewType,
87  typename workViewType,
88  EOperator opType,
89  ordinal_type numPtsEval>
90  struct Functor {
91  outputValueViewType _outputValues;
92  const inputPointViewType _inputPoints;
93  const vinvViewType _vinv;
94  workViewType _work;
95  const ordinal_type _opDn;
96 
97  KOKKOS_INLINE_FUNCTION
98  Functor( outputValueViewType outputValues_,
99  inputPointViewType inputPoints_,
100  vinvViewType vinv_,
101  workViewType work_,
102  const ordinal_type opDn_ = 0 )
103  : _outputValues(outputValues_), _inputPoints(inputPoints_),
104  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
105 
106  KOKKOS_INLINE_FUNCTION
107  void operator()(const size_type iter) const {
108  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
109  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
110 
111  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
112  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
113 
114  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
115 
116  auto vcprop = Kokkos::common_view_alloc_prop(_work);
117  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
118 
119  switch (opType) {
120  case OPERATOR_VALUE : {
121  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
122  Serial<opType>::getValues( output, input, work, _vinv );
123  break;
124  }
125  case OPERATOR_Dn : {
126  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
127  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
128  break;
129  }
130  default: {
131  INTREPID2_TEST_FOR_ABORT( true,
132  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Functor) operator is not supported");
133 
134  }
135  }
136  }
137  };
138  };
139  }
140 
141  template<typename DeviceType = void,
142  typename outputValueType = double,
143  typename pointValueType = double>
145  : public Basis<DeviceType,outputValueType,pointValueType> {
146  public:
148 
150 
151  using typename BasisBase::OrdinalTypeArray1DHost;
152  using typename BasisBase::OrdinalTypeArray2DHost;
153  using typename BasisBase::OrdinalTypeArray3DHost;
154 
155  using typename BasisBase::OutputViewType;
156  using typename BasisBase::PointViewType ;
157  using typename BasisBase::ScalarViewType;
158 
161  Basis_HVOL_LINE_Cn_FEM(const ordinal_type order,
162  const EPointType pointType = POINTTYPE_EQUISPACED);
163 
164  using BasisBase::getValues;
165 
166  virtual
167  void
168  getValues( OutputViewType outputValues,
169  const PointViewType inputPoints,
170  const EOperator operatorType = OPERATOR_VALUE ) const override {
171 #ifdef HAVE_INTREPID2_DEBUG
172  Intrepid2::getValues_HVOL_Args(outputValues,
173  inputPoints,
174  operatorType,
175  this->getBaseCellTopology(),
176  this->getCardinality() );
177 #endif
178  constexpr ordinal_type numPtsPerEval = 1;
179  Impl::Basis_HVOL_LINE_Cn_FEM::
180  getValues<DeviceType,numPtsPerEval>( outputValues,
181  inputPoints,
182  this->vinv_,
183  operatorType );
184  }
185 
186  virtual void
187  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
188  ordinal_type& perThreadSpaceSize,
189  const PointViewType inputPointsconst,
190  const EOperator operatorType = OPERATOR_VALUE) const override;
191 
192  KOKKOS_INLINE_FUNCTION
193  virtual void
194  getValues(
195  OutputViewType outputValues,
196  const PointViewType inputPoints,
197  const EOperator operatorType,
198  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
199  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
200  const ordinal_type subcellDim = -1,
201  const ordinal_type subcellOrdinal = -1) const override;
202 
203  virtual
204  void
205  getDofCoords( ScalarViewType dofCoords ) const override {
206 #ifdef HAVE_INTREPID2_DEBUG
207  // Verify rank of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
210  // Verify 0th dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
213  // Verify 1st dimension of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
216 #endif
217  Kokkos::deep_copy(dofCoords, this->dofCoords_);
218  }
219 
220  virtual
221  void
222  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
223 #ifdef HAVE_INTREPID2_DEBUG
224  // Verify rank of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
227  // Verify 0th dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
230 #endif
231  Kokkos::deep_copy(dofCoeffs, 1.0);
232  }
233 
234  void
235  getVandermondeInverse( ScalarViewType vinv ) const {
236  // has to be same rank and dimensions
237  Kokkos::deep_copy(vinv, this->vinv_);
238  }
239 
240  virtual
241  const char*
242  getName() const override {
243  return "Intrepid2_HVOL_LINE_Cn_FEM";
244  }
245 
246  virtual HostBasisPtr<outputValueType,pointValueType>
247  getHostBasis() const override{
249  }
250 
251  private:
252 
255  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
256  EPointType pointType_;
257  };
258 
259 }// namespace Intrepid2
260 
262 
263 #endif
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
small utility functions
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
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.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
virtual const char * getName() const override
Returns basis name.
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...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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 HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.