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"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 
59 namespace Intrepid2 {
60 
86  namespace Impl {
87 
92  public:
93  typedef struct Tetrahedron<4> cell_topology_type;
97  template<EOperator opType>
98  struct Serial {
99  template<typename outputValueViewType,
100  typename inputPointViewType,
101  typename workViewType,
102  typename vinvViewType>
103  KOKKOS_INLINE_FUNCTION
104  static void
105  getValues( outputValueViewType outputValues,
106  const inputPointViewType inputPoints,
107  workViewType work,
108  const vinvViewType vinv );
109  };
110 
111  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
112  typename outputValueValueType, class ...outputValueProperties,
113  typename inputPointValueType, class ...inputPointProperties,
114  typename vinvValueType, class ...vinvProperties>
115  static void
116  getValues( 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_GRAD :
164  case OPERATOR_D1 :
165  case OPERATOR_D2 :
166  //case OPERATOR_D3 :
167  {
168  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
169  Serial<opType>::getValues( output, input, work, _vinv );
170  break;
171  }
172  default: {
173  INTREPID2_TEST_FOR_ABORT( true,
174  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::Functor) operator is not supported");
175 
176  }
177  }
178  }
179  };
180  };
181  }
182 
183  template<typename ExecSpaceType = void,
184  typename outputValueType = double,
185  typename pointValueType = double>
187  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
188  public:
192 
196 
198 
199  private:
200 
203  Kokkos::DynRankView<scalarType,ExecSpaceType> vinv_;
204 
205  public:
206 
209  Basis_HGRAD_TET_Cn_FEM(const ordinal_type order,
210  const EPointType pointType = POINTTYPE_EQUISPACED);
211 
212 
213 
215 
216  virtual
217  void
218  getValues( OutputViewType outputValues,
219  const PointViewType inputPoints,
220  const EOperator operatorType = OPERATOR_VALUE) const {
221 #ifdef HAVE_INTREPID2_DEBUG
222  Intrepid2::getValues_HGRAD_Args(outputValues,
223  inputPoints,
224  operatorType,
225  this->getBaseCellTopology(),
226  this->getCardinality() );
227 #endif
228  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
229  Impl::Basis_HGRAD_TET_Cn_FEM::
230  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
231  inputPoints,
232  this->vinv_,
233  operatorType);
234  }
235 
236  virtual
237  void
238  getDofCoords( ScalarViewType dofCoords ) const {
239 #ifdef HAVE_INTREPID2_DEBUG
240  // Verify rank of output array.
241  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
242  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
243  // Verify 0th dimension of output array.
244  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
245  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
246  // Verify 1st dimension of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
249 #endif
250  Kokkos::deep_copy(dofCoords, this->dofCoords_);
251  }
252 
253  virtual
254  void
255  getDofCoeffs( ScalarViewType dofCoeffs ) const {
256 #ifdef HAVE_INTREPID2_DEBUG
257  // Verify rank of output array.
258  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
259  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
260  // Verify 0th dimension of output array.
261  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
262  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
263 #endif
264  Kokkos::deep_copy(dofCoeffs, 1.0);
265  }
266 
267 
268  void
269  getVandermondeInverse( ScalarViewType vinv ) const {
270  // has to be same rank and dimensions
271  Kokkos::deep_copy(vinv, this->vinv_);
272  }
273 
274  virtual
275  const char*
276  getName() const {
277  return "Intrepid2_HGRAD_TET_Cn_FEM";
278  }
279 
280  virtual
281  bool
283  return (this->basisDegree_ > 2);
284  }
285 
286  Kokkos::DynRankView<typename ScalarViewType::const_value_type,ExecSpaceType>
287  getVandermondeInverse() const {
288  return vinv_;
289  }
290 
291  ordinal_type
292  getWorkSizePerPoint(const EOperator operatorType) const {
293  auto cardinality = getPnCardinality<3>(this->basisDegree_);
294  switch (operatorType) {
295  case OPERATOR_GRAD:
296  case OPERATOR_CURL:
297  case OPERATOR_D1:
298  return 7*cardinality;
299  default:
300  return getDkCardinality(operatorType, 3)*cardinality;
301  }
302  }
303  };
304 
305 }// namespace Intrepid2
306 
308 
309 #endif
virtual const char * getName() const
Returns basis name.
small utility functions
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Basis_HGRAD_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual bool requireOrientation() const
True if orientation is required.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
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.