Intrepid2
Intrepid2_HVOL_TRI_Cn_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 
48 #ifndef __INTREPID2_HVOL_TRI_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_TRI_CN_FEM_HPP__
50 
51 #include "Intrepid2_Basis.hpp"
52 
53 #include "Intrepid2_PointTools.hpp"
54 #include "Teuchos_LAPACK.hpp"
55 
56 namespace Intrepid2 {
57 
72  namespace Impl {
73 
78  public:
79  typedef struct Triangle<3> cell_topology_type;
83  template<EOperator opType>
84  struct Serial {
85  template<typename outputValueViewType,
86  typename inputPointViewType,
87  typename workViewType,
88  typename vinvViewType>
89  KOKKOS_INLINE_FUNCTION
90  static void
91  getValues( outputValueViewType outputValues,
92  const inputPointViewType inputPoints,
93  workViewType work,
94  const vinvViewType vinv );
95 
96 
97  KOKKOS_INLINE_FUNCTION
98  static ordinal_type
99  getWorkSizePerPoint(ordinal_type order) {
100  auto cardinality = getPnCardinality<2>(order);
101  switch (opType) {
102  case OPERATOR_GRAD:
103  case OPERATOR_CURL:
104  case OPERATOR_D1:
105  return 5*cardinality;
106  default:
107  return getDkCardinality<opType,2>()*cardinality;
108  }
109  }
110  };
111 
112  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
113  typename outputValueValueType, class ...outputValueProperties,
114  typename inputPointValueType, class ...inputPointProperties,
115  typename vinvValueType, class ...vinvProperties>
116  static void
117  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
118  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
119  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
120  const EOperator operatorType);
121 
125  template<typename outputValueViewType,
126  typename inputPointViewType,
127  typename vinvViewType,
128  typename workViewType,
129  EOperator opType,
130  ordinal_type numPtsEval>
131  struct Functor {
132  outputValueViewType _outputValues;
133  const inputPointViewType _inputPoints;
134  const vinvViewType _vinv;
135  workViewType _work;
136 
137  KOKKOS_INLINE_FUNCTION
138  Functor( outputValueViewType outputValues_,
139  inputPointViewType inputPoints_,
140  vinvViewType vinv_,
141  workViewType work_)
142  : _outputValues(outputValues_), _inputPoints(inputPoints_),
143  _vinv(vinv_), _work(work_) {}
144 
145  KOKKOS_INLINE_FUNCTION
146  void operator()(const size_type iter) const {
147  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
148  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
149 
150  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
151  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
152 
153  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
154 
155  auto vcprop = Kokkos::common_view_alloc_prop(_work);
156  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
157 
158  switch (opType) {
159  case OPERATOR_VALUE : {
160  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
161  Serial<opType>::getValues( output, input, work, _vinv );
162  break;
163  }
164  case OPERATOR_D1:
165  case OPERATOR_D2: {
166  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
167  Serial<opType>::getValues( output, input, work, _vinv );
168  break;
169  }
170  default: {
171  INTREPID2_TEST_FOR_ABORT( true,
172  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::Functor) operator is not supported");
173 
174  }
175  }
176  }
177  };
178  };
179  }
180 
181  template<typename ExecSpaceType = void,
182  typename outputValueType = double,
183  typename pointValueType = double>
185  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
186  public:
190 
193  Basis_HVOL_TRI_Cn_FEM(const ordinal_type order,
194  const EPointType pointType = POINTTYPE_EQUISPACED);
195 
196 
200 
202 
204 
205  virtual
206  void
207  getValues( OutputViewType outputValues,
208  const PointViewType inputPoints,
209  const EOperator operatorType = OPERATOR_VALUE) const {
210 #ifdef HAVE_INTREPID2_DEBUG
211  Intrepid2::getValues_HVOL_Args(outputValues,
212  inputPoints,
213  operatorType,
214  this->getBaseCellTopology(),
215  this->getCardinality() );
216 #endif
217  constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
218  Impl::Basis_HVOL_TRI_Cn_FEM::
219  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
220  inputPoints,
221  this->vinv_,
222  operatorType);
223  }
224 
225  virtual
226  void
227  getDofCoords( ScalarViewType dofCoords ) const {
228 #ifdef HAVE_INTREPID2_DEBUG
229  // Verify rank of output array.
230  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
231  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
232  // Verify 0th dimension of output array.
233  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
234  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
235  // Verify 1st dimension of output array.
236  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
237  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
238 #endif
239  Kokkos::deep_copy(dofCoords, this->dofCoords_);
240  }
241 
242  virtual
243  void
244  getDofCoeffs( ScalarViewType dofCoeffs ) const {
245 #ifdef HAVE_INTREPID2_DEBUG
246  // Verify rank of output array.
247  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
248  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
249  // Verify 0th dimension of output array.
250  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
251  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
252 #endif
253  Kokkos::deep_copy(dofCoeffs, 1.0);
254  }
255 
256  void
257  getVandermondeInverse( ScalarViewType vinv ) const {
258  // has to be same rank and dimensions
259  Kokkos::deep_copy(vinv, this->vinv_);
260  }
261 
262  virtual
263  const char*
264  getName() const {
265  return "Intrepid2_HVOL_TRI_Cn_FEM";
266  }
267 
268  virtual
269  bool
271  return false;
272  }
273 
274  private:
275 
278  Kokkos::DynRankView<scalarType,ExecSpaceType> vinv_;
279 
280  };
281 
282 }// namespace Intrepid2
283 
285 
286 #endif
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
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.
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.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
virtual bool requireOrientation() const
True if orientation is required.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Definition file for FEM basis functions of degree n for H(vol) functions on TRI.
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.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual const char * getName() const
Returns basis name.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.