Intrepid2
Intrepid2_IntegratedLegendreBasis_HGRAD_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_IntegratedLegendreBasis_HGRAD_LINE_h
50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_LINE_h
51 
52 #include <Kokkos_View.hpp>
53 #include <Kokkos_DynRankView.hpp>
54 
55 #include <Intrepid2_config.h>
56 
58 #include "Intrepid2_Utils.hpp"
59 
60 namespace Intrepid2
61 {
67  template<class ExecutionSpace, class OutputScalar, class PointScalar,
68  class OutputFieldType, class InputPointsType>
70  {
71  using ScratchSpace = typename ExecutionSpace::scratch_memory_space;
72  using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
73  using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
74 
75  using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
76  using TeamMember = typename TeamPolicy::member_type;
77 
78  EOperator opType_; // OPERATOR_VALUE or OPERATOR_GRAD
79 
80  OutputFieldType output_; // F,P
81  InputPointsType inputPoints_; // P,D
82 
83  int polyOrder_;
84  bool defineVertexFunctions_;
85  int numFields_, numPoints_;
86 
87  size_t fad_size_output_;
88 
89  Hierarchical_HGRAD_LINE_Functor(EOperator opType, OutputFieldType output, InputPointsType inputPoints,
90  int polyOrder, bool defineVertexFunctions)
91  : opType_(opType), output_(output), inputPoints_(inputPoints),
92  polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
93  fad_size_output_(getScalarDimensionForView(output))
94  {
95  numFields_ = output.extent_int(0);
96  numPoints_ = output.extent_int(1);
97  INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, "point counts need to match!");
98  INTREPID2_TEST_FOR_EXCEPTION(numFields_ != polyOrder_+1, std::invalid_argument, "output field size does not match basis cardinality");
99  }
100 
101  KOKKOS_INLINE_FUNCTION
102  void operator()( const TeamMember & teamMember ) const
103  {
104  auto pointOrdinal = teamMember.league_rank();
105  OutputScratchView field_values_at_point;
106  if (fad_size_output_ > 0) {
107  field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_, fad_size_output_);
108  }
109  else {
110  field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_);
111  }
112 
113  const auto & input_x = inputPoints_(pointOrdinal,0);
114  const bool taking_derivative = (opType_ != OPERATOR_VALUE);
115  const bool callingShiftedScaledLegendre = (opType_ == OPERATOR_VALUE) || (opType_ == OPERATOR_GRAD) || (opType_ == OPERATOR_D1);
116 
117  // shiftedScaledIntegratedLegendreValues{_dx} expects x in [0,1]
118  const PointScalar x = callingShiftedScaledLegendre ? PointScalar((input_x + 1.0)/2.0) : PointScalar(input_x);
119  const double legendreScaling = 1.0;
120  const double outputScaling = taking_derivative ? 0.5 : 1.0; // output scaling -- 0.5 if we take derivatives, 1.0 otherwise
121 
122  switch (opType_)
123  {
124  case OPERATOR_VALUE:
125  // field values are integrated Legendre polynomials, except for the first and second field,
126  // which may be 1 and x or x and 1-x, depending on whether the vertex compatibility flag is set.
127  Polynomials::shiftedScaledIntegratedLegendreValues(field_values_at_point, polyOrder_, x, legendreScaling);
128 
129  // note that because shiftedScaledIntegratedLegendreValues determines field values recursively, there is not much
130  // opportunity at that level for further parallelism
131 
132  if (defineVertexFunctions_)
133  {
134  field_values_at_point(0) = 1. - x;
135  field_values_at_point(1) = x;
136  }
137  break;
138  case OPERATOR_GRAD:
139  case OPERATOR_D1:
140  // field values are Legendre polynomials, except for the first and second field,
141  // which may be 0 and 1 or -1 and 1, depending on whether the vertex compatibility flag is set.
142  Polynomials::shiftedScaledIntegratedLegendreValues_dx(field_values_at_point, polyOrder_, x, legendreScaling);
143 
144  // note that because shiftedScaledIntegratedLegendreValues_dx determines field values recursively, there is not much
145  // opportunity at that level for further parallelism
146 
147  if (defineVertexFunctions_)
148  {
149  field_values_at_point(0) = -1.0; // derivative of 1-x
150  field_values_at_point(1) = 1.0; // derivative of x
151  }
152  break;
153  case OPERATOR_D2:
154  case OPERATOR_D3:
155  case OPERATOR_D4:
156  case OPERATOR_D5:
157  case OPERATOR_D6:
158  case OPERATOR_D7:
159  case OPERATOR_D8:
160  case OPERATOR_D9:
161  case OPERATOR_D10:
162  {
163  auto derivativeOrder = getOperatorOrder(opType_) - 1;
164  Polynomials::legendreDerivativeValues(field_values_at_point, polyOrder_, x, derivativeOrder);
165 
166  // L_i is defined in terms of an integral of P_(i-1), so we need to shift the values by 1
167  if (numFields_ >= 3)
168  {
169  OutputScalar Pn_minus_one = field_values_at_point(1);
170  for (int fieldOrdinal=2; fieldOrdinal<numFields_; fieldOrdinal++)
171  {
172  OutputScalar Pn = field_values_at_point(fieldOrdinal);
173  field_values_at_point(fieldOrdinal) = Pn_minus_one;
174  Pn_minus_one = Pn;
175  }
176  }
177  if (numFields_ >= 1) field_values_at_point(0) = 0.0;
178  if (numFields_ >= 2) field_values_at_point(1) = 0.0;
179  // legendreDerivativeValues works on [-1,1], so no per-derivative scaling is necessary
180  // however, there is a factor of 0.5 that comes from the scaling of the Legendre polynomials prior to integration
181  // in the shiftedScaledIntegratedLegendreValues -- the first derivative of our integrated polynomials is 0.5 times the Legendre polynomial
182  break;
183  }
184  default:
185  // unsupported operator type
186  device_assert(false);
187  }
188 
189  // copy the values into the output container
190  for (int fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
191  {
192  // access() allows us to write one line that applies both to gradient (for which outputValues has rank 3, but third rank has only one entry) and to value (rank 2)
193  output_.access(fieldOrdinal,pointOrdinal,0) = outputScaling * field_values_at_point(fieldOrdinal);
194  }
195  }
196 
197  // Provide the shared memory capacity.
198  // This function takes the team_size as an argument,
199  // which allows team_size-dependent allocations.
200  size_t team_shmem_size (int team_size) const
201  {
202  // we want to use shared memory to create a fast buffer that we can use for basis computations
203  size_t shmem_size = 0;
204  if (fad_size_output_ > 0)
205  shmem_size += OutputScratchView::shmem_size(numFields_, fad_size_output_);
206  else
207  shmem_size += OutputScratchView::shmem_size(numFields_);
208 
209  return shmem_size;
210  }
211  };
212 
230  template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
231  typename OutputScalar = double,
232  typename PointScalar = double,
233  bool defineVertexFunctions = true, // if defineVertexFunctions is true, first and second basis functions are x and 1-x. Otherwise, they are 1 and x.
234  bool useMinusOneToOneReferenceElement = true> // if useMinusOneToOneReferenceElement is true, basis is define on [-1,1]. Otherwise, [0,1].
236  : public Basis<ExecutionSpace,OutputScalar,PointScalar>
237  {
238  public:
239  using OrdinalTypeArray1DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray1DHost;
240  using OrdinalTypeArray2DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray2DHost;
241 
245  protected:
246  int polyOrder_; // the maximum order of the polynomial
247  bool defineVertexFunctions_; // if true, first and second basis functions are x and 1-x. Otherwise, they are 1 and x.
248  public:
260  :
261  polyOrder_(polyOrder)
262  {
263  this->basisCardinality_ = polyOrder+1;
264  this->basisDegree_ = polyOrder;
265  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
266  this->basisType_ = BASIS_FEM_HIERARCHICAL;
267  this->basisCoordinates_ = COORDINATES_CARTESIAN;
268  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
269 
270  const int degreeLength = 1;
271  this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost("Integrated Legendre H(grad) line polynomial degree lookup", this->basisCardinality_, degreeLength);
272 
273  for (int i=0; i<this->basisCardinality_; i++)
274  {
275  // for H(grad) line, if defineVertexFunctions is false, first basis member is constant, second is first-degree, etc.
276  // if defineVertexFunctions is true, then the only difference is that the entry is also degree 1
277  this->fieldOrdinalPolynomialDegree_(i,0) = i;
278  }
279  if (defineVertexFunctions)
280  {
281  this->fieldOrdinalPolynomialDegree_(0,0) = 1;
282  }
283 
284  // initialize tags
285  {
286  const auto & cardinality = this->basisCardinality_;
287 
288  // Basis-dependent initializations
289  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
290  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
291  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
292  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
293 
294  OrdinalTypeArray1DHost tagView("tag view", cardinality*tagSize);
295 
296  if (defineVertexFunctions) {
297  {
298  const ordinal_type v0 = 0;
299  tagView(v0*tagSize+0) = 0; // vertex dof
300  tagView(v0*tagSize+1) = 0; // vertex id
301  tagView(v0*tagSize+2) = 0; // local dof id
302  tagView(v0*tagSize+3) = 1; // total number of dofs in this vertex
303 
304  const ordinal_type v1 = 1;
305  tagView(v1*tagSize+0) = 0; // vertex dof
306  tagView(v1*tagSize+1) = 1; // vertex id
307  tagView(v1*tagSize+2) = 0; // local dof id
308  tagView(v1*tagSize+3) = 1; // total number of dofs in this vertex
309 
310  const ordinal_type iend = cardinality - 2;
311  for (ordinal_type i=0;i<iend;++i) {
312  const auto e = i + 2;
313  tagView(e*tagSize+0) = 1; // edge dof
314  tagView(e*tagSize+1) = 0; // edge id
315  tagView(e*tagSize+2) = i; // local dof id
316  tagView(e*tagSize+3) = iend; // total number of dofs in this edge
317  }
318  }
319  } else {
320  for (ordinal_type i=0;i<cardinality;++i) {
321  tagView(i*tagSize+0) = 1; // edge dof
322  tagView(i*tagSize+1) = 0; // edge id
323  tagView(i*tagSize+2) = i; // local dof id
324  tagView(i*tagSize+3) = cardinality; // total number of dofs in this edge
325  }
326  }
327 
328  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
329  // tags are constructed on host
330  this->setOrdinalTagData(this->tagToOrdinal_,
331  this->ordinalToTag_,
332  tagView,
333  this->basisCardinality_,
334  tagSize,
335  posScDim,
336  posScOrd,
337  posDfOrd);
338  }
339  }
340 
345  const char* getName() const override {
346  return "Intrepid2_IntegratedLegendreBasis_HGRAD_LINE";
347  }
348 
351  virtual bool requireOrientation() const override {
352  return false;
353  }
354 
355  // since the getValues() below only overrides the FEM variant, we specify that
356  // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
357  // (It's an error to use the FVD variant on this basis.)
359 
378  virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
379  const EOperator operatorType = OPERATOR_VALUE ) const override
380  {
381  auto numPoints = inputPoints.extent_int(0);
382 
384 
385  FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
386 
387  const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
388  const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
389  const int vectorSize = std::max(outputVectorSize,pointVectorSize);
390  const int teamSize = 1; // because of the way the basis functions are computed, we don't have a second level of parallelism...
391 
392  auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
393  Kokkos::parallel_for( policy , functor, "Hierarchical_HGRAD_LINE_Functor");
394  }
395  };
396 } // end namespace Intrepid2
397 
398 #endif /* Intrepid2_IntegratedLegendreBasis_HGRAD_LINE_h */
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
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 (...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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.
Header function for Intrepid2::Util class and other utility functions.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
virtual bool requireOrientation() const override
True if orientation is required.
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.
Functor for computing values for the IntegratedLegendreBasis_HGRAD_LINE class.
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...