Intrepid2
Intrepid2_HVOL_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 
48 #ifndef __INTREPID2_HVOL_TET_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_TET_CN_FEM_HPP__
50 
51 #include "Intrepid2_Basis.hpp"
52 
53 #include "Intrepid2_PointTools.hpp"
54 #include "Teuchos_LAPACK.hpp"
55 
56 
57 namespace Intrepid2 {
58 
74  namespace Impl {
75 
80  public:
81  typedef struct Tetrahedron<4> cell_topology_type;
85  template<EOperator opType>
86  struct Serial {
87  template<typename outputValueViewType,
88  typename inputPointViewType,
89  typename workViewType,
90  typename vinvViewType>
91  KOKKOS_INLINE_FUNCTION
92  static void
93  getValues( outputValueViewType outputValues,
94  const inputPointViewType inputPoints,
95  workViewType work,
96  const vinvViewType vinv );
97 
98 
99  KOKKOS_INLINE_FUNCTION
100  static ordinal_type
101  getWorkSizePerPoint(ordinal_type order) {
102  auto cardinality = getPnCardinality<3>(order);
103  switch (opType) {
104  case OPERATOR_GRAD:
105  case OPERATOR_CURL:
106  case OPERATOR_D1:
107  return 7*cardinality;
108  default:
109  return getDkCardinality<opType,3>()*cardinality;
110  }
111  }
112  };
113 
114  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
115  typename outputValueValueType, class ...outputValueProperties,
116  typename inputPointValueType, class ...inputPointProperties,
117  typename vinvValueType, class ...vinvProperties>
118  static void
119  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
120  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
121  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
122  const EOperator operatorType);
123 
127  template<typename outputValueViewType,
128  typename inputPointViewType,
129  typename vinvViewType,
130  typename workViewType,
131  EOperator opType,
132  ordinal_type numPtsEval>
133  struct Functor {
134  outputValueViewType _outputValues;
135  const inputPointViewType _inputPoints;
136  const vinvViewType _vinv;
137  workViewType _work;
138 
139  KOKKOS_INLINE_FUNCTION
140  Functor( outputValueViewType outputValues_,
141  inputPointViewType inputPoints_,
142  vinvViewType vinv_,
143  workViewType work_)
144  : _outputValues(outputValues_), _inputPoints(inputPoints_),
145  _vinv(vinv_), _work(work_) {}
146 
147  KOKKOS_INLINE_FUNCTION
148  void operator()(const size_type iter) const {
149  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
150  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
151 
152  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
153  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
154 
155  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
156 
157  auto vcprop = Kokkos::common_view_alloc_prop(_work);
158  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
159 
160  switch (opType) {
161  case OPERATOR_VALUE : {
162  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
163  Serial<opType>::getValues( output, input, work, _vinv );
164  break;
165  }
166  case OPERATOR_GRAD :
167  case OPERATOR_D1 :
168  case OPERATOR_D2 :
169  //case OPERATOR_D3 :
170  {
171  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
172  Serial<opType>::getValues( output, input, work, _vinv );
173  break;
174  }
175  default: {
176  INTREPID2_TEST_FOR_ABORT( true,
177  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::Functor) operator is not supported");
178 
179  }
180  }
181  }
182  };
183  };
184  }
185 
186  template<typename ExecSpaceType = void,
187  typename outputValueType = double,
188  typename pointValueType = double>
190  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
191  public:
195 
198  Basis_HVOL_TET_Cn_FEM(const ordinal_type order,
199  const EPointType pointType = POINTTYPE_EQUISPACED);
200 
201 
206 
208 
209  virtual
210  void
211  getValues( outputViewType outputValues,
212  const pointViewType inputPoints,
213  const EOperator operatorType = OPERATOR_VALUE) const {
214 #ifdef HAVE_INTREPID2_DEBUG
215  Intrepid2::getValues_HVOL_Args(outputValues,
216  inputPoints,
217  operatorType,
218  this->getBaseCellTopology(),
219  this->getCardinality() );
220 #endif
221  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
222  Impl::Basis_HVOL_TET_Cn_FEM::
223  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
224  inputPoints,
225  this->vinv_,
226  operatorType);
227  }
228 
229  virtual
230  void
231  getDofCoords( scalarViewType dofCoords ) const {
232 #ifdef HAVE_INTREPID2_DEBUG
233  // Verify rank of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
236  // Verify 0th dimension of output array.
237  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
238  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
239  // Verify 1st dimension of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
242 #endif
243  Kokkos::deep_copy(dofCoords, this->dofCoords_);
244  }
245 
246  virtual
247  void
248  getDofCoeffs( scalarViewType dofCoeffs ) const {
249 #ifdef HAVE_INTREPID2_DEBUG
250  // Verify rank of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
253  // Verify 0th dimension of output array.
254  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
255  ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
256 #endif
257  Kokkos::deep_copy(dofCoeffs, 1.0);
258  }
259 
260  void
261  getVandermondeInverse( scalarViewType vinv ) const {
262  // has to be same rank and dimensions
263  Kokkos::deep_copy(vinv, this->vinv_);
264  }
265 
266  virtual
267  const char*
268  getName() const {
269  return "Intrepid2_HVOL_TET_Cn_FEM";
270  }
271 
272  virtual
273  bool
275  return false;
276  }
277 
278  private:
279 
282  Kokkos::DynRankView<scalarType,ExecSpaceType> vinv_;
283 
284  };
285 
286 }// namespace Intrepid2
287 
289 
290 #endif
small utility functions
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
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.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
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.
Basis_HVOL_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
virtual const char * getName() const
Returns basis name.
virtual bool requireOrientation() const
True if orientation is required.
Definition file for FEM basis functions of degree n for H(vol) functions on TET.
See Intrepid2::Basis_HVOL_TET_Cn_FEM.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Tetrahedron cell...
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.