]> gitweb.michael.orlitzky.com - sage.d.git/blobdiff - mjo/interpolation.py
COPYING,LICENSE: add (AGPL-3.0+)
[sage.d.git] / mjo / interpolation.py
index 0cafc283e7b445f0ec8af9487f3964195a291cc5..6e0ecb3d79c5129a2e6270177b648142085786d1 100644 (file)
@@ -1,5 +1,25 @@
 from sage.all import *
-from misc import product
+product = prod
+
+
+def lagrange_denominator(k, xs):
+    """
+    Return the denominator of the kth Lagrange coefficient.
+
+    INPUT:
+
+      - ``k`` -- The index of the coefficient.
+
+      - ``xs`` -- The list of points at which the function values are
+        known.
+
+    OUTPUT:
+
+    The product of all xs[j] with j != k.
+
+    """
+    return product( xs[k] - xs[j] for j in range(len(xs)) if j != k )
+
 
 def lagrange_coefficient(k, x, xs):
     """
@@ -12,9 +32,9 @@ def lagrange_coefficient(k, x, xs):
 
     INPUT:
 
-      - ``k`` -- the index of the coefficient.
+      - ``k`` -- The index of the coefficient.
 
-      - ``x`` -- the symbolic variable to use for the first argument
+      - ``x`` -- The symbolic variable to use for the first argument
         of l_{k}.
 
       - ``xs`` -- The list of points at which the function values are
@@ -22,17 +42,21 @@ def lagrange_coefficient(k, x, xs):
 
     OUTPUT:
 
-    A symbolic function of one variable.
+    A symbolic expression of one variable.
+
+    SETUP::
+
+        sage: from mjo.interpolation import lagrange_coefficient
 
     TESTS::
 
         sage: xs = [ -pi/2, -pi/6, 0, pi/6, pi/2 ]
         sage: lagrange_coefficient(0, x, xs)
-        1/8*(pi - 6*x)*(pi - 2*x)*(pi + 6*x)*x/pi^4
+        1/8*(pi + 6*x)*(pi - 2*x)*(pi - 6*x)*x/pi^4
 
     """
-    numerator = product([x - xs[j] for j in range(0, len(xs)) if j != k])
-    denominator = product([xs[k] - xs[j] for j in range(0, len(xs)) if j != k])
+    numerator = lagrange_psi(x, xs)/(x - xs[k])
+    denominator = lagrange_denominator(k, xs)
 
     return (numerator / denominator)
 
@@ -40,7 +64,7 @@ def lagrange_coefficient(k, x, xs):
 
 def lagrange_polynomial(x, xs, ys):
     """
-    Return the Lagrange form of the interpolation polynomial in `x` of
+    Return the Lagrange form of the interpolating polynomial in `x`
     at the points (xs[k], ys[k]).
 
     INPUT:
@@ -53,12 +77,16 @@ def lagrange_polynomial(x, xs, ys):
 
     OUTPUT:
 
-    A symbolic function (polynomial) interpolating each (xs[k], ys[k]).
+    A symbolic expression (polynomial) interpolating each (xs[k], ys[k]).
+
+    SETUP::
+
+        sage: from mjo.interpolation import lagrange_polynomial
 
     TESTS::
 
         sage: xs = [ -pi/2, -pi/6, 0, pi/6, pi/2 ]
-        sage: ys = map(sin, xs)
+        sage: ys = list(map(sin, xs))
         sage: L = lagrange_polynomial(x, xs, ys)
         sage: expected  = 27/16*(pi - 6*x)*(pi - 2*x)*(pi + 2*x)*x/pi^4
         sage: expected -= 1/8*(pi - 6*x)*(pi - 2*x)*(pi + 6*x)*x/pi^4
@@ -68,97 +96,273 @@ def lagrange_polynomial(x, xs, ys):
         True
 
     """
