Intrepid2
Intrepid2_HCURL_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_HCURL_TET_IN_FEM_HPP__
50 #define __INTREPID2_HCURL_TET_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 
93 #define CardinalityHCurlTet(order) (order*(order+2)*(order+3)/2)
94 
95 namespace Impl {
96 
101 public:
102  typedef struct Tetrahedron<4> cell_topology_type;
106  template<EOperator opType>
107  struct Serial {
108  template<typename outputValueViewType,
109  typename inputPointViewType,
110  typename workViewType,
111  typename vinvViewType>
112  KOKKOS_INLINE_FUNCTION
113  static void
114  getValues( outputValueViewType outputValues,
115  const inputPointViewType inputPoints,
116  workViewType work,
117  const vinvViewType vinv );
118 
119 
120  KOKKOS_INLINE_FUNCTION
121  static ordinal_type
122  getWorkSizePerPoint(ordinal_type order) {
123  auto cardinality = CardinalityHCurlTet(order);
124  switch (opType) {
125  case OPERATOR_DIV:
126  case OPERATOR_D1:
127  return 7*cardinality;
128  default:
129  return getDkCardinality<opType,3>()*cardinality;
130  }
131  }
132  };
133 
134  template<typename DeviceType, ordinal_type numPtsPerEval,
135  typename outputValueValueType, class ...outputValueProperties,
136  typename inputPointValueType, class ...inputPointProperties,
137  typename vinvValueType, class ...vinvProperties>
138  static void
139  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
140  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
141  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
142  const EOperator operatorType);
143 
147  template<typename outputValueViewType,
148  typename inputPointViewType,
149  typename vinvViewType,
150  typename workViewType,
151  EOperator opType,
152  ordinal_type numPtsEval>
153  struct Functor {
154  outputValueViewType _outputValues;
155  const inputPointViewType _inputPoints;
156  const vinvViewType _coeffs;
157  workViewType _work;
158 
159  KOKKOS_INLINE_FUNCTION
160  Functor( outputValueViewType outputValues_,
161  inputPointViewType inputPoints_,
162  vinvViewType coeffs_,
163  workViewType work_)
164  : _outputValues(outputValues_), _inputPoints(inputPoints_),
165  _coeffs(coeffs_), _work(work_) {}
166 
167  KOKKOS_INLINE_FUNCTION
168  void operator()(const size_type iter) const {
169  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
170  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
171 
172  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
173  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
174 
175  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
176 
177  auto vcprop = Kokkos::common_view_alloc_prop(_work);
178  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
179 
180  switch (opType) {
181  case OPERATOR_VALUE : {
182  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
183  Serial<opType>::getValues( output, input, work, _coeffs );
184  break;
185  }
186  case OPERATOR_CURL: {
187  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
188  Serial<opType>::getValues( output, input, work, _coeffs );
189  break;
190  }
191  default: {
192  INTREPID2_TEST_FOR_ABORT( true,
193  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::Functor) operator is not supported");
194 
195  }
196  }
197  }
198  };
199 };
200 }
201 
202 template<typename DeviceType = void,
203  typename outputValueType = double,
204  typename pointValueType = double>
206  : public Basis<DeviceType,outputValueType,pointValueType> {
207  public:
208  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost OrdinalTypeArray1DHost;
209  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost OrdinalTypeArray2DHost;
210  typedef typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost OrdinalTypeArray3DHost;
211 
214  Basis_HCURL_TET_In_FEM(const ordinal_type order,
215  const EPointType pointType = POINTTYPE_EQUISPACED);
216 
217 
221 
223 
225 
226  virtual
227  void
228  getValues( OutputViewType outputValues,
229  const PointViewType inputPoints,
230  const EOperator operatorType = OPERATOR_VALUE) const override {
231 #ifdef HAVE_INTREPID2_DEBUG
232  Intrepid2::getValues_HCURL_Args(outputValues,
233  inputPoints,
234  operatorType,
235  this->getBaseCellTopology(),
236  this->getCardinality() );
237 #endif
238  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
239  Impl::Basis_HCURL_TET_In_FEM::
240  getValues<DeviceType,numPtsPerEval>( outputValues,
241  inputPoints,
242  this->coeffs_,
243  operatorType);
244  }
245 
246  virtual
247  void
248  getDofCoords( ScalarViewType dofCoords ) const override {
249 #ifdef HAVE_INTREPID2_DEBUG
250  // Verify rank of output array.
251  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
252  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
253  // Verify 0th dimension of output array.
254  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
255  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
256  // Verify 1st dimension of output array.
257  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
258  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
259 #endif
260  Kokkos::deep_copy(dofCoords, this->dofCoords_);
261  }
262 
263  virtual
264  void
265  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
266 #ifdef HAVE_INTREPID2_DEBUG
267  // Verify rank of output array.
268  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
269  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
270  // Verify 0th dimension of output array.
271  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
272  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
273  // Verify 1st dimension of output array.
274  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
275  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
276 #endif
277  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
278  }
279 
280  void
281  getExpansionCoeffs( ScalarViewType coeffs ) const {
282  // has to be same rank and dimensions
283  Kokkos::deep_copy(coeffs, this->coeffs_);
284  }
285 
286  virtual
287  const char*
288  getName() const override {
289  return "Intrepid2_HCURL_TET_In_FEM";
290  }
291 
292  virtual
293  bool
294  requireOrientation() const override {
295  return true;
296  }
297 
308  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
309  if(subCellDim == 1) {
310  return Teuchos::rcp(new
312  (this->basisDegree_-1, pointType_));
313  } else if(subCellDim == 2) {
314  return Teuchos::rcp(new
316  (this->basisDegree_, pointType_));
317  }
318  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
319  }
320 
322  getHostBasis() const override{
324  }
325 
326  private:
327 
330  Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
331 
334 
335 };
336 
337 }// namespace Intrepid2
338 
340 
341 #endif
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
small utility functions
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
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...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type getCardinality() const
Returns cardinality of the basis.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
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.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
Definition file for FEM basis functions of degree n for H(curl) functions on TET. ...
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...
virtual bool requireOrientation() const override
True if orientation is required.
See Intrepid2::Basis_HCURL_TET_In_FEM.
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.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
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.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
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.