Intrepid2
Intrepid2_HDIV_QUAD_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_QUAD_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_QUAD_IN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 namespace Intrepid2 {
57 
58  namespace Impl {
59 
64  public:
65  typedef struct Quadrilateral<4> cell_topology_type;
69  template<EOperator opType>
70  struct Serial {
71  template<typename outputValueViewType,
72  typename inputPointViewType,
73  typename workViewType,
74  typename vinvViewType>
75  KOKKOS_INLINE_FUNCTION
76  static void
77  getValues( outputValueViewType outputValues,
78  const inputPointViewType inputPoints,
79  workViewType work,
80  const vinvViewType vinvLine,
81  const vinvViewType vinvBubble );
82 
83  KOKKOS_INLINE_FUNCTION
84  static ordinal_type
85  getWorkSizePerPoint(ordinal_type order) {
86  return 2*getPnCardinality<1>(order)+getPnCardinality<1>(order-1);
87  }
88  };
89 
90  template<typename DeviceType, ordinal_type numPtsPerEval,
91  typename outputValueValueType, class ...outputValueProperties,
92  typename inputPointValueType, class ...inputPointProperties,
93  typename vinvValueType, class ...vinvProperties>
94  static void
95  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
96  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
97  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
98  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
99  const EOperator operatorType );
100 
104  template<typename outputValueViewType,
105  typename inputPointViewType,
106  typename vinvViewType,
107  typename workViewType,
108  EOperator opType,
109  ordinal_type numPtsEval>
110  struct Functor {
111  outputValueViewType _outputValues;
112  const inputPointViewType _inputPoints;
113  const vinvViewType _vinvLine;
114  const vinvViewType _vinvBubble;
115  workViewType _work;
116 
117  KOKKOS_INLINE_FUNCTION
118  Functor( outputValueViewType outputValues_,
119  inputPointViewType inputPoints_,
120  vinvViewType vinvLine_,
121  vinvViewType vinvBubble_,
122  workViewType work_)
123  : _outputValues(outputValues_), _inputPoints(inputPoints_),
124  _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
125 
126  KOKKOS_INLINE_FUNCTION
127  void operator()(const size_type iter) const {
128  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
129  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
130 
131  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
132  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
133 
134  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
135 
136  auto vcprop = Kokkos::common_view_alloc_prop(_work);
137  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
138 
139  switch (opType) {
140  case OPERATOR_VALUE : {
141  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
142  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
143  break;
144  }
145  case OPERATOR_DIV : {
146  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
147  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
148  break;
149  }
150  default: {
151  INTREPID2_TEST_FOR_ABORT( true,
152  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::Functor) operator is not supported");
153 
154  }
155  }
156  }
157  };
158  };
159  }
160 
165  template<typename DeviceType = void,
166  typename outputValueType = double,
167  typename pointValueType = double>
169  : public Basis<DeviceType,outputValueType,pointValueType> {
170  public:
171  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
172  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
173  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
174 
177  Basis_HDIV_QUAD_In_FEM(const ordinal_type order,
178  const EPointType pointType = POINTTYPE_EQUISPACED);
179 
183 
185 
186  virtual
187  void
188  getValues( OutputViewType outputValues,
189  const PointViewType inputPoints,
190  const EOperator operatorType = OPERATOR_VALUE ) const override {
191 #ifdef HAVE_INTREPID2_DEBUG
192  Intrepid2::getValues_HDIV_Args(outputValues,
193  inputPoints,
194  operatorType,
195  this->getBaseCellTopology(),
196  this->getCardinality() );
197 #endif
198  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
199  Impl::Basis_HDIV_QUAD_In_FEM::
200  getValues<DeviceType,numPtsPerEval>( outputValues,
201  inputPoints,
202  this->vinvLine_,
203  this->vinvBubble_,
204  operatorType );
205  }
206  virtual
207  void
208  getDofCoords( ScalarViewType dofCoords ) const override {
209 #ifdef HAVE_INTREPID2_DEBUG
210  // Verify rank of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
213  // Verify 0th dimension of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
216  // Verify 1st dimension of output array.
217  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
218  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
219 #endif
220  Kokkos::deep_copy(dofCoords, this->dofCoords_);
221  }
222 
223  virtual
224  void
225  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
226 #ifdef HAVE_INTREPID2_DEBUG
227  // Verify rank of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
230  // Verify 0th dimension of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
233  // Verify 1st dimension of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
236 #endif
237  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
238  }
239 
240  virtual
241  const char*
242  getName() const override {
243  return "Intrepid2_HDIV_QUAD_In_FEM";
244  }
245 
246  virtual
247  bool
248  requireOrientation() const override {
249  return true;
250  }
251 
261  BasisPtr<DeviceType,outputValueType,pointValueType>
262  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
263  if(subCellDim == 1)
264  return Teuchos::rcp( new
266  ( this->basisDegree_ - 1, POINTTYPE_GAUSS ) );
267 
268  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
269  }
270 
271  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
272  getHostBasis() const override{
274  }
275  private:
276 
278  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinvLine_, vinvBubble_;
279  EPointType pointType_;
280  };
281 
282 }// namespace Intrepid2
283 
285 
286 #endif
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(div) functions on QUAD cells.
small utility functions
Implementation of the default H(div)-compatible FEM basis on Quadrilateral cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
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...
virtual bool requireOrientation() const override
True if orientation is required.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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...
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual const char * getName() const override
Returns basis name.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Basis_HDIV_QUAD_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HDIV_QUAD_In_FEM.
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.