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