Intrepid2
Intrepid2_Polynomials.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14
15#ifndef Intrepid2_Polynomials_h
16#define Intrepid2_Polynomials_h
17
18#include "Intrepid2_Polylib.hpp"
19#include "Intrepid2_Types.hpp"
20#include "Intrepid2_Utils.hpp"
21
22namespace Intrepid2
23{
24 namespace Polynomials
25 {
26 /*
27 These polynomials are supplemental to those defined in the Polylib class; there is some overlap.
28 We actually take advantage of the overlap in our verification tests, using the Polylib functions
29 to verify the ones defined here. Our interface here is a little simpler, and the functions are a little less
30 general than those in Polylib.
31
32 We define some polynomial functions that are useful in a variety of contexts.
33 In particular, the (integrated) Legendre and Jacobi polynomials are useful in defining
34 hierarchical bases. See in particular:
35
36 Federico Fuentes, Brendan Keith, Leszek Demkowicz, Sriram Nagaraj.
37 "Orientation embedded high order shape functions for the exact sequence elements of all shapes."
38 Computers & Mathematics with Applications, Volume 70, Issue 4, 2015, Pages 353-458, ISSN 0898-1221.
39 https://doi.org/10.1016/j.camwa.2015.04.027.
40
41 In this implementation, we take care to make minimal assumptions on both the containers
42 and the scalar type. The containers need to support a one-argument operator() for assignment and/or
43 lookup (as appropriate). The scalar type needs to support a cast from Intrepid2::ordinal_type, as well
44 as standard arithmetic operations. In particular, both 1-rank Kokkos::View and Kokkos::DynRankView are
45 supported, as are C++ floating point types and Sacado scalar types.
46 */
47
54 template<typename OutputValueViewType, typename ScalarType>
55 KOKKOS_INLINE_FUNCTION void legendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
56 {
57 if (n >= 0) outputValues(0) = 1.0;
58 if (n >= 1) outputValues(1) = x;
59 for (int i=2; i<=n; i++)
60 {
61 const ScalarType i_scalar = ScalarType(i);
62 outputValues(i) = (2. - 1. / i_scalar) * x * outputValues(i-1) - (1. - 1. / i_scalar) * outputValues(i-2);
63 }
64 }
65
73 template<typename OutputValueViewType, typename ScalarType>
74 KOKKOS_INLINE_FUNCTION void legendreDerivativeValues(OutputValueViewType outputValues, const OutputValueViewType legendreValues, Intrepid2::ordinal_type n, ScalarType x)
75 {
76 if (n >= 0) outputValues(0) = 0.0;
77 if (n >= 1) outputValues(1) = 1.0;
78 for (int i=2; i<=n; i++)
79 {
80 const ScalarType i_scalar = ScalarType(i);
81 outputValues(i) = outputValues(i-2) + (2. * i_scalar - 1.) * legendreValues(i-1);
82 }
83 }
84
85 // derivative values can be computed using the Legendre values
86 // n: number of Legendre polynomials' derivative values to compute. outputValues must have at least n+1 entries
87 // x: value in [-1, 1]
88 // dn: number of derivatives to take. Must be >= 1.
96 template<typename OutputValueViewType, typename PointScalarType>
97 KOKKOS_INLINE_FUNCTION void legendreDerivativeValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, PointScalarType x, Intrepid2::ordinal_type dn)
98 {
99 const OutputValueViewType nullOutputScalarView;
100 const double alpha = 0;
101 const double beta = 0;
102 const int numPoints = 1;
103
104 using Layout = typename NaturalLayoutForType<PointScalarType>::layout;
105
106 using UnmanagedPointScalarView = Kokkos::View<PointScalarType*, Layout, Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
107 UnmanagedPointScalarView pointView = UnmanagedPointScalarView(&x,numPoints);
108
109 for (int i=0; i<=n; i++)
110 {
111 auto jacobiValue = Kokkos::subview(outputValues,Kokkos::pair<Intrepid2::ordinal_type,Intrepid2::ordinal_type>(i,i+1));
112 jacobiValue(0) = 0.0;
113 Intrepid2::Polylib::Serial::JacobiPolynomial(numPoints, pointView, jacobiValue, nullOutputScalarView, i-dn, alpha+dn, beta+dn);
114
115 double scaleFactor = 1.0;
116 for (int j=1; j<=dn; j++)
117 {
118 scaleFactor *= 0.5 * (j+alpha+beta+i);
119 }
120
121 outputValues(i) = jacobiValue(0) * scaleFactor;
122 }
123 }
124
134 template<typename OutputValueViewType, typename ScalarType>
135 KOKKOS_INLINE_FUNCTION void shiftedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
136 {
137 legendreValues(outputValues, n, 2.*x-1.);
138 }
139
150 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
151 KOKKOS_INLINE_FUNCTION void shiftedScaledLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
152 {
153 using OutputScalar = typename OutputValueViewType::value_type;
154 OutputScalar two_x_minus_t = 2. * x - t;
155 OutputScalar t_squared = t * t;
156 if (n >= 0) outputValues(0) = 1.0;
157 if (n >= 1) outputValues(1) = two_x_minus_t;
158 for (int i=2; i<=n; i++)
159 {
160 const ScalarType one_over_i = 1.0 / ScalarType(i);
161 outputValues(i) = one_over_i * ( (2. *i - 1. ) * two_x_minus_t * outputValues(i-1) - (i - 1.) * t_squared * outputValues(i-2));
162 }
163 }
164
178 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
179 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
180 {
181 // reduced memory version: compute P_i in outputValues
182 shiftedScaledLegendreValues(outputValues,n,x,t);
183 // keep a copy of the last two P_i values around; update these before overwriting in outputValues
184 ScalarType P_i_minus_two, P_i_minus_one;
185 if (n >= 0) P_i_minus_two = outputValues(0);
186 if (n >= 1) P_i_minus_one = outputValues(1);
187
188 if (n >= 0) outputValues(0) = 1.0;
189 if (n >= 1) outputValues(1) = x;
190 for (int i=2; i<=n; i++)
191 {
192 const ScalarType & P_i = outputValues(i); // define as P_i just for clarity of the code below
193 const ScalarType i_scalar = ScalarType(i);
194 ScalarType L_i = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
195
196 // get the next values of P_{i-1} and P_{i-2} before overwriting the P_i value
197 P_i_minus_two = P_i_minus_one;
198 P_i_minus_one = P_i;
199
200 // overwrite P_i value
201 outputValues(i) = L_i;
202 }
203 }
204
219 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
220 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, const OutputValueViewType shiftedScaledLegendreValues,
221 Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
222 {
223 // reduced flops version: rely on previously computed P_i
224 if (n >= 0) outputValues(0) = 1.0;
225 if (n >= 1) outputValues(1) = x;
226 for (int i=2; i<=n; i++)
227 {
228 const ScalarType & P_i = shiftedScaledLegendreValues(i); // define as P_i just for clarity of the code below
229 const ScalarType & P_i_minus_two = shiftedScaledLegendreValues(i-2);
230 const ScalarType i_scalar = ScalarType(i);
231 outputValues(i) = (P_i - t * t * P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
232 }
233 }
234
235 // the below implementation is commented out for now to guard a certain confusion.
236 // the integratedLegendreValues() implementation below is implemented in such a way as to agree, modulo a coordinate transformation, with the
237 // shiftedScaledIntegratedLegendreValues() above. Since the latter really is an integral of shiftedScaledLegendreValues(), the former isn't
238 // actually an integral of the (unshifted, unscaled) Legendre polynomials.
239// template<typename OutputValueViewType, typename ScalarType>
240// KOKKOS_INLINE_FUNCTION void integratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
241// {
242// // to reduce memory requirements, compute P_i in outputValues
243// legendreValues(outputValues,n,x);
244// // keep a copy of the last two P_i values around; update these before overwriting in outputValues
245// ScalarType P_i_minus_two, P_i_minus_one;
246// if (n >= 0) P_i_minus_two = outputValues(0);
247// if (n >= 1) P_i_minus_one = outputValues(1);
248//
249// if (n >= 0) outputValues(0) = 1.0;
250// if (n >= 1) outputValues(1) = (x + 1.0) / 2.0;
251// for (int i=2; i<=n; i++)
252// {
253// const ScalarType & P_i = outputValues(i); // define as P_i just for clarity of the code below
254// const ScalarType i_scalar = ScalarType(i);
255// ScalarType L_i = (P_i - P_i_minus_two) /( 2. * (2. * i_scalar - 1.));
256//
257// // get the next values of P_{i-1} and P_{i-2} before overwriting the P_i value
258// P_i_minus_two = P_i_minus_one;
259// P_i_minus_one = P_i;
260//
261// // overwrite P_i value
262// outputValues(i) = L_i;
263// }
264// }
265
273 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
274 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dx(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
275 {
276 if (n >= 0) outputValues(0) = 0.0;
277 if (n >= 1) outputValues(1) = 1.0;
278 if (n >= 2) outputValues(2) = 2. * x - t;
279 for (int i=2; i<=n-1; i++)
280 {
281 const ScalarType one_over_i = 1.0 / ScalarType(i);
282 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);
283 }
284 }
285
295 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
296 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
297 {
298 // memory-conserving version -- place the Legendre values in the final output container
299 shiftedScaledLegendreValues(outputValues, n, x, t);
300
301 ScalarType P_i_minus_2 = outputValues(0);
302 ScalarType P_i_minus_1 = outputValues(1);
303
304 if (n >= 0) outputValues(0) = 0.0;
305 if (n >= 1) outputValues(1) = 0.0;
306
307 for (int i=2; i<=n; i++)
308 {
309 const ScalarType L_i_dt = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
310
311 P_i_minus_2 = P_i_minus_1;
312 P_i_minus_1 = outputValues(i);
313
314 outputValues(i) = L_i_dt;
315 }
316 }
317
328 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
329 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, const OutputValueViewType shiftedScaledLegendreValues,
330 Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
331 {
332 // reduced flops version: rely on previously computed P_i
333 if (n >= 0) outputValues(0) = 0.0;
334 if (n >= 1) outputValues(1) = 0.0;
335 for (int i=2; i<=n; i++)
336 {
337 const ScalarType & P_i_minus_1 = shiftedScaledLegendreValues(i-1); // define as P_i just for clarity of the code below
338 const ScalarType & P_i_minus_2 = shiftedScaledLegendreValues(i-2);
339 outputValues(i) = -0.5 * (P_i_minus_1 + t * P_i_minus_2);
340 }
341 }
342
357 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
358 KOKKOS_INLINE_FUNCTION void shiftedScaledJacobiValues(OutputValueViewType outputValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
359 {
360 ScalarType two_x_minus_t = 2. * x - t;
361 ScalarTypeForScaling alpha_squared_t = alpha * alpha * t;
362
363 if (n >= 0) outputValues(0) = 1.0;
364 if (n >= 1) outputValues(1) = two_x_minus_t + alpha * x;
365
366 for (int i=2; i<=n; i++)
367 {
368 const ScalarType & P_i_minus_one = outputValues(i-1); // define as P_i just for clarity of the code below
369 const ScalarType & P_i_minus_two = outputValues(i-2);
370
371 double a_i = (2. * i) * (i + alpha) * (2. * i + alpha - 2.);
372 double b_i = 2. * i + alpha - 1.;
373 double c_i = (2. * i + alpha) * (2. * i + alpha - 2.);
374 double d_i = 2. * (i + alpha - 1.) * (i - 1.) * (2. * i + alpha);
375
376 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;
377 }
378 }
379
397 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
398 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues(OutputValueViewType outputValues, const OutputValueViewType jacobiValues,
399 double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
400 {
401 // reduced flops version: rely on previously computed P_i
402 if (n >= 0) outputValues(0) = 1.0;
403 if (n >= 1) outputValues(1) = x;
404
405 ScalarType t_squared = t * t;
406 for (int i=2; i<=n; i++)
407 {
408 const ScalarType & P_i = jacobiValues(i); // define as P_i just for clarity of the code below
409 const ScalarType & P_i_minus_1 = jacobiValues(i-1);
410 const ScalarType & P_i_minus_2 = jacobiValues(i-2);
411
412 double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha ));
413 double b_i = alpha / ((2. * i + alpha - 2.) * (2. * i + alpha ));
414 double c_i = (i - 1.) / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
415
416 outputValues(i) = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
417 }
418 }
419
437 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
438 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues(OutputValueViewType outputValues,
439 double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
440 {
441 // memory-conserving version -- place the Jacobi values in the final output container
442 shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
443
444 ScalarType P_i_minus_2 = outputValues(0);
445 ScalarType P_i_minus_1 = outputValues(1);
446
447 if (n >= 0) outputValues(0) = 1.0;
448 if (n >= 1) outputValues(1) = x;
449
450 ScalarType t_squared = t * t;
451 for (int i=2; i<=n; i++)
452 {
453 const ScalarType & P_i = outputValues(i);
454
455 double a_i = (i + alpha) / ((2. * i + alpha - 1.) * (2. * i + alpha ));
456 double b_i = alpha / ((2. * i + alpha - 2.) * (2. * i + alpha ));
457 double c_i = (i - 1.) / ((2. * i + alpha - 2.) * (2. * i + alpha - 1.));
458
459 ScalarType L_i = a_i * P_i + b_i * t * P_i_minus_1 - c_i * t_squared * P_i_minus_2;
460
461 P_i_minus_2 = P_i_minus_1;
462 P_i_minus_1 = P_i;
463
464 outputValues(i) = L_i;
465 }
466 }
467
468 // x derivative of integrated Jacobi is just Jacobi
469 // only distinction is in the index -- outputValues indices are shifted by 1 relative to jacobiValues, above
477 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
478 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues_dx(OutputValueViewType outputValues,
479 double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
480 {
481 // rather than repeating the somewhat involved implementation of jacobiValues here,
482 // call with (n-1), and then move values accordingly
483 shiftedScaledJacobiValues(outputValues, alpha, n-1, x, t);
484
485 // forward implementation
486 ScalarType nextValue = 0.0;
487 ScalarType nextNextValue = 0.0;
488 for (int i=0; i<=n-1; i++)
489 {
490 nextNextValue = outputValues(i);
491 outputValues(i) = nextValue;
492 nextValue = nextNextValue;
493 }
494 outputValues(n-1) = nextValue;
495 }
496
507 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
508 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues_dt(OutputValueViewType outputValues, const OutputValueViewType jacobiValues,
509 double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
510 {
511 // reduced flops version: rely on previously computed P_i
512 if (n >= 0) outputValues(0) = 0.0;
513 if (n >= 1) outputValues(1) = 0.0;
514 for (int i=2; i<=n; i++)
515 {
516 const ScalarType & P_i_minus_1 = jacobiValues(i-1); // define as P_i just for clarity of the code below
517 const ScalarType & P_i_minus_2 = jacobiValues(i-2);
518 outputValues(i) = - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
519 }
520 }
521
532 template<typename OutputValueViewType, typename ScalarType, typename ScalarTypeForScaling>
533 KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues_dt(OutputValueViewType outputValues,
534 double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
535 {
536 // memory-conserving version -- place the Jacobi values in the final output container
537 shiftedScaledJacobiValues(outputValues, alpha, n, x, t);
538
539 ScalarType P_i_minus_2 = outputValues(0);
540 ScalarType P_i_minus_1 = outputValues(1);
541
542 if (n >= 0) outputValues(0) = 0.0;
543 if (n >= 1) outputValues(1) = 0.0;
544
545 for (int i=2; i<=n; i++)
546 {
547 const ScalarType L_i_dt = - (i-1.) / (2. * i - 2. + alpha) * (P_i_minus_1 + t * P_i_minus_2);
548
549 P_i_minus_2 = P_i_minus_1;
550 P_i_minus_1 = outputValues(i);
551
552 outputValues(i) = L_i_dt;
553 }
554 }
555 } // namespace Polynomials
556} // namespace Intrepid2
557
558#endif /* Intrepid2_Polynomials_h */
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation.
KOKKOS_INLINE_FUNCTION void shiftedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
Evaluate shifted Legendre polynomials up to order n at a specified point in [0,1].
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues_dx(OutputValueViewType outputValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
x derivative of integrated Jacobi polynomials L_i for i>=1, defined for x in [0,1].
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues(OutputValueViewType outputValues, const OutputValueViewType jacobiValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
Integrated Jacobi values, defined for x in [0,1].
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dt(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
t derivative of shifted, scaled integrated Legendre polynomials L_i for i>=1, defined for x in [0,...
KOKKOS_INLINE_FUNCTION void legendreDerivativeValues(OutputValueViewType outputValues, const OutputValueViewType legendreValues, Intrepid2::ordinal_type n, ScalarType x)
Evaluate first derivatives of Legendre polynomials up to order n at a specified point,...
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
Integrated Legendre polynomials L_i for i>=1, defined for x in [0,1].
KOKKOS_INLINE_FUNCTION void shiftedScaledLegendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
Evaluate shifted, scaled Legendre polynomials up to order n at a specified point in [0,...
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedLegendreValues_dx(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
x derivative of shifted, scaled integrated Legendre polynomials L_i for i>=1, defined for x in [0,...
KOKKOS_INLINE_FUNCTION void shiftedScaledIntegratedJacobiValues_dt(OutputValueViewType outputValues, const OutputValueViewType jacobiValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
t derivative of shifted, scaled integrated Jacobi polynomials L_i for i>=1, defined for x in [0,...
KOKKOS_INLINE_FUNCTION void shiftedScaledJacobiValues(OutputValueViewType outputValues, double alpha, Intrepid2::ordinal_type n, ScalarType x, ScalarTypeForScaling t)
Shifted, scaled Jacobi values, defined for x in [0,1].
KOKKOS_INLINE_FUNCTION void legendreValues(OutputValueViewType outputValues, Intrepid2::ordinal_type n, ScalarType x)
Evaluate Legendre polynomials up to order n at a specified point.
Contains definitions of custom data types in Intrepid2.
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, .