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 
222 
224 
225  virtual
226  void
227  getValues( OutputViewType outputValues,
228  const PointViewType inputPoints,
229  const EOperator operatorType = OPERATOR_VALUE) const {
230 #ifdef HAVE_INTREPID2_DEBUG
231  Intrepid2::getValues_HCURL_Args(outputValues,
232  inputPoints,
233  operatorType,
234  this->getBaseCellTopology(),
235  this->getCardinality() );
236 #endif
237  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
238  Impl::Basis_HCURL_TET_In_FEM::
239  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
240  inputPoints,
241  this->coeffs_,
242  operatorType);
243  }
244 
245  virtual
246  void
247  getDofCoords( ScalarViewType dofCoords ) const {
248 #ifdef HAVE_INTREPID2_DEBUG
249  // Verify rank of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_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_TET_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_TET_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 {
265 #ifdef HAVE_INTREPID2_DEBUG
266  // Verify rank of output array.
267  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
268  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_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_TET_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_TET_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 {
288  return "Intrepid2_HCURL_TET_In_FEM";
289  }
290 
291  virtual
292  bool
294  return true;
295  }
296 
297  private:
298 
301  Kokkos::DynRankView<scalarType,ExecSpaceType> coeffs_;
302 
303 };
304 
305 }// namespace Intrepid2
306 
308 
309 #endif
Header file for the abstract base class Intrepid2::Basis.
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....
Definition file for FEM basis functions of degree n for H(curl) functions on TET.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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
virtual bool requireOrientation() const
True if orientation is required.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, void > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, void > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
ordinal_type getCardinality() const
Returns cardinality of the basis.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
See Intrepid2::Basis_HCURL_TET_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions