Intrepid2
Intrepid2_LegendreBasis_HVOL_LINE.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),
38 // Mauro Perego (mperego@sandia.gov), or
39 // Nate Roberts (nvrober@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #ifndef Intrepid2_LegendreBasis_HVOL_LINE_h
50 #define Intrepid2_LegendreBasis_HVOL_LINE_h
51 
52 #include <Kokkos_View.hpp>
53 #include <Kokkos_DynRankView.hpp>
54 
55 #include <Intrepid2_config.h>
56 
57 // Sacado header that defines some fad sizing…
58 #ifdef HAVE_INTREPID2_SACADO
59 #include <KokkosExp_View_Fad.hpp>
60 #endif
61 
64 
65 namespace Intrepid2
66 {
72  template<class ExecutionSpace, class OutputScalar, class PointScalar,
73  class OutputFieldType, class InputPointsType>
75  {
76  using ScratchSpace = typename ExecutionSpace::scratch_memory_space;
77  using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
78  using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
79 
80  using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
81  using TeamMember = typename TeamPolicy::member_type;
82 
83  OutputFieldType output_; // F,P
84  InputPointsType inputPoints_; // P,D
85 
86  int polyOrder_;
87  int numFields_, numPoints_;
88 
89  EOperator op_;
90 
91  size_t fad_size_output_;
92 
93  Hierarchical_HVOL_LINE_Functor(OutputFieldType output, InputPointsType inputPoints,
94  int polyOrder, EOperator op)
95  : output_(output), inputPoints_(inputPoints), polyOrder_(polyOrder), op_(op),
96  fad_size_output_(getScalarDimensionForView(output))
97  {
98  numFields_ = output.extent_int(0);
99  numPoints_ = output.extent_int(1);
100  INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, "point counts need to match!");
101  INTREPID2_TEST_FOR_EXCEPTION(numFields_ != polyOrder_+1, std::invalid_argument, "output field size does not match basis cardinality");
102  }
103 
104  KOKKOS_INLINE_FUNCTION
105  void operator()( const TeamMember & teamMember ) const
106  {
107  auto pointOrdinal = teamMember.league_rank();
108  OutputScratchView field_values_at_point;
109  if (fad_size_output_ > 0) {
110  field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_, fad_size_output_);
111  }
112  else {
113  field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_);
114  }
115 
116  const PointScalar x = inputPoints_(pointOrdinal,0);
117 
118  switch (op_)
119  {
120  case OPERATOR_VALUE:
121  Polynomials::legendreValues(field_values_at_point, polyOrder_, x);
122 
123  // note that because legendreValues determines field values recursively, there is not much
124  // opportunity at that level for further parallelism
125  break;
126  case OPERATOR_GRAD:
127  case OPERATOR_D1:
128  case OPERATOR_D2:
129  case OPERATOR_D3:
130  case OPERATOR_D4:
131  case OPERATOR_D5:
132  case OPERATOR_D6:
133  case OPERATOR_D7:
134  case OPERATOR_D8:
135  case OPERATOR_D9:
136  case OPERATOR_D10:
137  {
138  auto derivativeOrder = getOperatorOrder(op_);
139  Polynomials::legendreDerivativeValues(field_values_at_point, polyOrder_, x, derivativeOrder);
140  break;
141  }
142  default:
143  // unsupported operator type
144  device_assert(false);
145  }
146  // copy the values into the output container
147  for (int fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
148  {
149  output_.access(fieldOrdinal,pointOrdinal,0) = field_values_at_point(fieldOrdinal);
150  }
151  }
152 
153  // Provide the shared memory capacity.
154  // This function takes the team_size as an argument,
155  // which allows team_size-dependent allocations.
156  size_t team_shmem_size (int team_size) const
157  {
158  // we want to use shared memory to create a fast buffer that we can use for basis computations
159  size_t shmem_size = 0;
160  if (fad_size_output_ > 0)
161  shmem_size += OutputScratchView::shmem_size(numFields_, fad_size_output_);
162  else
163  shmem_size += OutputScratchView::shmem_size(numFields_);
164 
165  return shmem_size;
166  }
167  };
168 
180  template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
181  typename OutputScalar = double,
182  typename PointScalar = double>
184  : public Basis<ExecutionSpace,OutputScalar,PointScalar>
185  {
186  public:
187  using OrdinalTypeArray1DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray1DHost;
188  using OrdinalTypeArray2DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray2DHost;
189 
193  protected:
194  int polyOrder_; // the maximum order of the polynomial
195  public:
204  LegendreBasis_HVOL_LINE(int polyOrder)
205  :
206  polyOrder_(polyOrder)
207  {
208  this->basisCardinality_ = polyOrder+1;
209  this->basisDegree_ = polyOrder;
210  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
211  this->basisType_ = BASIS_FEM_HIERARCHICAL;
212  this->basisCoordinates_ = COORDINATES_CARTESIAN;
213  this->functionSpace_ = FUNCTION_SPACE_HVOL;
214 
215  const int degreeLength = 1;
216  this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost("Integrated Legendre H(grad) line polynomial degree lookup", this->basisCardinality_, degreeLength);
217 
218  for (int i=0; i<this->basisCardinality_; i++)
219  {
220  // for H(vol) line, first basis member is constant, second is first-degree, etc.
221  this->fieldOrdinalPolynomialDegree_(i,0) = i;
222  }
223 
224  // initialize tags
225  {
226  const auto & cardinality = this->basisCardinality_;
227 
228  // Basis-dependent initializations
229  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
230  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
231  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
232  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
233 
234  OrdinalTypeArray1DHost tagView("tag view", cardinality*tagSize);
235 
236  for (ordinal_type i=0;i<cardinality;++i) {
237  tagView(i*tagSize+0) = 1; // edge dof
238  tagView(i*tagSize+1) = 0; // edge id
239  tagView(i*tagSize+2) = i; // local dof id
240  tagView(i*tagSize+3) = cardinality; // total number of dofs in this edge
241  }
242 
243  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
244  // tags are constructed on host
245  this->setOrdinalTagData(this->tagToOrdinal_,
246  this->ordinalToTag_,
247  tagView,
248  this->basisCardinality_,
249  tagSize,
250  posScDim,
251  posScOrd,
252  posDfOrd);
253  }
254  }
255 
256  // since the getValues() below only overrides the FEM variant, we specify that
257  // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
258  // (It's an error to use the FVD variant on this basis.)
260 
265  virtual
266  const char*
267  getName() const override {
268  return "Intrepid2_LegendreBasis_HVOL_LINE";
269  }
270 
289  virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
290  const EOperator operatorType = OPERATOR_VALUE ) const override
291  {
292  auto numPoints = inputPoints.extent_int(0);
293 
295 
296  FunctorType functor(outputValues, inputPoints, polyOrder_, operatorType);
297 
298  const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
299  const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
300  const int vectorSize = std::max(outputVectorSize,pointVectorSize);
301  const int teamSize = 1; // because of the way the basis functions are computed, we don't have a second level of parallelism...
302 
303  auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
304  Kokkos::parallel_for( policy , functor, "Hierarchical_HVOL_LINE_Functor");
305  }
306  };
307 } // end namespace Intrepid2
308 
309 #endif /* Intrepid2_LegendreBasis_HVOL_LINE_h */
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
Functor for computing values for the LegendreBasis_HVOL_LINE class.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
EFunctionSpace functionSpace_
The function space in which the basis is defined.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Free functions, callable from device code, that implement various polynomials useful in basis definit...
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Implementation of an assert that can safely be called from device code.
Basis defining Legendre basis on the line, a polynomial subspace of L^2 (a.k.a. H(vol)) on the line...
virtual const char * getName() const override
Returns basis name.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
void setOrdinalTagData(OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
OrdinalTypeArray2DHost ordinalToTag_
&quot;true&quot; if tagToOrdinal_ and ordinalToTag_ have been initialized
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...