Intrepid2
Intrepid2_HGRAD_LINE_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 
49 #ifndef __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
51 
52 #include "Intrepid2_Basis.hpp"
54 
55 #include "Intrepid2_PointTools.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 namespace Intrepid2 {
59 
76  namespace Impl {
77 
82  public:
83  typedef struct Line<2> cell_topology_type;
87  template<EOperator opType>
88  struct Serial {
89  template<typename outputValueViewType,
90  typename inputPointViewType,
91  typename workViewType,
92  typename vinvViewType>
93  KOKKOS_INLINE_FUNCTION
94  static void
95  getValues( outputValueViewType outputValues,
96  const inputPointViewType inputPoints,
97  workViewType work,
98  const vinvViewType vinv,
99  const ordinal_type operatorDn = 0 );
100  };
101 
102  template<typename ExecSpaceType, ordinal_type numPtsPerEval,
103  typename outputValueValueType, class ...outputValueProperties,
104  typename inputPointValueType, class ...inputPointProperties,
105  typename vinvValueType, class ...vinvProperties>
106  static void
107  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
108  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
109  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
110  const EOperator operatorType );
111 
115  template<typename outputValueViewType,
116  typename inputPointViewType,
117  typename vinvViewType,
118  typename workViewType,
119  EOperator opType,
120  ordinal_type numPtsEval>
121  struct Functor {
122  outputValueViewType _outputValues;
123  const inputPointViewType _inputPoints;
124  const vinvViewType _vinv;
125  workViewType _work;
126  const ordinal_type _opDn;
127 
128  KOKKOS_INLINE_FUNCTION
129  Functor( outputValueViewType outputValues_,
130  inputPointViewType inputPoints_,
131  vinvViewType vinv_,
132  workViewType work_,
133  const ordinal_type opDn_ = 0 )
134  : _outputValues(outputValues_), _inputPoints(inputPoints_),
135  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
136 
137  KOKKOS_INLINE_FUNCTION
138  void operator()(const size_type iter) const {
139  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
140  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
141 
142  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
143  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
144 
145  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
146 
147  auto vcprop = Kokkos::common_view_alloc_prop(_work);
148  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
149 
150  switch (opType) {
151  case OPERATOR_VALUE : {
152  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
153  Serial<opType>::getValues( output, input, work, _vinv );
154  break;
155  }
156  case OPERATOR_Dn : {
157  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
158  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
159  break;
160  }
161  default: {
162  INTREPID2_TEST_FOR_ABORT( true,
163  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Functor) operator is not supported");
164 
165  }
166  }
167  }
168  };
169  };
170  }
171 
172  template<typename ExecSpaceType = void,
173  typename outputValueType = double,
174  typename pointValueType = double>
176  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
177  public:
181 
185 
186  private:
187 
190  Kokkos::DynRankView<typename ScalarViewType::value_type,ExecSpaceType> vinv_;
191 
192  public:
195  Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order,
196  const EPointType pointType = POINTTYPE_EQUISPACED);
197 
199 
200  virtual
201  void
202  getValues( OutputViewType outputValues,
203  const PointViewType inputPoints,
204  const EOperator operatorType = OPERATOR_VALUE ) const {
205 #ifdef HAVE_INTREPID2_DEBUG
206  Intrepid2::getValues_HGRAD_Args(outputValues,
207  inputPoints,
208  operatorType,
209  this->getBaseCellTopology(),
210  this->getCardinality() );
211 #endif
212  constexpr ordinal_type numPtsPerEval = 1;
213  Impl::Basis_HGRAD_LINE_Cn_FEM::
214  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
215  inputPoints,
216  this->vinv_,
217  operatorType );
218  }
219 
220  virtual
221  void
222  getDofCoords( ScalarViewType dofCoords ) const {
223 #ifdef HAVE_INTREPID2_DEBUG
224  // Verify rank of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
227  // Verify 0th dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
230  // Verify 1st dimension of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
233 #endif
234  Kokkos::deep_copy(dofCoords, this->dofCoords_);
235  }
236 
237  virtual
238  void
239  getDofCoeffs( ScalarViewType dofCoeffs ) const {
240 #ifdef HAVE_INTREPID2_DEBUG
241  // Verify rank of output array.
242  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
243  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
244  // Verify 0th dimension of output array.
245  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
246  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
247 #endif
248  Kokkos::deep_copy(dofCoeffs, 1.0);
249  }
250 
251  virtual
252  const char*
253  getName() const {
254  return "Intrepid2_HGRAD_LINE_Cn_FEM";
255  }
256 
257  void
258  getVandermondeInverse( ScalarViewType vinv ) const {
259  // has to be same rank and dimensions
260  Kokkos::deep_copy(vinv, this->vinv_);
261  }
262 
263  Kokkos::DynRankView<typename ScalarViewType::const_value_type,ExecSpaceType>
264  getVandermondeInverse() const {
265  return vinv_;
266  }
267 
268  ordinal_type
269  getWorkSizePerPoint(const EOperator operatorType) const {
270  return getPnCardinality<1>(this->basisDegree_);
271  }
272 
273  };
274 
275 }// namespace Intrepid2
276 
278 
279 #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.
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.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
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< 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.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Definition file for FEM basis functions of degree n for H(grad) functions on LINE.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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.