Intrepid2
Intrepid2_HDIV_TET_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_TET_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TET_IN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 namespace Intrepid2 {
59 
86 #define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
87 
88 namespace Impl {
89 
94 public:
95 
99  template<EOperator opType>
100  struct Serial {
101  template<typename outputValueViewType,
102  typename inputPointViewType,
103  typename workViewType,
104  typename vinvViewType>
105  KOKKOS_INLINE_FUNCTION
106  static void
107  getValues( outputValueViewType outputValues,
108  const inputPointViewType inputPoints,
109  workViewType work,
110  const vinvViewType vinv );
111 
112 
113  KOKKOS_INLINE_FUNCTION
114  static ordinal_type
115  getWorkSizePerPoint(ordinal_type order) {
116  auto cardinality = CardinalityHDivTet(order);
117  switch (opType) {
118  case OPERATOR_DIV:
119  case OPERATOR_D1:
120  return 7*cardinality;
121  default:
122  return getDkCardinality<opType,3>()*cardinality;
123  }
124  }
125  };
126 
127  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
128  typename outputValueValueType, class ...outputValueProperties,
129  typename inputPointValueType, class ...inputPointProperties,
130  typename vinvValueType, class ...vinvProperties>
131  static void
132  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
133  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
134  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
135  const EOperator operatorType);
136 
140  template<typename outputValueViewType,
141  typename inputPointViewType,
142  typename vinvViewType,
143  typename workViewType,
144  EOperator opType,
145  ordinal_type numPtsEval>
146  struct Functor {
147  outputValueViewType _outputValues;
148  const inputPointViewType _inputPoints;
149  const vinvViewType _coeffs;
150  workViewType _work;
151 
152  KOKKOS_INLINE_FUNCTION
153  Functor( outputValueViewType outputValues_,
154  inputPointViewType inputPoints_,
155  vinvViewType coeffs_,
156  workViewType work_)
157  : _outputValues(outputValues_), _inputPoints(inputPoints_),
158  _coeffs(coeffs_), _work(work_) {}
159 
160  KOKKOS_INLINE_FUNCTION
161  void operator()(const size_type iter) const {
162  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
163  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
164 
165  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
166  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
167 
168  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
169 
170  auto vcprop = Kokkos::common_view_alloc_prop(_work);
171  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
172 
173  switch (opType) {
174  case OPERATOR_VALUE : {
175  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
176  Serial<opType>::getValues( output, input, work, _coeffs );
177  break;
178  }
179  case OPERATOR_DIV: {
180  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
181  Serial<opType>::getValues( output, input, work, _coeffs );
182  break;
183  }
184  default: {
185  INTREPID2_TEST_FOR_ABORT( true,
186  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
187 
188  }
189  }
190  }
191  };
192 };
193 }
194 
195 template<typename ExecSpaceType = void,
196  typename outputValueType = double,
197  typename pointValueType = double>
199  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
200  public:
204 
207  Basis_HDIV_TET_In_FEM(const ordinal_type order,
208  const EPointType pointType = POINTTYPE_EQUISPACED);
209 
210 
215 
217 
218  virtual
219  void
220  getValues( /* */ outputViewType outputValues,
221  const pointViewType inputPoints,
222  const EOperator operatorType = OPERATOR_VALUE) const {
223 #ifdef HAVE_INTREPID2_DEBUG
224  Intrepid2::getValues_HDIV_Args(outputValues,
225  inputPoints,
226  operatorType,
227  this->getBaseCellTopology(),
228  this->getCardinality() );
229 #endif
230 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
231 Impl::Basis_HDIV_TET_In_FEM::
232 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
233  inputPoints,
234  this->coeffs_,
235  operatorType);
236  }
237 
238  virtual
239  void
240  getDofCoords( scalarViewType dofCoords ) const {
241 #ifdef HAVE_INTREPID2_DEBUG
242  // Verify rank of output array.
243  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
244  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
245  // Verify 0th dimension of output array.
246  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
247  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
248  // Verify 1st dimension of output array.
249  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
250  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
251 #endif
252  Kokkos::deep_copy(dofCoords, this->dofCoords_);
253  }
254 
255  virtual
256  void
257  getDofCoeffs( scalarViewType dofCoeffs ) const {
258 #ifdef HAVE_INTREPID2_DEBUG
259  // Verify rank of output array.
260  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
261  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
262  // Verify 0th dimension of output array.
263  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
264  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
265  // Verify 1st dimension of output array.
266  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
267  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
268 #endif
269  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
270  }
271 
272  void
273  getExpansionCoeffs( scalarViewType coeffs ) const {
274  // has to be same rank and dimensions
275  Kokkos::deep_copy(coeffs, this->coeffs_);
276  }
277 
278  virtual
279  const char*
280  getName() const {
281  return "Intrepid2_HDIV_TET_In_FEM";
282  }
283 
284  virtual
285  bool
287  return true;
288  }
289 
290  private:
291 
294  Kokkos::DynRankView<scalarType,ExecSpaceType> coeffs_;
295 
296 };
297 
298 }// namespace Intrepid2
299 
301 
302 #endif
small utility functions
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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 void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
virtual const char * getName() const
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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_TET_In_FEM.
virtual bool requireOrientation() const
True if orientation is required.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.