Intrepid2
Intrepid2_HVOL_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 
48 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_HPP__
49 #define __INTREPID2_HVOL_LINE_CN_FEM_HPP__
50 
51 #include "Intrepid2_Basis.hpp"
53 
54 #include "Intrepid2_PointTools.hpp"
55 #include "Teuchos_LAPACK.hpp"
56 
57 namespace Intrepid2 {
58 
74  namespace Impl {
75 
80  public:
81  typedef struct Line<2> cell_topology_type;
85  template<EOperator opType>
86  struct Serial {
87  template<typename outputValueViewType,
88  typename inputPointViewType,
89  typename workViewType,
90  typename vinvViewType>
91  KOKKOS_INLINE_FUNCTION
92  static void
93  getValues( outputValueViewType outputValues,
94  const inputPointViewType inputPoints,
95  workViewType work,
96  const vinvViewType vinv,
97  const ordinal_type operatorDn = 0 );
98 
99  KOKKOS_INLINE_FUNCTION
100  static ordinal_type
101  getWorkSizePerPoint(ordinal_type order) {return getPnCardinality<1>(order);}
102  };
103 
104  template<typename DeviceType, ordinal_type numPtsPerEval,
105  typename outputValueValueType, class ...outputValueProperties,
106  typename inputPointValueType, class ...inputPointProperties,
107  typename vinvValueType, class ...vinvProperties>
108  static void
109  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
110  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
111  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
112  const EOperator operatorType );
113 
117  template<typename outputValueViewType,
118  typename inputPointViewType,
119  typename vinvViewType,
120  typename workViewType,
121  EOperator opType,
122  ordinal_type numPtsEval>
123  struct Functor {
124  outputValueViewType _outputValues;
125  const inputPointViewType _inputPoints;
126  const vinvViewType _vinv;
127  workViewType _work;
128  const ordinal_type _opDn;
129 
130  KOKKOS_INLINE_FUNCTION
131  Functor( outputValueViewType outputValues_,
132  inputPointViewType inputPoints_,
133  vinvViewType vinv_,
134  workViewType work_,
135  const ordinal_type opDn_ = 0 )
136  : _outputValues(outputValues_), _inputPoints(inputPoints_),
137  _vinv(vinv_), _work(work_), _opDn(opDn_) {}
138 
139  KOKKOS_INLINE_FUNCTION
140  void operator()(const size_type iter) const {
141  const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
142  const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
143 
144  const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
145  const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
146 
147  typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
148 
149  auto vcprop = Kokkos::common_view_alloc_prop(_work);
150  workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
151 
152  switch (opType) {
153  case OPERATOR_VALUE : {
154  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
155  Serial<opType>::getValues( output, input, work, _vinv );
156  break;
157  }
158  case OPERATOR_Dn : {
159  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
160  Serial<opType>::getValues( output, input, work, _vinv, _opDn );
161  break;
162  }
163  default: {
164  INTREPID2_TEST_FOR_ABORT( true,
165  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Functor) operator is not supported");
166 
167  }
168  }
169  }
170  };
171  };
172  }
173 
174  template<typename DeviceType = void,
175  typename outputValueType = double,
176  typename pointValueType = double>
178  : public Basis<DeviceType,outputValueType,pointValueType> {
179  public:
182 
183  using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
184  using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
185  using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
186 
187  using OutputViewType = typename BasisBase::OutputViewType;
188  using PointViewType = typename BasisBase::PointViewType ;
189  using ScalarViewType = typename BasisBase::ScalarViewType;
190 
193  Basis_HVOL_LINE_Cn_FEM(const ordinal_type order,
194  const EPointType pointType = POINTTYPE_EQUISPACED);
195 
196  using BasisBase::getValues;
197 
198  virtual
199  void
200  getValues( OutputViewType outputValues,
201  const PointViewType inputPoints,
202  const EOperator operatorType = OPERATOR_VALUE ) const override {
203 #ifdef HAVE_INTREPID2_DEBUG
204  Intrepid2::getValues_HVOL_Args(outputValues,
205  inputPoints,
206  operatorType,
207  this->getBaseCellTopology(),
208  this->getCardinality() );
209 #endif
210  constexpr ordinal_type numPtsPerEval = 1;
211  Impl::Basis_HVOL_LINE_Cn_FEM::
212  getValues<DeviceType,numPtsPerEval>( outputValues,
213  inputPoints,
214  this->vinv_,
215  operatorType );
216  }
217 
218  virtual
219  void
220  getDofCoords( ScalarViewType dofCoords ) const override {
221 #ifdef HAVE_INTREPID2_DEBUG
222  // Verify rank of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
225  // Verify 0th dimension of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
228  // Verify 1st dimension of output array.
229  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
230  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
231 #endif
232  Kokkos::deep_copy(dofCoords, this->dofCoords_);
233  }
234 
235  virtual
236  void
237  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
238 #ifdef HAVE_INTREPID2_DEBUG
239  // Verify rank of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
242  // Verify 0th dimension of output array.
243  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
244  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
245 #endif
246  Kokkos::deep_copy(dofCoeffs, 1.0);
247  }
248 
249  void
250  getVandermondeInverse( ScalarViewType vinv ) const {
251  // has to be same rank and dimensions
252  Kokkos::deep_copy(vinv, this->vinv_);
253  }
254 
255  virtual
256  const char*
257  getName() const override {
258  return "Intrepid2_HVOL_LINE_Cn_FEM";
259  }
260 
262  getHostBasis() const override{
264  }
265 
266  private:
267 
270  Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
271  EPointType pointType_;
272  };
273 
274 }// namespace Intrepid2
275 
277 
278 #endif
small utility functions
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type), allowing host access to input and output views, and ensuring host execution of basis evaluation.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
virtual const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
void getValues_HVOL_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 HVOL-conforming FEM basis...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
EPointType
Enumeration of types of point distributions in Intrepid.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
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.