-    ls = [ lagrange_coefficient(k, x, xs) for k in range(0, len(xs)) ]
-    sigma = sum([ ys[k] * ls[k] for k in range(0, len(xs)) ])
-    return sigma
+    ls = [ lagrange_coefficient(k, x, xs) for k in range(len(xs)) ]
+    return sum( ys[k] * ls[k] for k in range(len(xs)) )
+
 
 
-def divided_difference(f, xs):
+def lagrange_interpolate(f, x, xs):
     """
-    Return the Newton divided difference of `f` at the points
-    `xs`. Reference:
+    Interpolate the function ``f`` at the points ``xs`` using the
+    Lagrange form of the interpolating polynomial.
+
+    INPUT:
+
+      - ``f`` -- The function to interpolate.
+
+      - ``x`` -- The independent variable of the resulting polynomial.
+
+      - ``xs`` -- A list of points at which to interpolate ``f``.
+
+    OUTPUT:
+
+    A polynomial in ``x`` which interpolates ``f`` at ``xs``.
+
+    SETUP::
+
+        sage: from mjo.interpolation import lagrange_interpolate
+
+    EXAMPLES:
+
+    We're exact on polynomials of degree `n` if we use `n+1` points::
+
+        sage: t = SR.symbol('t', domain='real')
+        sage: lagrange_interpolate(x^2, t, [-1,0,1]).simplify_rational()
+        t^2
+
+    """
+    # f should be a function of one variable.
+    z = f.variables()[0]
+    # We're really just doing map(f, xs) here; the additional
+    # gymnastics are to avoid a warning when calling `f` with an
+    # unnamed argument.
+    ys = [ f({z: xk}) for xk in xs ]
+    return lagrange_polynomial(x, xs, ys)
+
+
+
+def divided_difference_coefficients(xs):
+    """
+    Assuming some function `f`, compute the coefficients of the
+    divided difference f[xs[0], ..., xs[n]].
+
+    SETUP::
+
+        sage: from mjo.interpolation import divided_difference_coefficients
+
+    TESTS::
+
+        sage: divided_difference_coefficients([0])
+        [1]
+        sage: divided_difference_coefficients([0, pi])
+        [-1/pi, 1/pi]
+        sage: divided_difference_coefficients([0, pi, 2*pi])
+        [1/2/pi^2, -1/pi^2, 1/2/pi^2]
+
+    """
+    return [ ~lagrange_denominator(k, xs) for k in range(len(xs)) ]
+
+
+def divided_difference(xs, ys):
+    """
+    Return the Newton divided difference of the points (xs[k],
+    ys[k]). Reference:
 
       http://en.wikipedia.org/wiki/Divided_differences
 
     INPUT:
 
-      - ``f`` -- The function whose divided difference we seek.
+      - ``xs`` -- The list of x-values.
 
-      - ``xs`` -- The list of points at which to compute `f`.
+      - ``ys`` -- The function values at `xs`.
 
     OUTPUT:
 
-    The divided difference of `f` at ``xs``.
+    The (possibly symbolic) divided difference function.
+
+    SETUP::
+
+        sage: from mjo.interpolation import divided_difference
 
     TESTS::
 
-        sage: divided_difference(sin, [0])
+        sage: xs = [0]
+        sage: ys = list(map(sin, xs))
+        sage: divided_difference(xs, ys)
         0
-        sage: divided_difference(sin, [0, pi])
+        sage: xs = [0, pi]
+        sage: ys = list(map(sin, xs))
+        sage: divided_difference(xs, ys)
         0
-        sage: divided_difference(sin, [0, pi, 2*pi])
+        sage: xs = [0, pi, 2*pi]
+        sage: ys = list(map(sin, xs))
+        sage: divided_difference(xs, ys)
         0
 
     We try something entirely symbolic::
 
-        sage: f = function('f'x)
-        sage: divided_difference(f, [x])
+        sage: f = function('f')(x)
+        sage: divided_difference([x], [f(x=x)])
         f(x)
-        sage: x1,x2 = var('x1,x2')
-        sage: divided_difference(f, [x1,x2])
-        (f(x1) - f(x2))/(x1 - x2)
+        sage: x1,x2 = SR.var('x1,x2')
+        sage: divided_difference([x1,x2], [f(x=x1),f(x=x2)])
+        f(x1)/(x1 - x2) - f(x2)/(x1 - x2)
 
     """
-    if (len(xs) == 1):
-        # Avoid that goddamned DeprecationWarning when we have a named
-        # argument but don't know what it is.
-        if len(f.variables()) == 0:
-            return f(xs[0])
-        else:
-            v = f.variables()[0]
-            return f({ v: xs[0] })
+    coeffs = divided_difference_coefficients(xs)
+    v_cs = vector(coeffs)
+    v_ys = vector(ys)
+    return v_cs.dot_product(v_ys)
 
-    # Use the recursive definition.
-    numerator  = divided_difference(f, xs[1:])
-    numerator -= divided_difference(f, xs[:-1])
-    return numerator / (xs[-1] - xs[0])
 
