Intrepid2
Intrepid2_HDIV_TET_In_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_HDIV_TET_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TET_IN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 #include "Intrepid2_PointTools.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 namespace Intrepid2 {
60 
87 #define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
88 
89 namespace Impl {
90 
95 public:
96 
100  template<EOperator opType>
101  struct Serial {
102  template<typename outputValueViewType,
103  typename inputPointViewType,
104  typename workViewType,
105  typename vinvViewType>
106  KOKKOS_INLINE_FUNCTION
107  static void
108  getValues( outputValueViewType outputValues,
109  const inputPointViewType inputPoints,
110  workViewType work,
111  const vinvViewType vinv );
112 
113 
114  KOKKOS_INLINE_FUNCTION
115  static ordinal_type
116  getWorkSizePerPoint(ordinal_type order) {
117  auto cardinality = CardinalityHDivTet(order);
118  switch (opType) {
119  case OPERATOR_DIV:
120  case OPERATOR_D1:
121  return 7*cardinality;
122  default:
123  return getDkCardinality<opType,3>()*cardinality;
124  }
125  }
126  };
127 
128  template<typename DeviceType, ordinal_type numPtsPerEval,
129  typename outputValueValueType, class ...outputValueProperties,
130  typename inputPointValueType, class ...inputPointProperties,
131  typename vinvValueType, class ...vinvProperties>
132  static void
133  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
134  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
135  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
136  const EOperator operatorType);
137 
141  template<typename outputValueViewType,
142  typename inputPointViewType,
143  typename vinvViewType,
144  typename workViewType,
145  EOperator opType,
146  ordinal_type numPtsEval>
147  struct Functor {
148  outputValueViewType _outputValues;
149  const inputPointViewType _inputPoints;
150  const vinvViewType _coeffs;
151  workViewType _work;
152 
153  KOKKOS_INLINE_FUNCTION
154  Functor( outputValueViewType outputValues_,
155  inputPointViewType inputPoints_,
156  vinvViewType coeffs_,
157  workViewType work_)
158  : _outputValues(outputValues_), _inputPoints(inputPoints_),
159  _coeffs(coeffs_), _work(work_) {}
160 
161  KOKKOS_INLINE_FUNCTION
162  void operator()(const size_type iter) const {
163  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
164  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
165 
166  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
167  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
168 
169  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
170 
171  auto vcprop = Kokkos::common_view_alloc_prop(_work);
172  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
173 
174  switch (opType) {
175  case OPERATOR_VALUE : {
176  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
177  Serial<opType>::getValues( output, input, work, _coeffs );
178  break;
179  }
180  case OPERATOR_DIV: {
181  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
182  Serial<opType>::getValues( output, input, work, _coeffs );
183  break;
184  }
185  default: {
186  INTREPID2_TEST_FOR_ABORT( true,
187  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
188 
189  }
190  }
191  }
192  };
193 };
194 }
195 
196 template<typename DeviceType = void,
197  typename outputValueType = double,
198  typename pointValueType = double>
200  : public Basis<DeviceType,outputValueType,pointValueType> {
201  public:
202  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost OrdinalTypeArray1DHost;
203  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost OrdinalTypeArray2DHost;
204  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost OrdinalTypeArray3DHost;
205 
208  Basis_HDIV_TET_In_FEM(const ordinal_type order,
209  const EPointType pointType = POINTTYPE_EQUISPACED);
210 
211 
215 
217 
219 
220  virtual
221  void
222  getValues( /* */ OutputViewType outputValues,
223  const PointViewType inputPoints,
224  const EOperator operatorType = OPERATOR_VALUE) const override {
225 #ifdef HAVE_INTREPID2_DEBUG
226  Intrepid2::getValues_HDIV_Args(outputValues,
227  inputPoints,
228  operatorType,
229  this->getBaseCellTopology(),
230  this->getCardinality() );
231 #endif
232 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
233 Impl::Basis_HDIV_TET_In_FEM::
234 getValues<DeviceType,numPtsPerEval>( outputValues,
235  inputPoints,
236  this->coeffs_,
237  operatorType);
238  }
239 
240  virtual
241  void
242  getDofCoords( ScalarViewType dofCoords ) const override {
243 #ifdef HAVE_INTREPID2_DEBUG
244  // Verify rank of output array.
245  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
246  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
247  // Verify 0th dimension of output array.
248  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
249  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
250  // Verify 1st dimension of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
253 #endif
254  Kokkos::deep_copy(dofCoords, this->dofCoords_);
255  }
256 
257  virtual
258  void
259  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
260 #ifdef HAVE_INTREPID2_DEBUG
261  // Verify rank of output array.
262  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
263  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
264  // Verify 0th dimension of output array.
265  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
266  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
267  // Verify 1st dimension of output array.
268  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
269  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
270 #endif
271  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
272  }
273 
274  void
275  getExpansionCoeffs( ScalarViewType coeffs ) const {
276  // has to be same rank and dimensions
277  Kokkos::deep_copy(coeffs, this->coeffs_);
278  }
279 
280  virtual
281  const char*
282  getName() const override {
283  return "Intrepid2_HDIV_TET_In_FEM";
284  }
285 
286  virtual
287  bool
288  requireOrientation() const override {
289  return true;
290  }
291 
302  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
303 
304  if(subCellDim == 2) {
305  return Teuchos::rcp(new
307  (this->basisDegree_-1, pointType_));
308  }
309  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
310  }
311 
313  getHostBasis() const override{
315  }
316  private:
317 
320  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
321 
324 
325 };
326 
327 }// namespace Intrepid2
328 
330 
331 #endif
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
small utility functions
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
void getValues_HDIV_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 HDIV-conforming FEM basis...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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...
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual bool requireOrientation() const override
True if orientation is required.
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...
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
EPointType
Enumeration of types of point distributions in Intrepid.
EPointType pointType_
type of lattice used for creating the DoF coordinates
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HDIV_TET_In_FEM.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
virtual const char * getName() const override
Returns basis name.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
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.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.