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"
54 
55 namespace Intrepid2 {
56 
57  namespace Impl {
58 
63  public:
64  typedef struct Quadrilateral<4> cell_topology_type;
68  template<EOperator opType>
69  struct Serial {
70  template<typename outputValueViewType,
71  typename inputPointViewType,
72  typename workViewType,
73  typename vinvViewType>
74  KOKKOS_INLINE_FUNCTION
75  static void
76  getValues( outputValueViewType outputValues,
77  const inputPointViewType inputPoints,
78  workViewType work,
79  const vinvViewType vinvLine,
80  const vinvViewType vinvBubble );
81 
82  KOKKOS_INLINE_FUNCTION
83  static ordinal_type
84  getWorkSizePerPoint(ordinal_type order) {
85  return 2*getPnCardinality<1>(order)+getPnCardinality<1>(order-1);
86  }
87  };
88 
89  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
90  typename outputValueValueType, class ...outputValueProperties,
91  typename inputPointValueType, class ...inputPointProperties,
92  typename vinvValueType, class ...vinvProperties>
93  static void
94  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
95  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
96  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
97  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
98  const EOperator operatorType );
99 
103  template<typename outputValueViewType,
104  typename inputPointViewType,
105  typename vinvViewType,
106  typename workViewType,
107  EOperator opType,
108  ordinal_type numPtsEval>
109  struct Functor {
110  outputValueViewType _outputValues;
111  const inputPointViewType _inputPoints;
112  const vinvViewType _vinvLine;
113  const vinvViewType _vinvBubble;
114  workViewType _work;
115 
116  KOKKOS_INLINE_FUNCTION
117  Functor( outputValueViewType outputValues_,
118  inputPointViewType inputPoints_,
119  vinvViewType vinvLine_,
120  vinvViewType vinvBubble_,
121  workViewType work_)
122  : _outputValues(outputValues_), _inputPoints(inputPoints_),
123  _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
124 
125  KOKKOS_INLINE_FUNCTION
126  void operator()(const size_type iter) const {
127  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
128  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
129 
130  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
131  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
132 
133  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
134 
135  auto vcprop = Kokkos::common_view_alloc_prop(_work);
136  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
137 
138  switch (opType) {
139  case OPERATOR_VALUE : {
140  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
141  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
142  break;
143  }
144  case OPERATOR_DIV : {
145  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
146  Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
147  break;
148  }
149  default: {
150  INTREPID2_TEST_FOR_ABORT( true,
151  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::Functor) operator is not supported");
152 
153  }
154  }
155  }
156  };
157  };
158  }
159 
164  template<typename ExecSpaceType = void,
165  typename outputValueType = double,
166  typename pointValueType = double>
168  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
169  public:
173 
176  Basis_HDIV_QUAD_In_FEM(const ordinal_type order,
177  const EPointType pointType = POINTTYPE_EQUISPACED);
178 
182 
184 
185  virtual
186  void
187  getValues( OutputViewType outputValues,
188  const PointViewType inputPoints,
189  const EOperator operatorType = OPERATOR_VALUE ) const {
190 #ifdef HAVE_INTREPID2_DEBUG
191  Intrepid2::getValues_HDIV_Args(outputValues,
192  inputPoints,
193  operatorType,
194  this->getBaseCellTopology(),
195  this->getCardinality() );
196 #endif
197  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
198  Impl::Basis_HDIV_QUAD_In_FEM::
199  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
200  inputPoints,
201  this->vinvLine_,
202  this->vinvBubble_,
203  operatorType );
204  }
205  virtual
206  void
207  getDofCoords( ScalarViewType dofCoords ) const {
208 #ifdef HAVE_INTREPID2_DEBUG
209  // Verify rank of output array.
210  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
211  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
212  // Verify 0th dimension of output array.
213  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
214  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
215  // Verify 1st dimension of output array.
216  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
217  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
218 #endif
219  Kokkos::deep_copy(dofCoords, this->dofCoords_);
220  }
221 
222  virtual
223  void
224  getDofCoeffs( ScalarViewType dofCoeffs ) const {
225 #ifdef HAVE_INTREPID2_DEBUG
226  // Verify rank of output array.
227  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
228  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
229  // Verify 0th dimension of output array.
230  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
231  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
232  // Verify 1st dimension of output array.
233  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
234  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
235 #endif
236  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
237  }
238 
239  virtual
240  const char*
241  getName() const {
242  return "Intrepid2_HDIV_QUAD_In_FEM";
243  }
244 
245  virtual
246  bool
248  return true;
249  }
250 
251  private:
252 
254  Kokkos::DynRankView<typename ScalarViewType::value_type,ExecSpaceType> vinvLine_, vinvBubble_;
255  };
256 
257 }// namespace Intrepid2
258 
260 
261 #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.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
virtual const char * getName() const
Returns basis name.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecSpaceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HDIV_QUAD_In_FEM.
Basis_HDIV_QUAD_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Header file for the abstract base class Intrepid2::Basis.