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"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 namespace Intrepid2 {
59 
92 #define CardinalityHCurlTet(order) (order*(order+2)*(order+3)/2)
93 
94 namespace Impl {
95 
100 public:
101  typedef struct Tetrahedron<4> cell_topology_type;
105  template<EOperator opType>
106  struct Serial {
107  template<typename outputValueViewType,
108  typename inputPointViewType,
109  typename workViewType,
110  typename vinvViewType>
111  KOKKOS_INLINE_FUNCTION
112  static void
113  getValues( outputValueViewType outputValues,
114  const inputPointViewType inputPoints,
115  workViewType work,
116  const vinvViewType vinv );
117 
118 
119  KOKKOS_INLINE_FUNCTION
120  static ordinal_type
121  getWorkSizePerPoint(ordinal_type order) {
122  auto cardinality = CardinalityHCurlTet(order);
123  switch (opType) {
124  case OPERATOR_DIV:
125  case OPERATOR_D1:
126  return 7*cardinality;
127  default:
128  return getDkCardinality<opType,3>()*cardinality;
129  }
130  }
131  };
132 
133  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
134  typename outputValueValueType, class ...outputValueProperties,
135  typename inputPointValueType, class ...inputPointProperties,
136  typename vinvValueType, class ...vinvProperties>
137  static void
138  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
139  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
140  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
141  const EOperator operatorType);
142 
146  template<typename outputValueViewType,
147  typename inputPointViewType,
148  typename vinvViewType,
149  typename workViewType,
150  EOperator opType,
151  ordinal_type numPtsEval>
152  struct Functor {
153  outputValueViewType _outputValues;
154  const inputPointViewType _inputPoints;
155  const vinvViewType _coeffs;
156  workViewType _work;
157 
158  KOKKOS_INLINE_FUNCTION
159  Functor( outputValueViewType outputValues_,
160  inputPointViewType inputPoints_,
161  vinvViewType coeffs_,
162  workViewType work_)
163  : _outputValues(outputValues_), _inputPoints(inputPoints_),
164  _coeffs(coeffs_), _work(work_) {}
165 
166  KOKKOS_INLINE_FUNCTION
167  void operator()(const size_type iter) const {
168  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
169  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
170 
171  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
172  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
173 
174  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
175 
176  auto vcprop = Kokkos::common_view_alloc_prop(_work);
177  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
178 
179  switch (opType) {
180  case OPERATOR_VALUE : {
181  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
182  Serial<opType>::getValues( output, input, work, _coeffs );
183  break;
184  }
185  case OPERATOR_CURL: {
186  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
187  Serial<opType>::getValues( output, input, work, _coeffs );
188  break;
189  }
190  default: {
191  INTREPID2_TEST_FOR_ABORT( true,
192  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::Functor) operator is not supported");
193 
194  }
195  }
196  }
197  };
198 };
199 }
200 
201 template<typename ExecSpaceType = void,
202  typename outputValueType = double,
203  typename pointValueType = double>
205  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
206  public:
210 
213  Basis_HCURL_TET_In_FEM(const ordinal_type order,
214  const EPointType pointType = POINTTYPE_EQUISPACED);
215 
216 
220 
221  using outputViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use OutputViewType instead","OutputViewType") = OutputViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use OutputViewType instead");
222  using pointViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use PointViewType instead","PointViewType") = PointViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use PointViewType instead");
223  using scalarViewType INTREPID2_DEPRECATED_TYPENAME_REPLACEMENT("use ScalarViewType instead","ScalarViewType") = ScalarViewType INTREPID2_DEPRECATED_TYPENAME_TRAILING_ATTRIBUTE("use ScalarViewType instead");
224 
226 
228 
229  virtual
230  void
231  getValues( OutputViewType outputValues,
232  const PointViewType inputPoints,
233  const EOperator operatorType = OPERATOR_VALUE) const {
234 #ifdef HAVE_INTREPID2_DEBUG
235  Intrepid2::getValues_HCURL_Args(outputValues,
236  inputPoints,
237  operatorType,
238  this->getBaseCellTopology(),
239  this->getCardinality() );
240 #endif
241  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
242  Impl::Basis_HCURL_TET_In_FEM::
243  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
244  inputPoints,
245  this->coeffs_,
246  operatorType);
247  }
248 
249  virtual
250  void
251  getDofCoords( ScalarViewType dofCoords ) const {
252 #ifdef HAVE_INTREPID2_DEBUG
253  // Verify rank of output array.
254  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
255  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
256  // Verify 0th dimension of output array.
257  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
258  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
259  // Verify 1st dimension of output array.
260  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
261  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
262 #endif
263  Kokkos::deep_copy(dofCoords, this->dofCoords_);
264  }
265 
266  virtual
267  void
268  getDofCoeffs( ScalarViewType dofCoeffs ) const {
269 #ifdef HAVE_INTREPID2_DEBUG
270  // Verify rank of output array.
271  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
272  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
273  // Verify 0th dimension of output array.
274  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
275  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
276  // Verify 1st dimension of output array.
277  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
278  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
279 #endif
280  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
281  }
282 
283  void
284  getExpansionCoeffs( ScalarViewType coeffs ) const {
285  // has to be same rank and dimensions
286  Kokkos::deep_copy(coeffs, this->coeffs_);
287  }
288 
289  virtual
290  const char*
291  getName() const {
292  return "Intrepid2_HCURL_TET_In_FEM";
293  }
294 
295  virtual
296  bool
298  return true;
299  }
300 
301  private:
302 
305  Kokkos::DynRankView<scalarType,ExecSpaceType> coeffs_;
306 
307 };
308 
309 }// namespace Intrepid2
310 
312 
313 #endif
small utility functions
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
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.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Definition file for FEM basis functions of degree n for H(curl) functions on TET. ...
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...
See Intrepid2::Basis_HCURL_TET_In_FEM.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
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.