Intrepid2
Intrepid2_HDIV_TRI_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_TRI_IN_FEM_HPP__
50 #define __INTREPID2_HDIV_TRI_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 
88 #define CardinalityHDivTri(order) (order*(order+2))
89 
90 namespace Impl {
91 
96 public:
97  typedef struct Triangle<3> cell_topology_type;
98 
102  template<EOperator opType>
103  struct Serial {
104  template<typename outputValueViewType,
105  typename inputPointViewType,
106  typename workViewType,
107  typename vinvViewType>
108  KOKKOS_INLINE_FUNCTION
109  static void
110  getValues( outputValueViewType outputValues,
111  const inputPointViewType inputPoints,
112  workViewType work,
113  const vinvViewType vinv );
114 
115  KOKKOS_INLINE_FUNCTION
116  static ordinal_type
117  getWorkSizePerPoint(ordinal_type order) {
118  auto cardinality = CardinalityHDivTri(order);
119  switch (opType) {
120  case OPERATOR_GRAD:
121  case OPERATOR_DIV:
122  case OPERATOR_D1:
123  return 5*cardinality;
124  default:
125  return getDkCardinality<opType,2>()*cardinality;
126  }
127  }
128  };
129 
130  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
131  typename outputValueValueType, class ...outputValueProperties,
132  typename inputPointValueType, class ...inputPointProperties,
133  typename vinvValueType, class ...vinvProperties>
134  static void
135  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
136  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
137  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
138  const EOperator operatorType);
139 
143  template<typename outputValueViewType,
144  typename inputPointViewType,
145  typename vinvViewType,
146  typename workViewType,
147  EOperator opType,
148  ordinal_type numPtsEval>
149  struct Functor {
150  outputValueViewType _outputValues;
151  const inputPointViewType _inputPoints;
152  const vinvViewType _coeffs;
153  workViewType _work;
154 
155  KOKKOS_INLINE_FUNCTION
156  Functor( outputValueViewType outputValues_,
157  inputPointViewType inputPoints_,
158  vinvViewType coeffs_,
159  workViewType work_)
160  : _outputValues(outputValues_), _inputPoints(inputPoints_),
161  _coeffs(coeffs_), _work(work_) {}
162 
163  KOKKOS_INLINE_FUNCTION
164  void operator()(const size_type iter) const {
165  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
166  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
167 
168  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
169  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
170 
171  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
172 
173  auto vcprop = Kokkos::common_view_alloc_prop(_work);
174  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
175 
176 
177  switch (opType) {
178  case OPERATOR_VALUE : {
179  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
180  Serial<opType>::getValues( output, input, work, _coeffs );
181  break;
182  }
183  case OPERATOR_DIV: {
184  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
185  Serial<opType>::getValues( output, input, work, _coeffs );
186  break;
187  }
188  default: {
189  INTREPID2_TEST_FOR_ABORT( true,
190  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::Functor) operator is not supported");
191 
192  }
193  }
194  }
195  };
196 };
197 }
198 
199 template<typename ExecSpaceType = void,
200  typename outputValueType = double,
201  typename pointValueType = double>
203  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
204  public:
208 
211  Basis_HDIV_TRI_In_FEM(const ordinal_type order,
212  const EPointType pointType = POINTTYPE_EQUISPACED);
213 
217 
219 
221 
222  virtual
223  void
224  getValues( /* */ OutputViewType outputValues,
225  const PointViewType inputPoints,
226  const EOperator operatorType = OPERATOR_VALUE) const {
227 #ifdef HAVE_INTREPID2_DEBUG
228  Intrepid2::getValues_HDIV_Args(outputValues,
229  inputPoints,
230  operatorType,
231  this->getBaseCellTopology(),
232  this->getCardinality() );
233 #endif
234  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
235  Impl::Basis_HDIV_TRI_In_FEM::
236  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
237  inputPoints,
238  this->coeffs_,
239  operatorType);
240  }
241 
242  virtual
243  void
244  getDofCoords( ScalarViewType dofCoords ) const {
245 #ifdef HAVE_INTREPID2_DEBUG
246  // Verify rank of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
249  // Verify 0th dimension of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
252  // Verify 1st dimension of output array.
253  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
254  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
255 #endif
256  Kokkos::deep_copy(dofCoords, this->dofCoords_);
257  }
258 
259  virtual
260  void
261  getDofCoeffs( ScalarViewType dofCoeffs ) const {
262 #ifdef HAVE_INTREPID2_DEBUG
263  // Verify rank of output array.
264  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
265  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
266  // Verify 0th dimension of output array.
267  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
268  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
269  // Verify 1st dimension of output array.
270  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
271  ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
272 #endif
273  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
274  }
275 
276  void
277  getExpansionCoeffs( ScalarViewType coeffs ) const {
278  // has to be same rank and dimensions
279  Kokkos::deep_copy(coeffs, this->coeffs_);
280  }
281 
282  virtual
283  const char*
284  getName() const {
285  return "Intrepid2_HDIV_TRI_In_FEM";
286  }
287 
288  virtual
289  bool
291  return true;
292  }
293 
294  private:
295 
298  Kokkos::DynRankView<scalarType,ExecSpaceType> coeffs_;
299 
300 };
301 
302 }// namespace Intrepid2
303 
305 
306 #endif
virtual bool requireOrientation() const
True if orientation is required.
See Intrepid2::Basis_HDIV_TRI_In_FEM.
small utility functions
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual const char * getName() const
Returns basis name.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Basis_HDIV_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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.
Definition file for FEM basis functions of degree n for H(div) functions on TRI cells.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.