-
-def newton_polynomial(f, x, xs):
+def newton_polynomial(x, xs, ys):
     """
-    Return the Newton form of the interpolating polynomial of `f` at
-    the points `xs` in the variable `x`.
+    Return the Newton form of the interpolating polynomial of the
+    points (xs[k], ys[k]) in the variable `x`.
 
     INPUT:
 
-      - ``f`` -- The function to interpolate.
-
       - ``x`` -- The independent variable to use for the interpolating
         polynomial.
 
-      - ``xs`` -- The list of points at which to interpolate `f`.
+      - ``xs`` -- The list of x-values.
+
+      - ``ys`` -- The function values at `xs`.
 
     OUTPUT:
 
-    A symbolic function.
+    A symbolic expression.
 
-    TESTS:
+    SETUP::
+
+        sage: from mjo.interpolation import lagrange_polynomial, newton_polynomial
+
+    TESTS::
 
         sage: xs = [ -pi/2, -pi/6, 0, pi/6, pi/2 ]
-        sage: ys = map(sin, xs)
+        sage: ys = list(map(sin, xs))
         sage: L = lagrange_polynomial(x, xs, ys)
-        sage: N = newton_polynomial(sin, x, xs)
+        sage: N = newton_polynomial(x, xs, ys)
         sage: bool(N == L)
         True
 
     """
-    degree = len(xs) - 1
+    return sum( divided_difference(xs[:k+1], ys[:k+1])*lagrange_psi(x, xs[:k])
+                for k in range(len(xs)) )
+
+
+def hermite_coefficient(k, x, xs):
+    """
+    Return the Hermite coefficient h_{k}(x). See Atkinson, p. 160.
+
+    INPUT:
+
+      - ``k`` -- The index of the coefficient.
+
+      - ``x`` -- The symbolic variable to use as the argument of h_{k}.
+
+      - ``xs`` -- The list of points at which the function values are
+        known.
+
+    OUTPUT:
+
+    A symbolic expression.
+
+    """
+    lk = lagrange_coefficient(k, x, xs)
+    return (1 - 2*lk.diff(x)(x=xs[k])*(x - xs[k]))*(lk**2)
+
+
+def hermite_deriv_coefficient(k, x, xs):
+    """
+    Return the Hermite derivative coefficient, \tilde{h}_{k}(x). See
+    Atkinson, p. 160.
+
+    INPUT:
+
+      - ``k`` -- The index of the coefficient.
+
+      - ``x`` -- The symbolic variable to use as the argument of h_{k}.
+
+      - ``xs`` -- The list of points at which the function values are
+        known.
+
+    OUTPUT:
+
+    A symbolic expression.
+
+    """
+    lk = lagrange_coefficient(k, x, xs)
+    return (x - xs[k])*(lk**2)
+
+
+def hermite_interpolant(x, xs, ys, y_primes):
+    """
+    Return the Hermite interpolant `H(x)` such that H(xs[k]) = ys[k]
+    and H'(xs[k]) = y_primes[k] for each k.
+
+    Reference: Atkinson, p. 160.
+
+    INPUT:
+
+      - ``x`` -- The symbolic variable to use as the argument of H(x).
+
+      - ``xs`` -- The list of points at which the function values are
+        known.
+
+      - ``ys`` -- The function values at the `xs`.
+
+      - ``y_primes`` -- The derivatives at the `xs`.
+
+    OUTPUT:
+
+    A symbolic expression.
+
+    SETUP::
+
+        sage: from mjo.interpolation import hermite_interpolant
+
+    TESTS::
+
+        sage: xs = [ 0, pi/6, pi/2 ]
+        sage: ys = list(map(sin, xs))
+        sage: y_primes = list(map(cos, xs))
+        sage: H = hermite_interpolant(x, xs, ys, y_primes)
+        sage: expected  = -27/4*(pi - 6*x)*(pi - 2*x)^2*sqrt(3)*x^2/pi^4
+        sage: expected += (5*(pi - 2*x)/pi + 1)*(pi - 6*x)^2*x^2/pi^4
+        sage: expected += 81/2*((pi - 6*x)/pi + 1)*(pi - 2*x)^2*x^2/pi^4
+        sage: expected += (pi - 6*x)^2*(pi - 2*x)^2*x/pi^4
+        sage: bool(H == expected)
+        True
+
+    """
+    s1 = sum( ys[k] * hermite_coefficient(k, x, xs)
+               for k in range(len(xs)) )
+
+    s2 = sum( y_primes[k] * hermite_deriv_coefficient(k, x, xs)
+               for k in range(len(xs)) )
+
+    return (s1 + s2)
+
+
+def lagrange_psi(x, xs):
+    """
+    The function,
 
-    N = SR(0)
+      Psi(x) = (x - xs[0])*(x - xs[1])* ... *(x - xs[-1])
 
-    for k in range(0, degree+1):
-        term  = divided_difference(f, xs[:k+1])
-        term *= product([ x - xk for xk in xs[:k]])
-        N += term
+    used in Lagrange and Hermite interpolation.
+
+    INPUT:
+
+      - ``x`` -- The independent variable of the resulting expression.
+
+      - ``xs`` -- A list of points.
+
+    OUTPUT:
+
+    A symbolic expression in one variable, `x`.
+
+    """
 
-    return N
+    return product( (x - xj) for xj in xs )