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 ExecSpaceType, 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 ExecSpaceType = void,
175  typename outputValueType = double,
176  typename pointValueType = double>
178  : public Basis<ExecSpaceType,outputValueType,pointValueType> {
179  public:
183 
186  Basis_HVOL_LINE_Cn_FEM(const ordinal_type order,
187  const EPointType pointType = POINTTYPE_EQUISPACED);
188 
192 
194 
195  virtual
196  void
197  getValues( outputViewType outputValues,
198  const pointViewType inputPoints,
199  const EOperator operatorType = OPERATOR_VALUE ) const {
200 #ifdef HAVE_INTREPID2_DEBUG
201  Intrepid2::getValues_HVOL_Args(outputValues,
202  inputPoints,
203  operatorType,
204  this->getBaseCellTopology(),
205  this->getCardinality() );
206 #endif
207  constexpr ordinal_type numPtsPerEval = 1;
208  Impl::Basis_HVOL_LINE_Cn_FEM::
209  getValues<ExecSpaceType,numPtsPerEval>( outputValues,
210  inputPoints,
211  this->vinv_,
212  operatorType );
213  }
214 
215  virtual
216  void
217  getDofCoords( scalarViewType dofCoords ) const {
218 #ifdef HAVE_INTREPID2_DEBUG
219  // Verify rank of output array.
220  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
222  // Verify 0th dimension of output array.
223  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
224  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225  // Verify 1st dimension of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228 #endif
229  Kokkos::deep_copy(dofCoords, this->dofCoords_);
230  }
231 
232  virtual
233  void
234  getDofCoeffs( scalarViewType dofCoeffs ) const {
235 #ifdef HAVE_INTREPID2_DEBUG
236  // Verify rank of output array.
237  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
239  // Verify 0th dimension of output array.
240  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241  ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242 #endif
243  Kokkos::deep_copy(dofCoeffs, 1.0);
244  }
245 
246  void
247  getVandermondeInverse( scalarViewType vinv ) const {
248  // has to be same rank and dimensions
249  Kokkos::deep_copy(vinv, this->vinv_);
250  }
251 
252  virtual
253  const char*
254  getName() const {
255  return "Intrepid2_HVOL_LINE_Cn_FEM";
256  }
257 
258  private:
259 
262  Kokkos::DynRankView<typename scalarViewType::value_type,ExecSpaceType> vinv_;
263  };
264 
265 }// namespace Intrepid2
266 
268 
269 #endif
small utility functions
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
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.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
Definition file for FEM basis functions of degree n for H(vol) functions on LINE. ...
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< typename scalarViewType::value_type, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.