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 DeviceType, ordinal_type numPtsPerEval,
87  typename outputValueValueType, class ...outputValueProperties,
88  typename inputPointValueType, class ...inputPointProperties,
89  typename vinvValueType, class ...vinvProperties>
90  static void
91  getValues( const typename DeviceType::execution_space& space,
92  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
93  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
94  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
95  const EOperator operatorType);
96 
97 
101  template<typename outputValueViewType,
102  typename inputPointViewType,
103  typename vinvViewType,
104  typename workViewType,
105  EOperator opType,
106  ordinal_type numPtsEval>
107  struct Functor {
108  outputValueViewType _outputValues;
109  const inputPointViewType _inputPoints;
110  const vinvViewType _vinv;
111  workViewType _work;
112  const ordinal_type _opDn;
113 
114  KOKKOS_INLINE_FUNCTION
115  Functor( outputValueViewType outputValues_,
116  inputPointViewType inputPoints_,
117  vinvViewType vinv_,
118  workViewType work_,
119  const ordinal_type opDn_ = 0 )
120  : _outputValues(outputValues_), _inputPoints(inputPoints_),
121  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
122 
123  KOKKOS_INLINE_FUNCTION
124  void operator()(const size_type iter) const {
125  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
126  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
127 
128  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
129  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
130 
131  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
132 
133  auto vcprop = Kokkos::common_view_alloc_prop(_work);
134  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
135 
136  switch (opType) {
137  case OPERATOR_VALUE : {
138  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
139  Serial<opType>::getValues( output, input, work, _vinv );
140  break;
141  }
142  case OPERATOR_CURL :
143  case OPERATOR_Dn : {
144  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
145  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
146  break;
147  }
148  default: {
149  INTREPID2_TEST_FOR_ABORT( true,
150  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::Functor) operator is not supported");
151 
152  }
153  }
154  }
155  };
156  };
157  }
158 
164  template<typename DeviceType = void,
165  typename outputValueType = double,
166  typename pointValueType = double>
168  : public Basis<DeviceType,outputValueType,pointValueType> {
169  public:
171  using typename BasisBase::ExecutionSpace;
172  using typename BasisBase::OrdinalTypeArray1DHost;
173  using typename BasisBase::OrdinalTypeArray2DHost;
174  using typename BasisBase::OrdinalTypeArray3DHost;
175 
176  using typename BasisBase::OutputViewType;
177  using typename BasisBase::PointViewType;
178  using typename BasisBase::ScalarViewType;
179 
180  private:
182  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
183 
185  EPointType pointType_;
186 
187  public:
190  Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order,
191  const EPointType pointType = POINTTYPE_EQUISPACED);
192 
193  using BasisBase::getValues;
194 
195  virtual
196  void
197  getValues( const ExecutionSpace& space,
198  OutputViewType outputValues,
199  const PointViewType inputPoints,
200  const EOperator operatorType = OPERATOR_VALUE ) const override {
201 #ifdef HAVE_INTREPID2_DEBUG
202  Intrepid2::getValues_HGRAD_Args(outputValues,
203  inputPoints,
204  operatorType,
205  this->getBaseCellTopology(),
206  this->getCardinality() );
207 #endif
208  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
209  Impl::Basis_HGRAD_QUAD_Cn_FEM::
210  getValues<DeviceType,numPtsPerEval>(space,
211  outputValues,
212  inputPoints,
213  this->vinv_,
214  operatorType);
215  }
216 
217  virtual
218  void
219  getDofCoords( ScalarViewType dofCoords ) const override {
220 #ifdef HAVE_INTREPID2_DEBUG
221  // Verify rank of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
224  // Verify 0th dimension of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
227  // Verify 1st dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
230 #endif
231  Kokkos::deep_copy(dofCoords, this->dofCoords_);
232  }
233 
234  virtual
235  void
236  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
237 #ifdef HAVE_INTREPID2_DEBUG
238  // Verify rank of output array.
239  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
240  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
241  // Verify 0th dimension of output array.
242  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
243  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
244 #endif
245  Kokkos::deep_copy(dofCoeffs, 1.0);
246  }
247 
248  virtual
249  const char*
250  getName() const override {
251  return "Intrepid2_HGRAD_QUAD_Cn_FEM";
252  }
253 
254  virtual
255  bool
256  requireOrientation() const override {
257  return (this->basisDegree_ > 2);
258  }
259 
260  Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
261  getVandermondeInverse() const {
262  return vinv_;
263  }
264 
265  ordinal_type
266  getWorkSizePerPoint(const EOperator operatorType) const {
267  return 3*getPnCardinality<1>(this->basisDegree_);
268  }
269 
278  BasisPtr<DeviceType,outputValueType,pointValueType>
279  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
280  if(subCellDim == 1) {
281  return Teuchos::rcp(new
283  (this->basisDegree_,pointType_));
284  }
285  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
286  }
287 
288  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
289  getHostBasis() const override{
291  }
292  };
293 
294 }// namespace Intrepid2
295 
297 
298 #endif
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
small utility functions
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...
ordinal_type getCardinality() const
Returns cardinality of the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
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...
Basis_HGRAD_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
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, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree n for H(grad) functions on QUAD cells...
See Intrepid2::Basis_HGRAD_QUAD_Cn_FEM work is a rank 1 view having the same value_type of inputPoint...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.