Intrepid2
Intrepid2_HGRAD_TET_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_TET_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_TET_CN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 #include "Intrepid2_PointTools.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 
60 namespace Intrepid2 {
61 
87  namespace Impl {
88 
93  public:
94  typedef struct Tetrahedron<4> cell_topology_type;
98  template<EOperator opType>
99  struct Serial {
100  template<typename outputValueViewType,
101  typename inputPointViewType,
102  typename workViewType,
103  typename vinvViewType>
104  KOKKOS_INLINE_FUNCTION
105  static void
106  getValues( outputValueViewType outputValues,
107  const inputPointViewType inputPoints,
108  workViewType work,
109  const vinvViewType vinv );
110  };
111 
112  template<typename DeviceType, ordinal_type numPtsPerEval,
113  typename outputValueValueType, class ...outputValueProperties,
114  typename inputPointValueType, class ...inputPointProperties,
115  typename vinvValueType, class ...vinvProperties>
116  static void
117  getValues( const typename DeviceType::execution_space& space,
118  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
119  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
120  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
121  const EOperator operatorType);
122 
126  template<typename outputValueViewType,
127  typename inputPointViewType,
128  typename vinvViewType,
129  typename workViewType,
130  EOperator opType,
131  ordinal_type numPtsEval>
132  struct Functor {
133  outputValueViewType _outputValues;
134  const inputPointViewType _inputPoints;
135  const vinvViewType _vinv;
136  workViewType _work;
137 
138  KOKKOS_INLINE_FUNCTION
139  Functor( outputValueViewType outputValues_,
140  inputPointViewType inputPoints_,
141  vinvViewType vinv_,
142  workViewType work_)
143  : _outputValues(outputValues_), _inputPoints(inputPoints_),
144  _vinv(vinv_), _work(work_) {}
145 
146  KOKKOS_INLINE_FUNCTION
147  void operator()(const size_type iter) const {
148  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
149  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
150 
151  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
152  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
153 
154  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
155 
156  auto vcprop = Kokkos::common_view_alloc_prop(_work);
157  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
158 
159  switch (opType) {
160  case OPERATOR_VALUE : {
161  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
162  Serial<opType>::getValues( output, input, work, _vinv );
163  break;
164  }
165  case OPERATOR_GRAD :
166  case OPERATOR_D1 :
167  case OPERATOR_D2 :
168  //case OPERATOR_D3 :
169  {
170  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
171  Serial<opType>::getValues( output, input, work, _vinv );
172  break;
173  }
174  default: {
175  INTREPID2_TEST_FOR_ABORT( true,
176  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::Functor) operator is not supported");
177 
178  }
179  }
180  }
181  };
182  };
183  }
184 
185  template<typename DeviceType = void,
186  typename outputValueType = double,
187  typename pointValueType = double>
189  : public Basis<DeviceType,outputValueType,pointValueType> {
190  public:
192  using typename BasisBase::ExecutionSpace;
193  using typename BasisBase::OrdinalTypeArray1DHost;
194  using typename BasisBase::OrdinalTypeArray2DHost;
195  using typename BasisBase::OrdinalTypeArray3DHost;
196 
197  using typename BasisBase::OutputViewType;
198  using typename BasisBase::PointViewType;
199  using typename BasisBase::ScalarViewType;
200 
201  using typename BasisBase::scalarType;
202 
203  private:
204 
207  Kokkos::DynRankView<scalarType,DeviceType> vinv_;
208 
210  EPointType pointType_;
211 
212  public:
213 
216  Basis_HGRAD_TET_Cn_FEM(const ordinal_type order,
217  const EPointType pointType = POINTTYPE_EQUISPACED);
218 
219  using BasisBase::getValues;
220 
221  virtual
222  void
223  getValues(const ExecutionSpace& space,
224  OutputViewType outputValues,
225  const PointViewType inputPoints,
226  const EOperator operatorType = OPERATOR_VALUE) const override {
227 #ifdef HAVE_INTREPID2_DEBUG
228  Intrepid2::getValues_HGRAD_Args(outputValues,
229  inputPoints,
230  operatorType,
231  this->getBaseCellTopology(),
232  this->getCardinality() );
233 #endif
234  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
235  Impl::Basis_HGRAD_TET_Cn_FEM::
236  getValues<DeviceType,numPtsPerEval>(space,
237  outputValues,
238  inputPoints,
239  this->vinv_,
240  operatorType);
241  }
242 
243  virtual
244  void
245  getDofCoords( ScalarViewType dofCoords ) const override {
246 #ifdef HAVE_INTREPID2_DEBUG
247  // Verify rank of output array.
248  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
249  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
250  // Verify 0th dimension of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
253  // Verify 1st dimension of output array.
254  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
255  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
256 #endif
257  Kokkos::deep_copy(dofCoords, this->dofCoords_);
258  }
259 
260  virtual
261  void
262  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
263 #ifdef HAVE_INTREPID2_DEBUG
264  // Verify rank of output array.
265  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
266  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
267  // Verify 0th dimension of output array.
268  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
269  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
270 #endif
271  Kokkos::deep_copy(dofCoeffs, 1.0);
272  }
273 
274 
275  void
276  getVandermondeInverse( ScalarViewType vinv ) const {
277  // has to be same rank and dimensions
278  Kokkos::deep_copy(vinv, this->vinv_);
279  }
280 
281  virtual
282  const char*
283  getName() const override {
284  return "Intrepid2_HGRAD_TET_Cn_FEM";
285  }
286 
287  virtual
288  bool
289  requireOrientation() const override {
290  return (this->basisDegree_ > 2);
291  }
292 
293  Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
294 
295  getVandermondeInverse() const {
296  return vinv_;
297  }
298 
299  ordinal_type
300  getWorkSizePerPoint(const EOperator operatorType) const {
301  auto cardinality = getPnCardinality<3>(this->basisDegree_);
302  switch (operatorType) {
303  case OPERATOR_GRAD:
304  case OPERATOR_CURL:
305  case OPERATOR_D1:
306  return 7*cardinality;
307  default:
308  return getDkCardinality(operatorType, 3)*cardinality;
309  }
310  }
311 
320  BasisPtr<DeviceType,outputValueType,pointValueType>
321  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
322  if(subCellDim == 1) {
323  return Teuchos::rcp(new
325  (this->basisDegree_, pointType_));
326  } else if(subCellDim == 2) {
327  return Teuchos::rcp(new
329  (this->basisDegree_, pointType_));
330  }
331 
332  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
333  }
334 
335  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
336  getHostBasis() const override{
338  }
339  };
340 
341 }// namespace Intrepid2
342 
344 
345 #endif
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
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< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
small utility functions
ordinal_type getCardinality() const
Returns cardinality of the basis.
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 void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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 bool requireOrientation() const override
True if orientation is required.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual const char * getName() const override
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
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...
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Basis_HGRAD_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
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.
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.