Intrepid2
Intrepid2_HVOL_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 
48 #ifndef __INTREPID2_HVOL_TRI_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_TRI_CN_FEM_HPP__
50 
51 #include "Intrepid2_Basis.hpp"
52 
53 #include "Intrepid2_PointTools.hpp"
54 #include "Teuchos_LAPACK.hpp"
55 
56 namespace Intrepid2 {
57 
72  namespace Impl {
73 
78  public:
79  typedef struct Triangle<3> cell_topology_type;
83  template<EOperator opType>
84  struct Serial {
85  template<typename outputValueViewType,
86  typename inputPointViewType,
87  typename workViewType,
88  typename vinvViewType>
89  KOKKOS_INLINE_FUNCTION
90  static void
91  getValues( outputValueViewType outputValues,
92  const inputPointViewType inputPoints,
93  workViewType work,
94  const vinvViewType vinv );
95 
96 
97  KOKKOS_INLINE_FUNCTION
98  static ordinal_type
99  getWorkSizePerPoint(ordinal_type order) {
100  auto cardinality = getPnCardinality<2>(order);
101  switch (opType) {
102  case OPERATOR_GRAD:
103  case OPERATOR_CURL:
104  case OPERATOR_D1:
105  return 5*cardinality;
106  default:
107  return getDkCardinality<opType,2>()*cardinality;
108  }
109  }
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( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
118  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
119  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
120  const EOperator operatorType);
121 
125  template<typename outputValueViewType,
126  typename inputPointViewType,
127  typename vinvViewType,
128  typename workViewType,
129  EOperator opType,
130  ordinal_type numPtsEval>
131  struct Functor {
132  outputValueViewType _outputValues;
133  const inputPointViewType _inputPoints;
134  const vinvViewType _vinv;
135  workViewType _work;
136 
137  KOKKOS_INLINE_FUNCTION
138  Functor( outputValueViewType outputValues_,
139  inputPointViewType inputPoints_,
140  vinvViewType vinv_,
141  workViewType work_)
142  : _outputValues(outputValues_), _inputPoints(inputPoints_),
143  _vinv(vinv_), _work(work_) {}
144 
145  KOKKOS_INLINE_FUNCTION
146  void operator()(const size_type iter) const {
147  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
148  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
149 
150  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
151  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
152 
153  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
154 
155  auto vcprop = Kokkos::common_view_alloc_prop(_work);
156  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
157 
158  switch (opType) {
159  case OPERATOR_VALUE : {
160  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
161  Serial<opType>::getValues( output, input, work, _vinv );
162  break;
163  }
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_HVOL_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:
188 
189  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
190  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
191  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
192 
195  Basis_HVOL_TRI_Cn_FEM(const ordinal_type order,
196  const EPointType pointType = POINTTYPE_EQUISPACED);
197 
198 
202 
204 
206 
207  virtual
208  void
209  getValues( OutputViewType outputValues,
210  const PointViewType inputPoints,
211  const EOperator operatorType = OPERATOR_VALUE) const override {
212 #ifdef HAVE_INTREPID2_DEBUG
213  Intrepid2::getValues_HVOL_Args(outputValues,
214  inputPoints,
215  operatorType,
216  this->getBaseCellTopology(),
217  this->getCardinality() );
218 #endif
219  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
220  Impl::Basis_HVOL_TRI_Cn_FEM::
221  getValues<DeviceType,numPtsPerEval>( outputValues,
222  inputPoints,
223  this->vinv_,
224  operatorType);
225  }
226 
227  virtual
228  void
229  getDofCoords( ScalarViewType dofCoords ) const override {
230 #ifdef HAVE_INTREPID2_DEBUG
231  // Verify rank of output array.
232  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
233  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
234  // Verify 0th dimension of output array.
235  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
236  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
237  // Verify 1st dimension of output array.
238  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
239  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
240 #endif
241  Kokkos::deep_copy(dofCoords, this->dofCoords_);
242  }
243 
244  virtual
245  void
246  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
247 #ifdef HAVE_INTREPID2_DEBUG
248  // Verify rank of output array.
249  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
250  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
251  // Verify 0th dimension of output array.
252  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
253  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
254 #endif
255  Kokkos::deep_copy(dofCoeffs, 1.0);
256  }
257 
258  void
259  getVandermondeInverse( ScalarViewType vinv ) const {
260  // has to be same rank and dimensions
261  Kokkos::deep_copy(vinv, this->vinv_);
262  }
263 
264  virtual
265  const char*
266  getName() const override {
267  return "Intrepid2_HVOL_TRI_Cn_FEM";
268  }
269 
270  virtual
271  bool
272  requireOrientation() const override {
273  return false;
274  }
275 
277  getHostBasis() const override{
279  }
280 
281  private:
282 
285  Kokkos::DynRankView<scalarType,DeviceType> vinv_;
286  EPointType pointType_;
287 
288  };
289 
290 }// namespace Intrepid2
291 
293 
294 #endif
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
small utility functions
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
ordinal_type getCardinality() const
Returns cardinality of the basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type), allowing host access to input and output views, and ensuring host execution of basis evaluation.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
virtual const char * getName() const override
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
void getValues_HVOL_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HVOL-conforming FEM basis...
EPointType
Enumeration of types of point distributions in Intrepid.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Definition file for FEM basis functions of degree n for H(vol) functions on TRI.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.
virtual bool requireOrientation() const override
True if orientation is required.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.