Intrepid2
Intrepid2_HGRAD_TRI_Cn_FEM.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 #include "Intrepid2_PointTools.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 namespace Intrepid2 {
60 
82  namespace Impl {
83 
88 
89  public:
90  typedef struct Triangle<3> cell_topology_type;
96  template<EOperator opType>
97  struct Serial {
98  template<typename outputValueViewType,
99  typename inputPointViewType,
100  typename workViewType,
101  typename vinvViewType>
102  KOKKOS_INLINE_FUNCTION
103  static void
104  getValues( outputValueViewType outputValues,
105  const inputPointViewType inputPoints,
106  workViewType work,
107  const vinvViewType vinv );
108  };
109 
110  template<typename DeviceType, ordinal_type numPtsPerEval,
111  typename outputValueValueType, class ...outputValueProperties,
112  typename inputPointValueType, class ...inputPointProperties,
113  typename vinvValueType, class ...vinvProperties>
114  static void
115  getValues(const typename DeviceType::execution_space& space,
116  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
117  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
118  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
119  const EOperator operatorType);
120 
124  template<typename outputValueViewType,
125  typename inputPointViewType,
126  typename vinvViewType,
127  typename workViewType,
128  EOperator opType,
129  ordinal_type numPtsEval>
130  struct Functor {
131  outputValueViewType _outputValues;
132  const inputPointViewType _inputPoints;
133  const vinvViewType _vinv;
134  workViewType _work;
135 
136  KOKKOS_INLINE_FUNCTION
137  Functor( outputValueViewType outputValues_,
138  inputPointViewType inputPoints_,
139  vinvViewType vinv_,
140  workViewType work_)
141  : _outputValues(outputValues_), _inputPoints(inputPoints_),
142  _vinv(vinv_), _work(work_) {}
143 
144  KOKKOS_INLINE_FUNCTION
145  void operator()(const size_type iter) const {
146  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
147  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
148 
149  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
150  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
151 
152  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
153 
154  auto vcprop = Kokkos::common_view_alloc_prop(_work);
155  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
156 
157  switch (opType) {
158  case OPERATOR_VALUE : {
159  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
160  Serial<opType>::getValues( output, input, work, _vinv );
161  break;
162  }
163  case OPERATOR_CURL:
164  case OPERATOR_D1:
165  case OPERATOR_D2: {
166  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
167  Serial<opType>::getValues( output, input, work, _vinv );
168  break;
169  }
170  default: {
171  INTREPID2_TEST_FOR_ABORT( true,
172  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
173 
174  }
175  }
176  }
177  };
178  };
179  }
180 
181  template<typename DeviceType = void,
182  typename outputValueType = double,
183  typename pointValueType = double>
185  : public Basis<DeviceType,outputValueType,pointValueType> {
186  public:
189  using typename BasisBase::ExecutionSpace;
190  using typename BasisBase::OrdinalTypeArray1DHost;
191  using typename BasisBase::OrdinalTypeArray2DHost;
192  using typename BasisBase::OrdinalTypeArray3DHost;
193 
194  using typename BasisBase::OutputViewType;
195  using typename BasisBase::PointViewType;
196  using typename BasisBase::ScalarViewType;
197 
198  using typename BasisBase::scalarType;
199 
200  private:
201 
204  Kokkos::DynRankView<scalarType,DeviceType> vinv_;
205 
207  EPointType pointType_;
208 
209  public:
212  Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order,
213  const EPointType pointType = POINTTYPE_EQUISPACED);
214 
215  using BasisBase::getValues;
216 
217  virtual
218  void
219  getValues( const ExecutionSpace& space,
220  OutputViewType outputValues,
221  const PointViewType inputPoints,
222  const EOperator operatorType = OPERATOR_VALUE) const override {
223 #ifdef HAVE_INTREPID2_DEBUG
224  Intrepid2::getValues_HGRAD_Args(outputValues,
225  inputPoints,
226  operatorType,
227  this->getBaseCellTopology(),
228  this->getCardinality() );
229 #endif
230  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
231  Impl::Basis_HGRAD_TRI_Cn_FEM::
232  getValues<DeviceType,numPtsPerEval>(space,
233  outputValues,
234  inputPoints,
235  this->vinv_,
236  operatorType);
237  }
238 
239  virtual
240  void
241  getDofCoords( ScalarViewType dofCoords ) const override {
242 #ifdef HAVE_INTREPID2_DEBUG
243  // Verify rank of output array.
244  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
245  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
246  // Verify 0th dimension of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
249  // Verify 1st dimension of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
252 #endif
253  Kokkos::deep_copy(dofCoords, this->dofCoords_);
254  }
255 
256  virtual
257  void
258  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
259 #ifdef HAVE_INTREPID2_DEBUG
260  // Verify rank of output array.
261  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
262  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
263  // Verify 0th dimension of output array.
264  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
265  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
266 #endif
267  Kokkos::deep_copy(dofCoeffs, 1.0);
268  }
269 
270 
271  virtual
272  const char*
273  getName() const override {
274  return "Intrepid2_HGRAD_TRI_Cn_FEM";
275  }
276 
277  virtual
278  bool
279  requireOrientation() const override {
280  return (this->basisDegree_ > 2);
281  }
282 
283  void
284  getVandermondeInverse( ScalarViewType vinv ) const {
285  // has to be same rank and dimensions
286  Kokkos::deep_copy(vinv, this->vinv_);
287  }
288 
289  Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
290  getVandermondeInverse() const {
291  return vinv_;
292  }
293 
294  ordinal_type
295  getWorkSizePerPoint(const EOperator operatorType) const {
296  auto cardinality = getPnCardinality<2>(this->basisDegree_);
297  switch (operatorType) {
298  case OPERATOR_GRAD:
299  case OPERATOR_CURL:
300  case OPERATOR_D1:
301  return 5*cardinality;
302  default:
303  return getDkCardinality(operatorType, 2)*cardinality;
304  }
305  }
306 
315  BasisPtr<DeviceType,outputValueType,pointValueType>
316  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
317  if(subCellDim == 1) {
318  return Teuchos::rcp(new
320  (this->basisDegree_, pointType_));
321  }
322  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
323  }
324 
325  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
326  getHostBasis() const override{
328  }
329  };
330 
331 }// namespace Intrepid2
332 
334 
335 #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
virtual bool requireOrientation() const override
True if orientation is required.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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 (...
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.
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.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
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.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Definition file for FEM basis functions of degree n for H(grad) functions on TRI cells.
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...
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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...
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.
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 abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.