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"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 namespace Intrepid2 {
59 
81  namespace Impl {
82 
87  public:
88  typedef struct Triangle<3> cell_topology_type;
94  template<EOperator opType>
95  struct Serial {
96  template<typename outputValueViewType,
97  typename inputPointViewType,
98  typename workViewType,
99  typename vinvViewType>
100  KOKKOS_INLINE_FUNCTION
101  static void
102  getValues( outputValueViewType outputValues,
103  const inputPointViewType inputPoints,
104  workViewType work,
105  const vinvViewType vinv );
106  };
107 
108  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
109  typename outputValueValueType, class ...outputValueProperties,
110  typename inputPointValueType, class ...inputPointProperties,
111  typename vinvValueType, class ...vinvProperties>
112  static void
113  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
114  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
115  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
116  const EOperator operatorType);
117 
121  template<typename outputValueViewType,
122  typename inputPointViewType,
123  typename vinvViewType,
124  typename workViewType,
125  EOperator opType,
126  ordinal_type numPtsEval>
127  struct Functor {
128  outputValueViewType _outputValues;
129  const inputPointViewType _inputPoints;
130  const vinvViewType _vinv;
131  workViewType _work;
132 
133  KOKKOS_INLINE_FUNCTION
134  Functor( outputValueViewType outputValues_,
135  inputPointViewType inputPoints_,
136  vinvViewType vinv_,
137  workViewType work_)
138  : _outputValues(outputValues_), _inputPoints(inputPoints_),
139  _vinv(vinv_), _work(work_) {}
140 
141  KOKKOS_INLINE_FUNCTION
142  void operator()(const size_type iter) const {
143  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
144  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
145 
146  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
147  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
148 
149  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
150 
151  auto vcprop = Kokkos::common_view_alloc_prop(_work);
152  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
153 
154  switch (opType) {
155  case OPERATOR_VALUE : {
156  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
157  Serial<opType>::getValues( output, input, work, _vinv );
158  break;
159  }
160  case OPERATOR_CURL:
161  case OPERATOR_D1:
162  case OPERATOR_D2: {
163  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
164  Serial<opType>::getValues( output, input, work, _vinv );
165  break;
166  }
167  default: {
168  INTREPID2_TEST_FOR_ABORT( true,
169  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
170 
171  }
172  }
173  }
174  };
175  };
176  }
177 
178  template<typename ExecSpaceType = void,
179  typename outputValueType = double,
180  typename pointValueType = double>
182  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
183  public:
187 
188  using ordinal_type_array_1d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray1DHost instead","OrdinalTypeArray1DHost") = OrdinalTypeArray1DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray1DHost instead");
189  using ordinal_type_array_2d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray2DHost instead","OrdinalTypeArray2DHost") = OrdinalTypeArray2DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray2DHost instead");
190  using ordinal_type_array_3d_host INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OrdinalTypeArray3DHost instead","OrdinalTypeArray3DHost") = OrdinalTypeArray3DHost INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OrdinalTypeArray3DHost instead");
191 
195 
196  using outputViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OutputViewType instead","OutputViewType") = OutputViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OutputViewType instead");
197  using pointViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use PointViewType instead","PointViewType") = PointViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use PointViewType instead");
198  using scalarViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use ScalarViewType instead","ScalarViewType") = ScalarViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use ScalarViewType instead");
199 
201 
202  private:
203 
206  Kokkos::DynRankView<scalarType,ExecSpaceType> vinv_;
207 
208  public:
211  Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order,
212  const EPointType pointType = POINTTYPE_EQUISPACED);
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_TRI_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_TRI_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_TRI_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_TRI_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_TRI_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_TRI_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  virtual
269  const char*
270  getName() const {
271  return "Intrepid2_HGRAD_TRI_Cn_FEM";
272  }
273 
274  virtual
275  bool
277  return (this->basisDegree_ > 2);
278  }
279 
280  void
281  getVandermondeInverse( ScalarViewType vinv ) const {
282  // has to be same rank and dimensions
283  Kokkos::deep_copy(vinv, this->vinv_);
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<2>(this->basisDegree_);
294  switch (operatorType) {
295  case OPERATOR_GRAD:
296  case OPERATOR_CURL:
297  case OPERATOR_D1:
298  return 5*cardinality;
299  default:
300  return getDkCardinality(operatorType, 2)*cardinality;
301  }
302  }
303 
304  };
305 
306 }// namespace Intrepid2
307 
309 
310 #endif
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.
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...
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual const char * getName() const
Returns basis name.
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual bool requireOrientation() const
True if orientation is required.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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.