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