-- Hoogle documentation, generated by Haddock -- See Hoogle, http://www.haskell.org/hoogle/ -- | Polynomials -- -- A type for representing polynomials, several functions for -- manipulating and evaluating them, and several interesting polynomial -- sequences. @package polynomial @version 0.6 module Math.Polynomial data Endianness -- | Big-Endian (head is highest-order term) BE :: Endianness -- | Little-Endian (head is const term) LE :: Endianness data Poly a -- | Make a Poly from a list of coefficients using the specified -- coefficient order. poly :: Num a => Endianness -> [a] -> Poly a -- | Get the coefficients of a a Poly in the specified order. polyCoeffs :: Num a => Endianness -> Poly a -> [a] polyIsZero :: Num a => Poly a -> Bool polyIsOne :: Num a => Poly a -> Bool -- | The polynomial "0" zero :: Num a => Poly a -- | The polynomial "1" one :: Num a => Poly a -- | Given some constant k, construct the polynomial whose value -- is constantly k. constPoly :: Num a => a -> Poly a -- | The polynomial (in x) "x" x :: Num a => Poly a -- | Given some scalar s and a polynomial f, computes the -- polynomial g such that: -- --
-- evalPoly g x = s * evalPoly f x --scalePoly :: Num a => a -> Poly a -> Poly a -- | Given some polynomial f, computes the polynomial g -- such that: -- --
-- evalPoly g x = negate (evalPoly f x) --negatePoly :: Num a => Poly a -> Poly a -- | composePoly f g constructs the polynomial h such -- that: -- --
-- evalPoly h = evalPoly f . evalPoly g ---- -- This is a very expensive operation and, in general, returns a -- polynomial that is quite a bit more expensive to evaluate than -- f and g together (because it is of a much higher -- order than either). Unless your polynomials are quite small or you are -- quite certain you need the coefficients of the composed polynomial, it -- is recommended that you simply evaluate f and g and -- explicitly compose the resulting functions. This will usually be much -- more efficient. composePoly :: Num a => Poly a -> Poly a -> Poly a -- | Given polynomials f and g, computes the polynomial -- h such that: -- --
-- evalPoly h x = evalPoly f x + evalPoly g x --addPoly :: Num a => Poly a -> Poly a -> Poly a sumPolys :: Num a => [Poly a] -> Poly a -- | Given polynomials f and g, computes the polynomial -- h such that: -- --
-- evalPoly h x = evalPoly f x * evalPoly g x --multPoly :: Num a => Poly a -> Poly a -> Poly a -- | Given a polynomial f and exponent n, computes the -- polynomial g such that: -- --
-- evalPoly g x = evalPoly f x ^ n --powPoly :: (Num a, Integral b) => Poly a -> b -> Poly a -- | Given polynomials a and b, with b not -- zero, computes polynomials q and r such that: -- --
-- addPoly (multPoly q b) r == a --quotRemPoly :: Fractional a => Poly a -> Poly a -> (Poly a, Poly a) quotPoly :: Fractional a => Poly a -> Poly a -> Poly a remPoly :: Fractional a => Poly a -> Poly a -> Poly a -- | Evaluate a polynomial at a point or, equivalently, convert a -- polynomial to the function it represents. For example, evalPoly -- x = id and evalPoly (constPoly k) = -- const k. evalPoly :: Num a => Poly a -> a -> a -- | Evaluate a polynomial and its derivative (respectively) at a point. evalPolyDeriv :: Num a => Poly a -> a -> (a, a) -- | Evaluate a polynomial and all of its nonzero derivatives at a point. -- This is roughly equivalent to: -- --
-- evalPolyDerivs p x = map (`evalPoly` x) (takeWhile (not . polyIsZero) (iterate polyDeriv p)) --evalPolyDerivs :: Num a => Poly a -> a -> [a] -- | "Contract" a polynomial by attempting to divide out a root. -- -- contractPoly p a returns (q,r) such that q*(x-a) -- + r == p contractPoly :: Num a => Poly a -> a -> (Poly a, a) -- | gcdPoly a b computes the highest order monic polynomial that -- is a divisor of both a and b. If both a and -- b are zero, the result is undefined. gcdPoly :: Fractional a => Poly a -> Poly a -> Poly a -- | Separate a nonzero polynomial into a set of factors none of which have -- multiple roots, and the product of which is the original polynomial. -- Note that if division is not exact, it may fail to separate roots. -- Rational coefficients is a good idea. -- -- Useful when applicable as a way to simplify root-finding problems. separateRoots :: Fractional a => Poly a -> [Poly a] -- | Compute the derivative of a polynomial. polyDeriv :: Num a => Poly a -> Poly a -- | Compute the definite integral (from 0 to x) of a polynomial. polyIntegral :: Fractional a => Poly a -> Poly a module Math.Polynomial.Bernstein -- | The Bernstein basis polynomials. The nth inner list is a -- basis for the polynomials of order n or lower. The -- nth basis consists of n polynomials of order -- n which sum to 1, and have roots of varying -- multiplicities at 0 and 1. bernstein :: [[Poly Integer]] -- | evalBernstein n v x evaluates the v'th Bernstein -- polynomial of order n at the point x. evalBernstein :: (Integral a, Num b) => a -> a -> b -> b -- | bernsteinFit n f: Approximate a function f as a -- linear combination of Bernstein polynomials of order n. This -- approximation converges slowly but uniformly to f on the -- interval [0,1]. bernsteinFit :: (Fractional b, Integral a) => a -> (b -> b) -> [b] -- | Evaluate a polynomial given as a list of n coefficients for -- the nth Bernstein basis. Roughly: -- --
-- evalBernsteinSeries cs = sum (zipWith scalePoly cs (bernstein !! (length cs - 1))) --evalBernsteinSeries :: Num a => [a] -> a -> a -- | de Casteljau's algorithm, returning the whole tableau. Used both for -- evaluating and splitting polynomials in Bernstein form. deCasteljau :: Num a => [a] -> a -> [[a]] -- | Given a polynomial in Bernstein form (that is, a list of coefficients -- for a basis set from bernstein, such as is returned by -- bernsteinFit) and a parameter value x, split the -- polynomial into two halves, mapping [0,x] and [x,1] -- respectively onto [0,1]. -- -- A typical use for this operation would be to split a Bezier curve -- (inserting a new knot at x). splitBernsteinSeries :: Num a => [a] -> a -> ([a], [a]) module Math.Polynomial.Chebyshev -- | The Chebyshev polynomials of the first kind with Integer -- coefficients. ts :: [Poly Integer] us :: [Poly Integer] -- | Compute the coefficients of the n'th Chebyshev polynomial of the first -- kind. t :: Num a => Int -> Poly a -- | Compute the coefficients of the n'th Chebyshev polynomial of the -- second kind. u :: Num a => Int -> Poly a -- | Evaluate the n'th Chebyshev polynomial of the first kind at a point X. -- Both more efficient and more numerically stable than computing the -- coefficients and evaluating the polynomial. evalT :: Num a => Int -> a -> a -- | Evaluate all the Chebyshev polynomials of the first kind at a point X. evalTs :: Num a => a -> [a] -- | Evaluate the n'th Chebyshev polynomial of the second kind at a point -- X. Both more efficient and more numerically stable than computing the -- coefficients and evaluating the polynomial. evalU :: Num a => Int -> a -> a -- | Evaluate all the Chebyshev polynomials of the second kind at a point -- X. evalUs :: Num a => a -> [a] -- | Evaluate the n'th Chebyshev polynomials of both kinds at a point X. evalTU :: Num a => Int -> a -> (a, a) -- | Evaluate all the Chebyshev polynomials of both kinds at a point X. evalTsUs :: Num a => a -> ([a], [a]) -- | Compute the roots of the n'th Chebyshev polynomial of the first kind. tRoots :: Floating a => Int -> [a] -- | Compute the extreme points of the n'th Chebyshev polynomial of the -- first kind. tExtrema :: Floating a => Int -> [a] -- | chebyshevFit n f returns a list of N coefficients cs -- such that f x ~= sum (zipWith (*) cs (evalTs x)) on -- the interval -1 < x < 1. -- -- The N roots of the N'th Chebyshev polynomial are the fitting points at -- which the function will be evaluated and at which the approximation -- will be exact. These points always lie within the interval -1 < x -- < 1. Outside this interval, the approximation will diverge quickly. -- -- This function deviates from most chebyshev-fit implementations in that -- it returns the first coefficient pre-scaled so that the series -- evaluation operation is a simple inner product, since in most other -- algorithms operating on chebyshev series, that factor is almost always -- a nuissance. chebyshevFit :: Floating a => Int -> (a -> a) -> [a] -- | Evaluate a Chebyshev series expansion with a finite number of terms. -- -- Note that this function expects the first coefficient to be pre-scaled -- by 1/2, which is what is produced by chebyshevFit. Thus, this -- computes a simple inner product of the given list with a -- matching-length sequence of chebyshev polynomials. evalChebyshevSeries :: Num a => [a] -> a -> a module Math.Polynomial.Lagrange -- | Returns the Lagrange basis set of polynomials associated with a set of -- points. This is the set of polynomials each of which is 1 at -- its corresponding point in the input list and 0 at all -- others. -- -- These polynomials are especially convenient, mathematically, for -- interpolation. The interpolating polynomial for a set of points -- (x,y) is given by using the ys as coefficients for -- the basis given by lagrangeBasis xs. Computationally, this is -- not an especially stable procedure though. -- Math.Polynomial.Interpolation.lagrangePolyFit implements a -- slightly better algorithm based on the same idea. -- -- Generally it is better to not compute the coefficients at all. -- Math.Polynomial.Interpolation.polyInterp evaluates the -- interpolating polynomial directly, and is both quicker and more stable -- than any method I know of that computes the coefficients. lagrangeBasis :: Fractional a => [a] -> [Poly a] -- | Construct the Lagrange master polynomial for the Lagrange -- barycentric form: That is, the monic polynomial with a root at each -- point in the input list. lagrange :: Num a => [a] -> Poly a -- | Compute the weights associated with each abscissa in the Lagrange -- barycentric form. lagrangeWeights :: Fractional a => [a] -> [a] module Math.Polynomial.Interpolation -- | Evaluate a polynomial passing through the specified set of points. The -- order of the interpolating polynomial will (at most) be one less than -- the number of points given. polyInterp :: Fractional a => [(a, a)] -> a -> a -- | Computes the tableau generated by Neville's algorithm. Each successive -- row of the table is a list of interpolants one order higher than the -- previous, using a range of input points starting at the same position -- in the input list as the interpolant's position in the output list. neville :: Fractional a => [(a, a)] -> a -> [[a]] -- | Computes the tableau generated by a modified form of Neville's -- algorithm described in Numerical Recipes, Ch. 3, Sec. 2, which records -- the differences between interpolants at each level. Each pair (c,d) is -- the amount to add to the previous level's interpolant at either the -- same or the subsequent position (respectively) in order to obtain the -- new level's interpolant. Mathematically, either sum yields the same -- value, but due to numerical errors they may differ slightly, and some -- "paths" through the table may yield more accurate final results than -- others. nevilleDiffs :: Fractional a => [(a, a)] -> a -> [[(a, a)]] -- | Fit a polynomial to a set of points by iteratively evaluating the -- interpolated polynomial (using polyInterp) at 0 to establish -- the constant coefficient and reducing the polynomial by subtracting -- that coefficient from all y's and dividing by their corresponding x's. -- -- Slower than lagrangePolyFit but stable under different sets of -- conditions. -- -- Note that computing the coefficients of a fitting polynomial is an -- inherently ill-conditioned problem. In most cases it is both faster -- and more accurate to use polyInterp or nevilleDiffs -- instead of evaluating a fitted polynomial. iterativePolyFit :: Fractional a => [(a, a)] -> Poly a -- | Fit a polynomial to a set of points using barycentric Lagrange -- polynomials. -- -- Note that computing the coefficients of a fitting polynomial is an -- inherently ill-conditioned problem. In most cases it is both faster -- and more accurate to use polyInterp or nevilleDiffs -- instead of evaluating a fitted polynomial. lagrangePolyFit :: Fractional a => [(a, a)] -> Poly a module Math.Polynomial.Legendre -- | The Legendre polynomials with Rational coefficients. These -- polynomials form an orthogonal basis of the space of all polynomials, -- relative to the L2 inner product on [-1,1] (which is given by -- integrating the product of 2 polynomials over that range). legendres :: [Poly Rational] -- | Compute the coefficients of the n'th Legendre polynomial. legendre :: Fractional a => Int -> Poly a -- | Evaluate the n'th Legendre polynomial at a point X. Both more -- efficient and more numerically stable than computing the coefficients -- and evaluating the polynomial. evalLegendre :: Fractional a => Int -> a -> a -- | Evaluate all the Legendre polynomials at a point X. evalLegendres :: Fractional a => a -> [a] -- | Evaluate the n'th Legendre polynomial and its derivative at a point X. -- Both more efficient and more numerically stable than computing the -- coefficients and evaluating the polynomial. evalLegendreDeriv :: Fractional a => Int -> a -> (a, a) -- | Zeroes of the n'th Legendre polynomial. legendreRoots :: (Fractional b, Ord b) => Int -> b -> [b] module Math.Polynomial.Newton -- | Returns the Newton basis set of polynomials associated with a set of -- abscissas. This is the set of monic polynomials each of which is -- 0 at all previous points in the input list. newtonBasis :: Num a => [a] -> [Poly a] -- | This module exports a Num instance for the Poly type. -- This instance does not implement all operations, because abs -- and signum are simply not definable, so I have placed it into a -- separate module so that I can make people read this caveat ;). -- -- Use at your own risk. module Math.Polynomial.NumInstance instance Num a => Num (Poly a)