Intrepid2
Intrepid2_HDIV_TRI_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_TRI_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TRI_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 
89 #define CardinalityHDivTri(order) (order*(order+2))
90 
91 namespace Impl {
92 
97 public:
98  typedef struct Triangle<3> cell_topology_type;
99 
103  template<EOperator opType>
104  struct Serial {
105  template<typename outputValueViewType,
106  typename inputPointViewType,
107  typename workViewType,
108  typename vinvViewType>
109  KOKKOS_INLINE_FUNCTION
110  static void
111  getValues( outputValueViewType outputValues,
112  const inputPointViewType inputPoints,
113  workViewType work,
114  const vinvViewType vinv );
115 
116  KOKKOS_INLINE_FUNCTION
117  static ordinal_type
118  getWorkSizePerPoint(ordinal_type order) {
119  auto cardinality = CardinalityHDivTri(order);
120  switch (opType) {
121  case OPERATOR_GRAD:
122  case OPERATOR_DIV:
123  case OPERATOR_D1:
124  return 5*cardinality;
125  default:
126  return getDkCardinality<opType,2>()*cardinality;
127  }
128  }
129  };
130 
131  template<typename DeviceType, ordinal_type numPtsPerEval,
132  typename outputValueValueType, class ...outputValueProperties,
133  typename inputPointValueType, class ...inputPointProperties,
134  typename vinvValueType, class ...vinvProperties>
135  static void
136  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
137  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
138  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
139  const EOperator operatorType);
140 
144  template<typename outputValueViewType,
145  typename inputPointViewType,
146  typename vinvViewType,
147  typename workViewType,
148  EOperator opType,
149  ordinal_type numPtsEval>
150  struct Functor {
151  outputValueViewType _outputValues;
152  const inputPointViewType _inputPoints;
153  const vinvViewType _coeffs;
154  workViewType _work;
155 
156  KOKKOS_INLINE_FUNCTION
157  Functor( outputValueViewType outputValues_,
158  inputPointViewType inputPoints_,
159  vinvViewType coeffs_,
160  workViewType work_)
161  : _outputValues(outputValues_), _inputPoints(inputPoints_),
162  _coeffs(coeffs_), _work(work_) {}
163 
164  KOKKOS_INLINE_FUNCTION
165  void operator()(const size_type iter) const {
166  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
167  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
168 
169  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
170  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
171 
172  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
173 
174  auto vcprop = Kokkos::common_view_alloc_prop(_work);
175  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
176 
177 
178  switch (opType) {
179  case OPERATOR_VALUE : {
180  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
181  Serial<opType>::getValues( output, input, work, _coeffs );
182  break;
183  }
184  case OPERATOR_DIV: {
185  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
186  Serial<opType>::getValues( output, input, work, _coeffs );
187  break;
188  }
189  default: {
190  INTREPID2_TEST_FOR_ABORT( true,
191  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
192 
193  }
194  }
195  }
196  };
197 };
198 }
199 
200 template<typename DeviceType = void,
201  typename outputValueType = double,
202  typename pointValueType = double>
204  : public Basis<DeviceType,outputValueType,pointValueType> {
205  public:
206  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost OrdinalTypeArray1DHost;
207  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost OrdinalTypeArray2DHost;
208  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost OrdinalTypeArray3DHost;
209 
212  Basis_HDIV_TRI_In_FEM(const ordinal_type order,
213  const EPointType pointType = POINTTYPE_EQUISPACED);
214 
218 
220 
222 
223  virtual
224  void
225  getValues( /* */ OutputViewType outputValues,
226  const PointViewType inputPoints,
227  const EOperator operatorType = OPERATOR_VALUE) const override {
228 #ifdef HAVE_INTREPID2_DEBUG
229  Intrepid2::getValues_HDIV_Args(outputValues,
230  inputPoints,
231  operatorType,
232  this->getBaseCellTopology(),
233  this->getCardinality() );
234 #endif
235  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
236  Impl::Basis_HDIV_TRI_In_FEM::
237  getValues<DeviceType,numPtsPerEval>( outputValues,
238  inputPoints,
239  this->coeffs_,
240  operatorType);
241  }
242 
243  virtual
244  void
245  getDofCoords( ScalarViewType dofCoords ) const override {
246 #ifdef HAVE_INTREPID2_DEBUG
247  // Verify rank of output array.
248  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
249  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
250  // Verify 0th dimension of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
253  // Verify 1st dimension of output array.
254  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
255  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
256 #endif
257  Kokkos::deep_copy(dofCoords, this->dofCoords_);
258  }
259 
260  virtual
261  void
262  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
263 #ifdef HAVE_INTREPID2_DEBUG
264  // Verify rank of output array.
265  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
266  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
267  // Verify 0th dimension of output array.
268  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
269  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
270  // Verify 1st dimension of output array.
271  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
272  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
273 #endif
274  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
275  }
276 
277  void
278  getExpansionCoeffs( ScalarViewType coeffs ) const {
279  // has to be same rank and dimensions
280  Kokkos::deep_copy(coeffs, this->coeffs_);
281  }
282 
283  virtual
284  const char*
285  getName() const override {
286  return "Intrepid2_HDIV_TRI_In_FEM";
287  }
288 
289  virtual
290  bool
291  requireOrientation() const override {
292  return true;
293  }
294 
305  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
306  if(subCellDim == 1) {
307  return Teuchos::rcp(new
309  (this->basisDegree_-1, pointType_));
310  }
311  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
312  }
313 
315  getHostBasis() const override{
317  }
318  private:
319 
322  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
323 
325  EPointType pointType_;
326 
327 };
328 
329 }// namespace Intrepid2
330 
332 
333 #endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Definition file for FEM basis functions of degree n for H(div) functions on TRI cells.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
virtual const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Basis_HDIV_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual bool requireOrientation() const override
True if orientation is required.
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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions