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"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 namespace Intrepid2 {
59 
80 #define CardinalityHCurlTri(order) (order*(order+2))
81 
82 namespace Impl {
83 
88 public:
89  typedef struct Triangle<3> cell_topology_type;
90 
94  template<EOperator opType>
95  struct Serial {
96  template<typename outputValueViewType,
97  typename inputPointViewType,
98  typename workViewType,
99  typename vinvViewType>
100  KOKKOS_INLINE_FUNCTION
101  static void
102  getValues( outputValueViewType outputValues,
103  const inputPointViewType inputPoints,
104  workViewType work,
105  const vinvViewType vinv );
106 
107  KOKKOS_INLINE_FUNCTION
108  static ordinal_type
109  getWorkSizePerPoint(ordinal_type order) {
110  auto cardinality = CardinalityHCurlTri(order);
111  switch (opType) {
112  case OPERATOR_GRAD:
113  case OPERATOR_CURL:
114  case OPERATOR_D1:
115  return 5*cardinality;
116  default:
117  return getDkCardinality<opType,2>()*cardinality;
118  }
119  }
120  };
121 
122  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
123  typename outputValueValueType, class ...outputValueProperties,
124  typename inputPointValueType, class ...inputPointProperties,
125  typename vinvValueType, class ...vinvProperties>
126  static void
127  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
130  const EOperator operatorType);
131 
135  template<typename outputValueViewType,
136  typename inputPointViewType,
137  typename vinvViewType,
138  typename workViewType,
139  EOperator opType,
140  ordinal_type numPtsEval>
141  struct Functor {
142  outputValueViewType _outputValues;
143  const inputPointViewType _inputPoints;
144  const vinvViewType _coeffs;
145  workViewType _work;
146 
147  KOKKOS_INLINE_FUNCTION
148  Functor( outputValueViewType outputValues_,
149  inputPointViewType inputPoints_,
150  vinvViewType coeffs_,
151  workViewType work_)
152  : _outputValues(outputValues_), _inputPoints(inputPoints_),
153  _coeffs(coeffs_), _work(work_) {}
154 
155  KOKKOS_INLINE_FUNCTION
156  void operator()(const size_type iter) const {
157  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
158  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
159 
160  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
161  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
162 
163 
164  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
165 
166  auto vcprop = Kokkos::common_view_alloc_prop(_work);
167  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
168 
169  switch (opType) {
170  case OPERATOR_VALUE : {
171  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
172  Serial<opType>::getValues( output, input, work, _coeffs );
173  break;
174  }
175  case OPERATOR_CURL: {
176  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
177  Serial<opType>::getValues( output, input, work, _coeffs );
178  break;
179  }
180  default: {
181  INTREPID2_TEST_FOR_ABORT( true,
182  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
183 
184  }
185  }
186  }
187  };
188 };
189 }
190 
191 template<typename ExecSpaceType = void,
192  typename outputValueType = double,
193  typename pointValueType = double>
195  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
196  public:
200 
203  Basis_HCURL_TRI_In_FEM(const ordinal_type order,
204  const EPointType pointType = POINTTYPE_EQUISPACED);
205 
206 
210 
212 
214 
215  virtual
216  void
217  getValues( OutputViewType outputValues,
218  const PointViewType inputPoints,
219  const EOperator operatorType = OPERATOR_VALUE) const {
220 #ifdef HAVE_INTREPID2_DEBUG
221  Intrepid2::getValues_HCURL_Args(outputValues,
222  inputPoints,
223  operatorType,
224  this->getBaseCellTopology(),
225  this->getCardinality() );
226 #endif
227  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
228  Impl::Basis_HCURL_TRI_In_FEM::
229  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
230  inputPoints,
231  this->coeffs_,
232  operatorType);
233  }
234 
235  virtual
236  void
237  getDofCoords( ScalarViewType dofCoords ) const {
238 #ifdef HAVE_INTREPID2_DEBUG
239  // Verify rank of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
242  // Verify 0th dimension of output array.
243  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
244  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
245  // Verify 1st dimension of output array.
246  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
247  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
248 #endif
249  Kokkos::deep_copy(dofCoords, this->dofCoords_);
250  }
251 
252  virtual
253  void
254  getDofCoeffs( ScalarViewType dofCoeffs ) const {
255 #ifdef HAVE_INTREPID2_DEBUG
256  // Verify rank of output array.
257  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
258  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
259  // Verify 0th dimension of output array.
260  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
261  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
262  // Verify 1st dimension of output array.
263  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
264  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
265 #endif
266  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
267  }
268 
269  void
270  getExpansionCoeffs( ScalarViewType coeffs ) const {
271  // has to be same rank and dimensions
272  Kokkos::deep_copy(coeffs, this->coeffs_);
273  }
274 
275  virtual
276  const char*
277  getName() const {
278  return "Intrepid2_HCURL_TRI_In_FEM";
279  }
280 
281  virtual
282  bool
284  return true;
285  }
286 
287  private:
288 
291  Kokkos::DynRankView<scalarType,ExecSpaceType> coeffs_;
292 
293 };
294 
295 }// namespace Intrepid2
296 
298 
299 #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 getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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 getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual const char * getName() const
Returns basis name.
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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< 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...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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.