Intrepid2
Intrepid2_Polynomials.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_Polynomials_h
50 #define Intrepid2_Polynomials_h
51 
52 #include "Intrepid2_Polylib.hpp"
53 #include "Intrepid2_Types.hpp"
54 #include "Intrepid2_Utils.hpp"
55 
56 namespace Intrepid2
57 {
58  namespace Polynomials
59  {
60  /*
61  These polynomials are supplemental to those defined in the Polylib class; there is some overlap.
62  We actually take advantage of the overlap in our verification tests, using the Polylib functions
63  to verify the ones defined here. Our interface here is a little simpler, and the functions are a little less
64  general than those in Polylib.
65 
66  We define some polynomial functions that are useful in a variety of contexts.
67  In particular, the (integrated) Legendre and Jacobi polynomials are useful in defining
68  hierarchical bases. See in particular:
69 
70  Federico Fuentes, Brendan Keith, Leszek Demkowicz, Sriram Nagaraj.
71  "Orientation embedded high order shape functions for the exact sequence elements of all shapes."
72  Computers & Mathematics with Applications, Volume 70, Issue 4, 2015, Pages 353-458, ISSN 0898-1221.
73  https://doi.org/10.1016/j.camwa.2015.04.027.
74 
75  In this implementation, we take care to make minimal assumptions on both the containers
76  and the scalar type. The containers need to support a one-argument operator() for assignment and/or
77  lookup (as appropriate). The scalar type needs to support a cast from Intrepid2::ordinal_type, as well
78  as standard arithmetic operations. In particular, both 1-rank Kokkos::View and Kokkos::DynRankView are
79  supported, as are C++ floating point types and Sacado scalar types.
80  */
81 
88  template<typename OutputValueViewType, typename ScalarType>
89  KOKKOS_INLINE_FUNCTION void legendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
90  {
91  if (n >= 0) outputValues(0) = 1.0;
92  if (n >= 1) outputValues(1) = x;
93  for (int i=2; i<=n; i++)
94  {
95  const ScalarType i_scalar = ScalarType(i);
96  outputValues(i) = (2. - 1. / i_scalar) * x * outputValues(i-1) - (1. - 1. / i_scalar) * outputValues(i-2);
97  }
98  }
99 
107  template<typename OutputValueViewType, typename ScalarType>
108  KOKKOS_INLINE_FUNCTION void legendreDerivativeValues(OutputValueViewType outputValues, const OutputValueViewType legendreValues, Intrepid2::ordinal_type n, ScalarType x)
109  {
110  if (n >= 0) outputValues(0) = 0.0;
111  if (n >= 1) outputValues(1) = 1.0;
112  for (int i=2; i<=n; i++)
113  {
114  const ScalarType i_scalar = ScalarType(i);
115  outputValues(i) = outputValues(i-2) + (2. * i_scalar - 1.) * legendreValues(i-1);
116  }
117  }
118 
119  // derivative values can be computed using the Legendre values
120  // n: number of Legendre polynomials' derivative values to compute. outputValues must have at least n+1 entries
121  // x: value in [-1, 1]
122  // dn: number of derivatives to take. Must be >= 1.
130  template<typename OutputValueViewType, typename PointScalarType>
131  KOKKOS_INLINE_FUNCTION void legendreDerivativeValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, PointScalarType x, Intrepid2::ordinal_type dn)
132  {
133  const OutputValueViewType nullOutputScalarView;
134  const double alpha = 0;
135  const double beta = 0;
136  const int numPoints = 1;
137 
138  using Layout = typename NaturalLayoutForType<PointScalarType>::layout;
139 
140  using UnmanagedPointScalarView = Kokkos::View<PointScalarType*, Layout, Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
141  UnmanagedPointScalarView pointView = UnmanagedPointScalarView(&x,numPoints);
142 
143  for (int i=0; i<=n; i++)
144  {
145  auto jacobiValue = Kokkos::subview(outputValues,Kokkos::pair<Intrepid2::ordinal_type,Intrepid2::ordinal_type>(i,i+1));
146  jacobiValue(0) = 0.0;
147  Intrepid2::Polylib::Serial::JacobiPolynomial(numPoints, pointView, jacobiValue, nullOutputScalarView, i-dn, alpha+dn, beta+dn);
148 
149  double scaleFactor = 1.0;
150  for (int j=1; j<=dn; j++)
151  {
152  scaleFactor *= 0.5 * (j+alpha+beta+i);
153  }
154 
155  outputValues(i) = jacobiValue(0) * scaleFactor;
156  }
157  }
158 
168  template<typename OutputValueViewType, typename ScalarType>
169  KOKKOS_INLINE_FUNCTION void shiftedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
170  {
171  legendreValues(outputValues, n, 2.*x-1.);
172  }
173 
184  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
185  KOKKOS_INLINE_FUNCTION void shiftedScaledLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
186  {
187  using OutputScalar = typename OutputValueViewType::value_type;
188  OutputScalar two_x_minus_t = 2. * x - t;
189  OutputScalar t_squared = t * t;
190  if (n >= 0) outputValues(0) = 1.0;
191  if (n >= 1) outputValues(1) = two_x_minus_t;
192  for (int i=2; i<=n; i++)
193  {
194  const ScalarType one_over_i = 1.0 / ScalarType(i);
195  outputValues(i) = one_over_i * ( (2. *i - 1. ) * two_x_minus_t * outputValues(i-1) - (i - 1.) * t_squared * outputValues(i-2));
196  }
197  }
198 
212  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
213  KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
214  {
215  // reduced memory version: compute P_i in outputValues
216  shiftedScaledLegendreValues(outputValues,n,x,t);
217  // keep a copy of the last two P_i values around; update these before overwriting in outputValues
218  ScalarType P_i_minus_two, P_i_minus_one;
219  if (n >= 0) P_i_minus_two = outputValues(0);
220  if (n >= 1) P_i_minus_one = outputValues(1);
221 
222  if (n >= 0) outputValues(0) = 1.0;
223  if (n >= 1) outputValues(1) = x;
224  for (int i=2; i<=n; i++)
225  {
226  const ScalarType & P_i = outputValues(i); // define as P_i just for clarity of the code below
227  const ScalarType i_scalar = ScalarType(i);
228  ScalarType L_i = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
229 
230  // get the next values of P_{i-1} and P_{i-2} before overwriting the P_i value
231  P_i_minus_two = P_i_minus_one;
232  P_i_minus_one = P_i;
233 
234  // overwrite P_i value
235  outputValues(i) = L_i;
236  }
237  }
238 
253  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
254  KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, const OutputValueViewType shiftedScaledLegendreValues,
255  Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
256  {
257  // reduced flops version: rely on previously computed P_i
258  if (n >= 0) outputValues(0) = 1.0;
259  if (n >= 1) outputValues(1) = x;
260  for (int i=2; i<=n; i++)
261  {
262  const ScalarType & P_i = shiftedScaledLegendreValues(i); // define as P_i just for clarity of the code below
263  const ScalarType & P_i_minus_two = shiftedScaledLegendreValues(i-2);
264  const ScalarType i_scalar = ScalarType(i);
265  outputValues(i) = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
266  }
267  }
268 
269  // the below implementation is commented out for now to guard a certain confusion.
270  // the integratedLegendreValues() implementation below is implemented in such a way as to agree, modulo a coordinate transformation, with the
271  // shiftedScaledIntegratedLegendreValues() above. Since the latter really is an integral of shiftedScaledLegendreValues(), the former isn't
272  // actually an integral of the (unshifted, unscaled) Legendre polynomials.
273 // template<typename OutputValueViewType, typename ScalarType>
274 // KOKKOS_INLINE_FUNCTION void integratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
275 // {
276 // // to reduce memory requirements, compute P_i in outputValues
277 // legendreValues(outputValues,n,x);
278 // // keep a copy of the last two P_i values around; update these before overwriting in outputValues
279 // ScalarType P_i_minus_two, P_i_minus_one;
280 // if (n >= 0) P_i_minus_two = outputValues(0);
281 // if (n >= 1) P_i_minus_one = outputValues(1);
282 //
283 // if (n >= 0) outputValues(0) = 1.0;
284 // if (n >= 1) outputValues(1) = (x + 1.0) / 2.0;
285 // for (int i=2; i<=n; i++)
286 // {
287 // const ScalarType & P_i = outputValues(i); // define as P_i just for clarity of the code below
288 // const ScalarType i_scalar = ScalarType(i);
289 // ScalarType L_i = (P_i - P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
290 //
291 // // get the next values of P_{i-1} and P_{i-2} before overwriting the P_i value
292 // P_i_minus_two = P_i_minus_one;
293 // P_i_minus_one = P_i;
294 //
295 // // overwrite P_i value
296 // outputValues(i) = L_i;
297 // }
298 // }
299 
307  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
308  KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dx(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
309  {
310  if (n >= 0) outputValues(0) = 0.0;
311  if (n >= 1) outputValues(1) = 1.0;
312  if (n >= 2) outputValues(2) = 2. * x - t;
313  for (int i=2; i<=n-1; i++)
314  {
315  const ScalarType one_over_i = 1.0 / ScalarType(i);
316  outputValues(i+1) = one_over_i * (2. * i - 1.) * (2. * x - t) * outputValues(i) - one_over_i * (i - 1.0) * t * t * outputValues(i-1);
317  }
318  }
319 
329  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
330  KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
331  {
332  // memory-conserving version -- place the Legendre values in the final output container
333  shiftedScaledLegendreValues(outputValues, n, x, t);
334 
335  ScalarType P_i_minus_2 = outputValues(0);
336  ScalarType P_i_minus_1 = outputValues(1);
337 
338  if (n >= 0) outputValues(0) = 0.0;
339  if (n >= 1) outputValues(1) = 0.0;
340 
341  for (int i=2; i<=n; i++)
342  {
343  const ScalarType L_i_dt = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
344 
345  P_i_minus_2 = P_i_minus_1;
346  P_i_minus_1 = outputValues(i);
347 
348  outputValues(i) = L_i_dt;
349  }
350  }
351 
362  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
363  KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, const OutputValueViewType shiftedScaledLegendreValues,
364  Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
365  {
366  // reduced flops version: rely on previously computed P_i
367  if (n >= 0) outputValues(0) = 0.0;
368  if (n >= 1) outputValues(1) = 0.0;
369  for (int i=2; i<=n; i++)
370  {
371  const ScalarType & P_i_minus_1 = shiftedScaledLegendreValues(i-1); // define as P_i just for clarity of the code below
372  const ScalarType & P_i_minus_2 = shiftedScaledLegendreValues(i-2);
373  outputValues(i) = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
374  }
375  }
376 
391  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
392  KOKKOS_INLINE_FUNCTION void shiftedScaledJacobiValues(OutputValueViewType outputValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
393  {
394  ScalarType two_x_minus_t = 2. * x - t;
395  ScalarTypeForScaling alpha_squared_t = alpha * alpha * t;
396 
397  if (n >= 0) outputValues(0) = 1.0;
398  if (n >= 1) outputValues(1) = two_x_minus_t + alpha * x;
399 
400  for (int i=2; i<=n; i++)
401  {
402  const ScalarType & P_i_minus_one = outputValues(i-1); // define as P_i just for clarity of the code below
403  const ScalarType & P_i_minus_two = outputValues(i-2);
404 
405  double a_i = (2. * i) * (i + alpha) * (2. * i + alpha - 2.);
406  double b_i = 2. * i + alpha - 1.;
407  double c_i = (2. * i + alpha) * (2. * i + alpha - 2.);
408  double d_i = 2. * (i + alpha - 1.) * (i - 1.) * (2. * i + alpha);
409 
410  outputValues(i) = (b_i / a_i) * (c_i * two_x_minus_t + alpha_squared_t) * P_i_minus_one - (d_i / a_i) * t * t * P_i_minus_two;
411  }
412  }
413 
431  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
432  KOKKOS_INLINE_FUNCTION void integratedJacobiValues(OutputValueViewType outputValues, const OutputValueViewType jacobiValues,
433  double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
434  {
435  // reduced flops version: rely on previously computed P_i
436  if (n >= 0) outputValues(0) = 1.0;
437  if (n >= 1) outputValues(1) = x;
438 
439  ScalarType t_squared = t * t;
440  for (int i=2; i<=n; i++)
441  {
442  const ScalarType & P_i = jacobiValues(i); // define as P_i just for clarity of the code below
443  const ScalarType & P_i_minus_1 = jacobiValues(i-1);
444  const ScalarType & P_i_minus_2 = jacobiValues(i-2);
445 
446  double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha ));
447  double b_i = alpha / ((2. * i + alpha - 2.) * (2. * i + alpha ));
448  double c_i = (i - 1.) / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
449 
450  outputValues(i) = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
451  }
452  }
453 
471  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
472  KOKKOS_INLINE_FUNCTION void integratedJacobiValues(OutputValueViewType outputValues,
473  double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
474  {
475  // memory-conserving version -- place the Jacobi values in the final output container
476  shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
477 
478  ScalarType P_i_minus_2 = outputValues(0);
479  ScalarType P_i_minus_1 = outputValues(1);
480 
481  if (n >= 0) outputValues(0) = 1.0;
482  if (n >= 1) outputValues(1) = x;
483 
484  ScalarType t_squared = t * t;
485  for (int i=2; i<=n; i++)
486  {
487  const ScalarType & P_i = outputValues(i);
488 
489  double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha ));
490  double b_i = alpha / ((2. * i + alpha - 2.) * (2. * i + alpha ));
491  double c_i = (i - 1.) / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
492 
493  ScalarType L_i = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
494 
495  P_i_minus_2 = P_i_minus_1;
496  P_i_minus_1 = P_i;
497 
498  outputValues(i) = L_i;
499  }
500  }
501 
502  // x derivative of integrated Jacobi is just Jacobi
503  // only distinction is in the index -- outputValues indices are shifted by 1 relative to jacobiValues, above
511  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
512  KOKKOS_INLINE_FUNCTION void integratedJacobiValues_dx(OutputValueViewType outputValues,
513  double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
514  {
515  // rather than repeating the somewhat involved implementation of jacobiValues here,
516  // call with (n-1), and then move values accordingly
517  shiftedScaledJacobiValues(outputValues, alpha, n-1, x, t);
518 
519  // forward implementation
520  ScalarType nextValue = 0.0;
521  ScalarType nextNextValue = 0.0;
522  for (int i=0; i<=n-1; i++)
523  {
524  nextNextValue = outputValues(i);
525  outputValues(i) = nextValue;
526  nextValue = nextNextValue;
527  }
528  outputValues(n-1) = nextValue;
529  }
530 
541  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
542  KOKKOS_INLINE_FUNCTION void integratedJacobiValues_dt(OutputValueViewType outputValues, const OutputValueViewType jacobiValues,
543  double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
544  {
545  // reduced flops version: rely on previously computed P_i
546  if (n >= 0) outputValues(0) = 0.0;
547  if (n >= 1) outputValues(1) = 0.0;
548  for (int i=2; i<=n; i++)
549  {
550  const ScalarType & P_i_minus_1 = jacobiValues(i-1); // define as P_i just for clarity of the code below
551  const ScalarType & P_i_minus_2 = jacobiValues(i-2);
552  outputValues(i) = - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
553  }
554  }
555 
566  template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
567  KOKKOS_INLINE_FUNCTION void integratedJacobiValues_dt(OutputValueViewType outputValues,
568  double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
569  {
570  // memory-conserving version -- place the Jacobi values in the final output container
571  shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
572 
573  ScalarType P_i_minus_2 = outputValues(0);
574  ScalarType P_i_minus_1 = outputValues(1);
575 
576  if (n >= 0) outputValues(0) = 0.0;
577  if (n >= 1) outputValues(1) = 0.0;
578 
579  for (int i=2; i<=n; i++)
580  {
581  const ScalarType L_i_dt = - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
582 
583  P_i_minus_2 = P_i_minus_1;
584  P_i_minus_1 = outputValues(i);
585 
586  outputValues(i) = L_i_dt;
587  }
588  }
589  } // namespace Polynomials
590 } // namespace Intrepid2
591 
592 #endif /* Intrepid2_Polynomials_h */
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation...
Header function for Intrepid2::Util class and other utility functions.
static KOKKOS_INLINE_FUNCTION void JacobiPolynomial(const ordinal_type np, const zViewType z, polyiViewType poly_in, polydViewType polyd, const ordinal_type n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Contains definitions of custom data types in Intrepid2.