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( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
131  const EOperator operatorType);
132 
136  template<typename outputValueViewType,
137  typename inputPointViewType,
138  typename vinvViewType,
139  typename workViewType,
140  EOperator opType,
141  ordinal_type numPtsEval>
142  struct Functor {
143  outputValueViewType _outputValues;
144  const inputPointViewType _inputPoints;
145  const vinvViewType _coeffs;
146  workViewType _work;
147 
148  KOKKOS_INLINE_FUNCTION
149  Functor( outputValueViewType outputValues_,
150  inputPointViewType inputPoints_,
151  vinvViewType coeffs_,
152  workViewType work_)
153  : _outputValues(outputValues_), _inputPoints(inputPoints_),
154  _coeffs(coeffs_), _work(work_) {}
155 
156  KOKKOS_INLINE_FUNCTION
157  void operator()(const size_type iter) const {
158  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
159  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
160 
161  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
162  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
163 
164 
165  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
166 
167  auto vcprop = Kokkos::common_view_alloc_prop(_work);
168  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
169 
170  switch (opType) {
171  case OPERATOR_VALUE : {
172  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
173  Serial<opType>::getValues( output, input, work, _coeffs );
174  break;
175  }
176  case OPERATOR_CURL: {
177  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
178  Serial<opType>::getValues( output, input, work, _coeffs );
179  break;
180  }
181  default: {
182  INTREPID2_TEST_FOR_ABORT( true,
183  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
184 
185  }
186  }
187  }
188  };
189 };
190 }
191 
192 template<typename DeviceType = void,
193  typename outputValueType = double,
194  typename pointValueType = double>
196  : public Basis<DeviceType,outputValueType,pointValueType> {
197  public:
198  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost OrdinalTypeArray1DHost;
199  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost OrdinalTypeArray2DHost;
200  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost OrdinalTypeArray3DHost;
201 
204  Basis_HCURL_TRI_In_FEM(const ordinal_type order,
205  const EPointType pointType = POINTTYPE_EQUISPACED);
206 
207 
211 
213 
215 
216  virtual
217  void
218  getValues( OutputViewType outputValues,
219  const PointViewType inputPoints,
220  const EOperator operatorType = OPERATOR_VALUE) const override {
221 #ifdef HAVE_INTREPID2_DEBUG
222  Intrepid2::getValues_HCURL_Args(outputValues,
223  inputPoints,
224  operatorType,
225  this->getBaseCellTopology(),
226  this->getCardinality() );
227 #endif
228  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
229  Impl::Basis_HCURL_TRI_In_FEM::
230  getValues<DeviceType,numPtsPerEval>( outputValues,
231  inputPoints,
232  this->coeffs_,
233  operatorType);
234  }
235 
236  virtual
237  void
238  getDofCoords( ScalarViewType dofCoords ) const override {
239 #ifdef HAVE_INTREPID2_DEBUG
240  // Verify rank of output array.
241  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
242  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
243  // Verify 0th dimension of output array.
244  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
245  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
246  // Verify 1st dimension of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
249 #endif
250  Kokkos::deep_copy(dofCoords, this->dofCoords_);
251  }
252 
253  virtual
254  void
255  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
256 #ifdef HAVE_INTREPID2_DEBUG
257  // Verify rank of output array.
258  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
259  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
260  // Verify 0th dimension of output array.
261  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
262  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
263  // Verify 1st dimension of output array.
264  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
265  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
266 #endif
267  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
268  }
269 
270  void
271  getExpansionCoeffs( ScalarViewType coeffs ) const {
272  // has to be same rank and dimensions
273  Kokkos::deep_copy(coeffs, this->coeffs_);
274  }
275 
276  virtual
277  const char*
278  getName() const override {
279  return "Intrepid2_HCURL_TRI_In_FEM";
280  }
281 
282  virtual
283  bool
284  requireOrientation() const override {
285  return true;
286  }
287 
298  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
299  if(subCellDim == 1) {
300  return Teuchos::rcp(new
302  ( this->basisDegree_ - 1, pointType_) );
303  }
304  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
305  }
306 
308  getHostBasis() const override{
310  }
311  private:
312 
315  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
316 
319 
320 };
321 
322 }// namespace Intrepid2
323 
325 
326 #endif
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
small utility functions
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
void getValues_HCURL_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HCURL-conforming FEM basis...
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 (...
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual bool requireOrientation() const override
True if orientation is required.
ordinal_type getCardinality() const
Returns cardinality of the basis.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
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
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
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.
EPointType
Enumeration of types of point distributions in Intrepid.
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...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
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.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.