Intrepid2
Intrepid2_HGRAD_QUAD_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 
49 #ifndef __INTREPID2_HGRAD_QUAD_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_QUAD_CN_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;
65 
71  template<EOperator opType>
72  struct Serial {
73  template<typename outputValueViewType,
74  typename inputPointViewType,
75  typename workViewType,
76  typename vinvViewType>
77  KOKKOS_INLINE_FUNCTION
78  static void
79  getValues( outputValueViewType outputValues,
80  const inputPointViewType inputPoints,
81  workViewType work,
82  const vinvViewType vinv,
83  const ordinal_type operatorDn = 0 );
84  };
85 
86  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
87  typename outputValueValueType, class ...outputValueProperties,
88  typename inputPointValueType, class ...inputPointProperties,
89  typename vinvValueType, class ...vinvProperties>
90  static void
91  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
94  const EOperator operatorType);
95 
96 
100  template<typename outputValueViewType,
101  typename inputPointViewType,
102  typename vinvViewType,
103  typename workViewType,
104  EOperator opType,
105  ordinal_type numPtsEval>
106  struct Functor {
107  outputValueViewType _outputValues;
108  const inputPointViewType _inputPoints;
109  const vinvViewType _vinv;
110  workViewType _work;
111  const ordinal_type _opDn;
112 
113  KOKKOS_INLINE_FUNCTION
114  Functor( outputValueViewType outputValues_,
115  inputPointViewType inputPoints_,
116  vinvViewType vinv_,
117  workViewType work_,
118  const ordinal_type opDn_ = 0 )
119  : _outputValues(outputValues_), _inputPoints(inputPoints_),
120  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
121 
122  KOKKOS_INLINE_FUNCTION
123  void operator()(const size_type iter) const {
124  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
125  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
126 
127  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
128  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
129 
130  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
131 
132  auto vcprop = Kokkos::common_view_alloc_prop(_work);
133  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
134 
135  switch (opType) {
136  case OPERATOR_VALUE : {
137  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
138  Serial<opType>::getValues( output, input, work, _vinv );
139  break;
140  }
141  case OPERATOR_CURL :
142  case OPERATOR_Dn : {
143  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
144  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
145  break;
146  }
147  default: {
148  INTREPID2_TEST_FOR_ABORT( true,
149  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
150 
151  }
152  }
153  }
154  };
155  };
156  }
157 
163  template<typename ExecSpaceType = void,
164  typename outputValueType = double,
165  typename pointValueType = double>
167  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
168  public:
172 
176 
177  private:
179  Kokkos::DynRankView<typename ScalarViewType::value_type,ExecSpaceType> vinv_;
180 
181  public:
184  Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order,
185  const EPointType pointType = POINTTYPE_EQUISPACED);
186 
188 
189  virtual
190  void
191  getValues( OutputViewType outputValues,
192  const PointViewType inputPoints,
193  const EOperator operatorType = OPERATOR_VALUE ) const {
194 #ifdef HAVE_INTREPID2_DEBUG
195  Intrepid2::getValues_HGRAD_Args(outputValues,
196  inputPoints,
197  operatorType,
198  this->getBaseCellTopology(),
199  this->getCardinality() );
200 #endif
201  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
202  Impl::Basis_HGRAD_QUAD_Cn_FEM::
203  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
204  inputPoints,
205  this->vinv_,
206  operatorType);
207  }
208 
209  virtual
210  void
211  getDofCoords( ScalarViewType dofCoords ) const {
212 #ifdef HAVE_INTREPID2_DEBUG
213  // Verify rank of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
216  // Verify 0th dimension of output array.
217  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
219  // Verify 1st dimension of output array.
220  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
222 #endif
223  Kokkos::deep_copy(dofCoords, this->dofCoords_);
224  }
225 
226  virtual
227  void
228  getDofCoeffs( ScalarViewType dofCoeffs ) const {
229 #ifdef HAVE_INTREPID2_DEBUG
230  // Verify rank of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
233  // Verify 0th dimension of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
236 #endif
237  Kokkos::deep_copy(dofCoeffs, 1.0);
238  }
239 
240  virtual
241  const char*
242  getName() const {
243  return "Intrepid2_HGRAD_QUAD_Cn_FEM";
244  }
245 
246  virtual
247  bool
249  return (this->basisDegree_ > 2);
250  }
251 
252  Kokkos::DynRankView<typename ScalarViewType::const_value_type,ExecSpaceType>
253  getVandermondeInverse() const {
254  return vinv_;
255  }
256 
257  ordinal_type
258  getWorkSizePerPoint(const EOperator operatorType) const {
259  return 3*getPnCardinality<1>(this->basisDegree_);
260  }
261 
262  };
263 
264 }// namespace Intrepid2
265 
267 
268 #endif
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
small utility functions
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.
virtual bool requireOrientation() const
True if orientation is required.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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 const char * getName() const
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree n for H(grad) functions on QUAD cells...
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM work is a rank 1 view having the same value_type of inputPoint...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
Header file for the abstract base class Intrepid2::Basis.