Intrepid2
Intrepid2_HCURL_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_HCURL_TRI_IN_FEM_HPP__
50 #define __INTREPID2_HCURL_TRI_IN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
55 
56 #include "Intrepid2_PointTools.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 namespace Intrepid2 {
60 
81 #define CardinalityHCurlTri(order) (order*(order+2))
82 
83 namespace Impl {
84 
89 public:
90  typedef struct Triangle<3> cell_topology_type;
91 
95  template<EOperator opType>
96  struct Serial {
97  template<typename outputValueViewType,
98  typename inputPointViewType,
99  typename workViewType,
100  typename vinvViewType>
101  KOKKOS_INLINE_FUNCTION
102  static void
103  getValues( outputValueViewType outputValues,
104  const inputPointViewType inputPoints,
105  workViewType work,
106  const vinvViewType vinv );
107 
108  KOKKOS_INLINE_FUNCTION
109  static ordinal_type
110  getWorkSizePerPoint(ordinal_type order) {
111  auto cardinality = CardinalityHCurlTri(order);
112  switch (opType) {
113  case OPERATOR_GRAD:
114  case OPERATOR_CURL:
115  case OPERATOR_D1:
116  return 5*cardinality;
117  default:
118  return getDkCardinality<opType,2>()*cardinality;
119  }
120  }
121  };
122 
123  template<typename DeviceType, ordinal_type numPtsPerEval,
124  typename outputValueValueType, class ...outputValueProperties,
125  typename inputPointValueType, class ...inputPointProperties,
126  typename vinvValueType, class ...vinvProperties>
127  static void
128  getValues(
129  const typename DeviceType::execution_space& space,
130  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
131  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
132  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
133  const EOperator operatorType);
134 
138  template<typename outputValueViewType,
139  typename inputPointViewType,
140  typename vinvViewType,
141  typename workViewType,
142  EOperator opType,
143  ordinal_type numPtsEval>
144  struct Functor {
145  outputValueViewType _outputValues;
146  const inputPointViewType _inputPoints;
147  const vinvViewType _coeffs;
148  workViewType _work;
149 
150  KOKKOS_INLINE_FUNCTION
151  Functor( outputValueViewType outputValues_,
152  inputPointViewType inputPoints_,
153  vinvViewType coeffs_,
154  workViewType work_)
155  : _outputValues(outputValues_), _inputPoints(inputPoints_),
156  _coeffs(coeffs_), _work(work_) {}
157 
158  KOKKOS_INLINE_FUNCTION
159  void operator()(const size_type iter) const {
160  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
161  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
162 
163  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
164  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
165 
166 
167  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
168 
169  auto vcprop = Kokkos::common_view_alloc_prop(_work);
170  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
171 
172  switch (opType) {
173  case OPERATOR_VALUE : {
174  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
175  Serial<opType>::getValues( output, input, work, _coeffs );
176  break;
177  }
178  case OPERATOR_CURL: {
179  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
180  Serial<opType>::getValues( output, input, work, _coeffs );
181  break;
182  }
183  default: {
184  INTREPID2_TEST_FOR_ABORT( true,
185  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
186 
187  }
188  }
189  }
190  };
191 };
192 }
193 
194 template<typename DeviceType = void,
195  typename outputValueType = double,
196  typename pointValueType = double>
198  : public Basis<DeviceType,outputValueType,pointValueType> {
199  public:
201  using typename BasisBase::ExecutionSpace;
202  using typename BasisBase::OrdinalTypeArray1DHost;
203  using typename BasisBase::OrdinalTypeArray2DHost;
204  using typename BasisBase::OrdinalTypeArray3DHost;
205 
208  Basis_HCURL_TRI_In_FEM(const ordinal_type order,
209  const EPointType pointType = POINTTYPE_EQUISPACED);
210 
212 
213  using typename BasisBase::OutputViewType;
214  using typename BasisBase::PointViewType;
215  using typename BasisBase::ScalarViewType;
216 
217  using typename BasisBase::scalarType;
218 
219  using BasisBase::getValues;
220 
221  virtual
222  void
224  const ExecutionSpace& space,
225  OutputViewType outputValues,
226  const PointViewType inputPoints,
227  const EOperator operatorType = OPERATOR_VALUE) const override {
228 #ifdef HAVE_INTREPID2_DEBUG
229  Intrepid2::getValues_HCURL_Args(outputValues,
230  inputPoints,
231  operatorType,
232  this->getBaseCellTopology(),
233  this->getCardinality() );
234 #endif
235  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
236  Impl::Basis_HCURL_TRI_In_FEM::
237  getValues<DeviceType,numPtsPerEval>(
238  space,
239  outputValues,
240  inputPoints,
241  this->coeffs_,
242  operatorType);
243  }
244 
245  virtual
246  void
247  getDofCoords( ScalarViewType dofCoords ) const override {
248 #ifdef HAVE_INTREPID2_DEBUG
249  // Verify rank of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
252  // Verify 0th dimension of output array.
253  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
254  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
255  // Verify 1st dimension of output array.
256  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
257  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
258 #endif
259  Kokkos::deep_copy(dofCoords, this->dofCoords_);
260  }
261 
262  virtual
263  void
264  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
265 #ifdef HAVE_INTREPID2_DEBUG
266  // Verify rank of output array.
267  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
268  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
269  // Verify 0th dimension of output array.
270  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
271  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
272  // Verify 1st dimension of output array.
273  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
274  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
275 #endif
276  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
277  }
278 
279  void
280  getExpansionCoeffs( ScalarViewType coeffs ) const {
281  // has to be same rank and dimensions
282  Kokkos::deep_copy(coeffs, this->coeffs_);
283  }
284 
285  virtual
286  const char*
287  getName() const override {
288  return "Intrepid2_HCURL_TRI_In_FEM";
289  }
290 
291  virtual
292  bool
293  requireOrientation() const override {
294  return true;
295  }
296 
306  BasisPtr<DeviceType,outputValueType,pointValueType>
307  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
308  if(subCellDim == 1) {
309  return Teuchos::rcp(new
311  ( this->basisDegree_ - 1, pointType_) );
312  }
313  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
314  }
315 
316  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
317  getHostBasis() const override{
319  }
320  private:
321 
324  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
325 
327  EPointType pointType_;
328 
329 };
330 
331 }// namespace Intrepid2
332 
334 
335 #endif
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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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...
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.
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
virtual bool requireOrientation() const override
True if orientation is required.
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 basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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 HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
Definition file for FEM basis functions of degree n for H(curl) functions on TRI. ...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual const char * getName() const override
Returns basis name.
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.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.