-- Hoogle documentation, generated by Haddock
-- See Hoogle, http://www.haskell.org/hoogle/


-- | Numerical computation
--   
--   Purely functional interface to selected numerical computations,
--   internally implemented using GSL.
@package hmatrix-gsl
@version 0.19.0.1


-- | Numerical integration routines.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Numerical-Integration.html#Numerical-Integration</a>
module Numeric.GSL.Integration

-- | Numerical integration using <i>gsl_integration_qng</i> (useful for
--   fast integration of smooth functions). For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateQNG 1E-6
--   
--   &gt;&gt;&gt; quad (\x -&gt; 4/(1+x*x)) 0 1
--   (3.141592653589793,3.487868498008632e-14)
--   </pre>
integrateQNG :: Double -> (Double -> Double) -> Double -> Double -> (Double, Double)

-- | Numerical integration using <i>gsl_integration_qags</i> (adaptive
--   integration with singularities). For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateQAGS 1E-9 1000
--   
--   &gt;&gt;&gt; let f a x = x**(-0.5) * log (a*x)
--   
--   &gt;&gt;&gt; quad (f 1) 0 1
--   (-3.999999999999974,4.871658632055187e-13)
--   </pre>
integrateQAGS :: Double -> Int -> (Double -> Double) -> Double -> Double -> (Double, Double)

-- | Numerical integration using <i>gsl_integration_qagi</i> (integration
--   over the infinite integral -Inf..Inf using QAGS). For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateQAGI 1E-9 1000
--   
--   &gt;&gt;&gt; let f a x = exp(-a * x^2)
--   
--   &gt;&gt;&gt; quad (f 0.5)
--   (2.5066282746310002,6.229215880648858e-11)
--   </pre>
integrateQAGI :: Double -> Int -> (Double -> Double) -> (Double, Double)

-- | Numerical integration using <i>gsl_integration_qagiu</i> (integration
--   over the semi-infinite integral a..Inf). For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateQAGIU 1E-9 1000
--   
--   &gt;&gt;&gt; let f a x = exp(-a * x^2)
--   
--   &gt;&gt;&gt; quad (f 0.5) 0
--   (1.2533141373155001,3.114607940324429e-11)
--   </pre>
integrateQAGIU :: Double -> Int -> (Double -> Double) -> Double -> (Double, Double)

-- | Numerical integration using <i>gsl_integration_qagil</i> (integration
--   over the semi-infinite integral -Inf..b). For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateQAGIL 1E-9 1000
--   
--   &gt;&gt;&gt; let f a x = exp(-a * x^2)
--   
--   &gt;&gt;&gt; quad (f 0.5) 0
--   (1.2533141373155001,3.114607940324429e-11)
--   </pre>
integrateQAGIL :: Double -> Int -> (Double -> Double) -> Double -> (Double, Double)

-- | Numerical integration using <i>gsl_integration_cquad</i> (quadrature
--   for general integrands). From the GSL manual:
--   
--   <pre>
--   CQUAD is a new doubly-adaptive general-purpose quadrature routine
--   which can handle most types of singularities, non-numerical function
--   values such as Inf or NaN, as well as some divergent integrals. It
--   generally requires more function evaluations than the integration
--   routines in QUADPACK, yet fails less often for difficult integrands.
--   </pre>
--   
--   For example:
--   
--   <pre>
--   &gt;&gt;&gt; let quad = integrateCQUAD 1E-12 1000
--   
--   &gt;&gt;&gt; let f a x = exp(-a * x^2)
--   
--   &gt;&gt;&gt; quad (f 0.5) 2 5
--   (5.7025405463957006e-2,9.678874441303705e-16,95)
--   </pre>
--   
--   Unlike other quadrature methods, integrateCQUAD also returns the
--   number of function evaluations required.
integrateCQUAD :: Double -> Int -> (Double -> Double) -> Double -> Double -> (Double, Double, Int)


-- | Fourier Transform.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Fast-Fourier-Transforms.html#Fast-Fourier-Transforms</a>
module Numeric.GSL.Fourier

-- | Fast 1D Fourier transform of a <a>Vector</a> <tt>(</tt><a>Complex</a>
--   <a>Double</a><tt>)</tt> using <i>gsl_fft_complex_forward</i>. It uses
--   the same scaling conventions as GNU Octave.
--   
--   <pre>
--   &gt;&gt;&gt; fft (fromList [1,2,3,4])
--   fromList [10.0 :+ 0.0,(-2.0) :+ 2.0,(-2.0) :+ 0.0,(-2.0) :+ (-2.0)]
--   </pre>
fft :: Vector (Complex Double) -> Vector (Complex Double)

-- | The inverse of <a>fft</a>, using <i>gsl_fft_complex_inverse</i>.
ifft :: Vector (Complex Double) -> Vector (Complex Double)


-- | Nonlinear Least-Squares Fitting
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Nonlinear-Least_002dSquares-Fitting.html</a>
--   
--   The example program in the GSL manual (see examples/fitting.hs):
--   
--   <pre>
--   dat = [
--    ([0.0],([6.0133918608118675],0.1)),
--    ([1.0],([5.5153769909966535],0.1)),
--    ([2.0],([5.261094606015287],0.1)),
--    ...
--    ([39.0],([1.0619821710802808],0.1))]
--   
--   expModel [a,lambda,b] [t] = [a * exp (-lambda * t) + b]
--   
--   expModelDer [a,lambda,b] [t] = [[exp (-lambda * t), -t * a * exp(-lambda*t) , 1]]
--   
--   (sol,path) = fitModelScaled 1E-4 1E-4 20 (expModel, expModelDer) dat [1,0,0]
--   </pre>
--   
--   <pre>
--   &gt;&gt;&gt; path
--   (6&gt;&lt;5)
--    [ 1.0,  76.45780563978782, 1.6465931240727802, 1.8147715267618197e-2, 0.6465931240727797
--    , 2.0, 37.683816318260355,  2.858760367632973,  8.092094813253975e-2, 1.4479636296208662
--    , 3.0,    9.5807893736187,  4.948995119561291,   0.11942927999921617, 1.0945766509238248
--    , 4.0,  5.630494933603935,  5.021755718065913,   0.10287787128056883, 1.0338835440862608
--    , 5.0,  5.443976278682909,  5.045204331329302,   0.10405523433131504,  1.019416067207375
--    , 6.0, 5.4439736648994685,  5.045357818922331,   0.10404905846029407, 1.0192487112786812 ]
--   
--   &gt;&gt;&gt; sol
--   [(5.045357818922331,6.027976702418132e-2),
--   (0.10404905846029407,3.157045047172834e-3),
--   (1.0192487112786812,3.782067731353722e-2)]
--   </pre>
module Numeric.GSL.Fitting

-- | Nonlinear multidimensional least-squares fitting.
nlFitting :: FittingMethod -> Double -> Double -> Int -> (Vector Double -> Vector Double) -> (Vector Double -> Matrix Double) -> Vector Double -> (Vector Double, Matrix Double)
data FittingMethod

-- | Interface to gsl_multifit_fdfsolver_lmsder. This is a robust and
--   efficient version of the Levenberg-Marquardt algorithm as implemented
--   in the scaled lmder routine in minpack. Minpack was written by Jorge
--   J. More, Burton S. Garbow and Kenneth E. Hillstrom.
LevenbergMarquardtScaled :: FittingMethod

-- | This is an unscaled version of the lmder algorithm. The elements of
--   the diagonal scaling matrix D are set to 1. This algorithm may be
--   useful in circumstances where the scaled version of lmder converges
--   too slowly, or the function is already scaled appropriately.
LevenbergMarquardt :: FittingMethod

-- | Higher level interface to <a>nlFitting</a>
--   <a>LevenbergMarquardtScaled</a>. The optimization function and
--   Jacobian are automatically built from a model f vs x = y and its
--   derivatives, and a list of instances (x, (y,sigma)) to be fitted.
fitModelScaled :: Double -> Double -> Int -> ([Double] -> x -> [Double], [Double] -> x -> [[Double]]) -> [(x, ([Double], Double))] -> [Double] -> ([(Double, Double)], Matrix Double)

-- | Higher level interface to <a>nlFitting</a> <a>LevenbergMarquardt</a>.
--   The optimization function and Jacobian are automatically built from a
--   model f vs x = y and its derivatives, and a list of instances (x,y) to
--   be fitted.
fitModel :: Double -> Double -> Int -> ([Double] -> x -> [Double], [Double] -> x -> [[Double]]) -> [(x, [Double])] -> [Double] -> ([Double], Matrix Double)
instance GHC.Internal.Enum.Bounded Numeric.GSL.Fitting.FittingMethod
instance GHC.Internal.Enum.Enum Numeric.GSL.Fitting.FittingMethod
instance GHC.Classes.Eq Numeric.GSL.Fitting.FittingMethod
instance GHC.Internal.Show.Show Numeric.GSL.Fitting.FittingMethod


-- | Numerical differentiation.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Numerical-Differentiation.html#Numerical-Differentiation</a>
--   
--   From the GSL manual: "The functions described in this chapter compute
--   numerical derivatives by finite differencing. An adaptive algorithm is
--   used to find the best choice of finite difference and to estimate the
--   error in the derivative."
module Numeric.GSL.Differentiation

-- | Adaptive central difference algorithm, <i>gsl_deriv_central</i>. For
--   example:
--   
--   <pre>
--   &gt;&gt;&gt; let deriv = derivCentral 0.01
--   
--   &gt;&gt;&gt; deriv sin (pi/4)
--   (0.7071067812000676,1.0600063101654055e-10)
--   
--   &gt;&gt;&gt; cos (pi/4)
--   0.7071067811865476
--   </pre>
derivCentral :: Double -> (Double -> Double) -> Double -> (Double, Double)

-- | Adaptive forward difference algorithm, <i>gsl_deriv_forward</i>. The
--   function is evaluated only at points greater than x, and never at x
--   itself. The derivative is returned in result and an estimate of its
--   absolute error is returned in abserr. This function should be used if
--   f(x) has a discontinuity at x, or is undefined for values less than x.
--   A backward derivative can be obtained using a negative step.
derivForward :: Double -> (Double -> Double) -> Double -> (Double, Double)

-- | Adaptive backward difference algorithm, <i>gsl_deriv_backward</i>.
derivBackward :: Double -> (Double -> Double) -> Double -> (Double, Double)


-- | Interpolation routines.
--   
--   
--   <a>https://www.gnu.org/software/gsl/manual/html_node/Interpolation.html#Interpolation</a>
--   
--   The GSL routines <tt>gsl_spline_eval</tt> and friends are used, but in
--   spite of the names, they are not restricted to spline interpolation.
--   The functions in this module will work for any
--   <a>InterpolationMethod</a>.
module Numeric.GSL.Interpolation
data InterpolationMethod
Linear :: InterpolationMethod
Polynomial :: InterpolationMethod
CSpline :: InterpolationMethod
CSplinePeriodic :: InterpolationMethod
Akima :: InterpolationMethod
AkimaPeriodic :: InterpolationMethod

-- | Evaluate a function by interpolating within the given dataset. For
--   example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = [1..10]
--   
--   &gt;&gt;&gt; let ys map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluate Akima (zip xs ys) 2.2
--   4.840000000000001
--   </pre>
--   
--   To successfully <tt>evaluate points x</tt>, the domain (<tt>x</tt>)
--   values in <tt>points</tt> must be monotonically increasing. The
--   evaluation point <tt>x</tt> must lie between the smallest and largest
--   values in the sampled domain.
evaluate :: InterpolationMethod -> [(Double, Double)] -> Double -> Double

-- | Evaluate a function by interpolating within the given dataset. For
--   example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = vector [1..10]
--   
--   &gt;&gt;&gt; let ys = vector $ map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateV CSpline xs ys 2.2
--   4.818867924528303
--   </pre>
--   
--   To successfully <tt>evaluateV xs ys x</tt>, the vectors of
--   corresponding domain-range values <tt>xs</tt> and <tt>ys</tt> must
--   have identical lengths, and <tt>xs</tt> must be monotonically
--   increasing. The evaluation point <tt>x</tt> must lie between the
--   smallest and largest values in <tt>xs</tt>.
evaluateV :: InterpolationMethod -> Vector Double -> Vector Double -> Double -> Double

-- | Evaluate the derivative of a function by interpolating within the
--   given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = [1..10]
--   
--   &gt;&gt;&gt; let ys map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateDerivative Akima (zip xs ys) 2.2
--   4.4
--   </pre>
--   
--   To successfully <tt>evaluateDerivative points x</tt>, the domain
--   (<tt>x</tt>) values in <tt>points</tt> must be monotonically
--   increasing. The evaluation point <tt>x</tt> must lie between the
--   smallest and largest values in the sampled domain.
evaluateDerivative :: InterpolationMethod -> [(Double, Double)] -> Double -> Double

-- | Evaluate the second derivative of a function by interpolating within
--   the given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = [1..10]
--   
--   &gt;&gt;&gt; let ys map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateDerivative2 Akima (zip xs ys) 2.2
--   2.0
--   </pre>
--   
--   To successfully <tt>evaluateDerivative2 points x</tt>, the domain
--   (<tt>x</tt>) values in <tt>points</tt> must be monotonically
--   increasing. The evaluation point <tt>x</tt> must lie between the
--   smallest and largest values in the sampled domain.
evaluateDerivative2 :: InterpolationMethod -> [(Double, Double)] -> Double -> Double

-- | Evaluate the derivative of a function by interpolating within the
--   given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = vector [1..10]
--   
--   &gt;&gt;&gt; let ys = vector $ map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateDerivativeV CSpline xs ys 2.2
--   4.338867924528302
--   </pre>
--   
--   To successfully <tt>evaluateDerivativeV xs ys x</tt>, the vectors of
--   corresponding domain-range values <tt>xs</tt> and <tt>ys</tt> must
--   have identical lengths, and <tt>xs</tt> must be monotonically
--   increasing. The interpolation point <tt>x</tt> must lie between the
--   smallest and largest values in <tt>xs</tt>.
evaluateDerivativeV :: InterpolationMethod -> Vector Double -> Vector Double -> Double -> Double

-- | Evaluate the second derivative of a function by interpolating within
--   the given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = vector [1..10]
--   
--   &gt;&gt;&gt; let ys = vector $ map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateDerivative2V CSpline xs ys 2.2
--   2.4
--   </pre>
--   
--   To successfully <tt>evaluateDerivative2V xs ys x</tt>, the vectors
--   <tt>xs</tt> and <tt>ys</tt> must have identical lengths, and
--   <tt>xs</tt> must be monotonically increasing. The evaluation point
--   <tt>x</tt> must lie between the smallest and largest values in
--   <tt>xs</tt>.
evaluateDerivative2V :: InterpolationMethod -> Vector Double -> Vector Double -> Double -> Double

-- | Evaluate the definite integral of a function by interpolating within
--   the given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = [1..10]
--   
--   &gt;&gt;&gt; let ys = map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateIntegralV CSpline (zip xs ys) (2.2, 5.5)
--   51.909
--   </pre>
--   
--   To successfully <tt>evaluateIntegral points (a, b)</tt>, the domain
--   (<tt>x</tt>) values of <tt>points</tt> must be monotonically
--   increasing. The integration bounds <tt>a</tt> and <tt>b</tt> must lie
--   between the smallest and largest values in the sampled domain..
evaluateIntegral :: InterpolationMethod -> [(Double, Double)] -> (Double, Double) -> Double

-- | Evaluate the definite integral of a function by interpolating within
--   the given dataset. For example:
--   
--   <pre>
--   &gt;&gt;&gt; let xs = vector [1..10]
--   
--   &gt;&gt;&gt; let ys = vector $ map (**2) [1..10]
--   
--   &gt;&gt;&gt; evaluateIntegralV CSpline xs ys 2.2 5.5
--   51.89853207547169
--   </pre>
--   
--   To successfully <tt>evaluateIntegralV xs ys a b</tt>, the vectors
--   <tt>xs</tt> and <tt>ys</tt> must have identical lengths, and
--   <tt>xs</tt> must be monotonically increasing. The integration bounds
--   <tt>a</tt> and <tt>b</tt> must lie between the smallest and largest
--   values in <tt>xs</tt>.
evaluateIntegralV :: InterpolationMethod -> Vector Double -> Vector Double -> Double -> Double -> Double
instance GHC.Classes.Eq Numeric.GSL.Interpolation.InterpolationMethod
instance GHC.Internal.Read.Read Numeric.GSL.Interpolation.InterpolationMethod
instance GHC.Internal.Show.Show Numeric.GSL.Interpolation.InterpolationMethod


module Numeric.GSL.LinearAlgebra
data RandDist

-- | uniform distribution in [0,1)
Uniform :: RandDist

-- | normal distribution with mean zero and standard deviation one
Gaussian :: RandDist

-- | Obtains a vector of pseudorandom elements from the the mt19937
--   generator in GSL, with a given seed. Use randomIO to get a random
--   seed.
randomVector :: Int -> RandDist -> Int -> Vector Double

-- | Saves a matrix as 2D ASCII table.
saveMatrix :: FilePath -> String -> Matrix Double -> IO ()

-- | Saves the elements of a vector to a binary file.
fwriteVector :: FilePath -> Vector Double -> IO ()

-- | Loads a vector from a binary file (the number of elements must be
--   known in advance).
freadVector :: FilePath -> Int -> IO (Vector Double)

-- | Saves the elements of a vector, with a given format (%f, %e, %g), to
--   an ASCII file.
fprintfVector :: FilePath -> String -> Vector Double -> IO ()

-- | Loads a vector from an ASCII file (the number of elements must be
--   known in advance).
fscanfVector :: FilePath -> Int -> IO (Vector Double)

-- | obtains the number of rows and columns in an ASCII data file
--   (provisionally using unix's wc).
fileDimensions :: FilePath -> IO (Int, Int)

-- | Loads a matrix from an ASCII file formatted as a 2D table.
loadMatrix :: FilePath -> IO (Matrix Double)

-- | Loads a matrix from an ASCII file (the number of rows and columns must
--   be known in advance).
fromFile :: FilePath -> (Int, Int) -> IO (Matrix Double)
instance GHC.Internal.Enum.Enum Numeric.GSL.LinearAlgebra.RandDist


-- | Minimization of a multidimensional function using some of the
--   algorithms described in:
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Multidimensional-Minimization.html</a>
--   
--   The example in the GSL manual:
--   
--   <pre>
--   f [x,y] = 10*(x-1)^2 + 20*(y-2)^2 + 30
--   
--   main = do
--       let (s,p) = minimize NMSimplex2 1E-2 30 [1,1] f [5,7]
--       print s
--       print p
--   </pre>
--   
--   <pre>
--   &gt;&gt;&gt; main
--   [0.9920430849306288,1.9969168063253182]
--    0.000  512.500  1.130  6.500  5.000
--    1.000  290.625  1.409  5.250  4.000
--    2.000  290.625  1.409  5.250  4.000
--    3.000  252.500  1.409  5.500  1.000
--    ...
--   22.000   30.001  0.013  0.992  1.997
--   23.000   30.001  0.008  0.992  1.997
--   </pre>
--   
--   The path to the solution can be graphically shown by means of:
--   
--   <pre>
--   <a>mplot</a> $ drop 3 (<a>toColumns</a> p)
--   </pre>
--   
--   Taken from the GSL manual:
--   
--   The vector Broyden-Fletcher-Goldfarb-Shanno (BFGS) algorithm is a
--   quasi-Newton method which builds up an approximation to the second
--   derivatives of the function f using the difference between successive
--   gradient vectors. By combining the first and second derivatives the
--   algorithm is able to take Newton-type steps towards the function
--   minimum, assuming quadratic behavior in that region.
--   
--   The bfgs2 version of this minimizer is the most efficient version
--   available, and is a faithful implementation of the line minimization
--   scheme described in Fletcher's Practical Methods of Optimization,
--   Algorithms 2.6.2 and 2.6.4. It supercedes the original bfgs routine
--   and requires substantially fewer function and gradient evaluations.
--   The user-supplied tolerance tol corresponds to the parameter sigma
--   used by Fletcher. A value of 0.1 is recommended for typical use
--   (larger values correspond to less accurate line searches).
--   
--   The nmsimplex2 version is a new O(N) implementation of the earlier
--   O(N^2) nmsimplex minimiser. It calculates the size of simplex as the
--   rms distance of each vertex from the center rather than the mean
--   distance, which has the advantage of allowing a linear update.
module Numeric.GSL.Minimization

-- | Minimization without derivatives
minimize :: MinimizeMethod -> Double -> Int -> [Double] -> ([Double] -> Double) -> [Double] -> ([Double], Matrix Double)

-- | Minimization without derivatives (vector version)
minimizeV :: MinimizeMethod -> Double -> Int -> Vector Double -> (Vector Double -> Double) -> Vector Double -> (Vector Double, Matrix Double)
data MinimizeMethod
NMSimplex :: MinimizeMethod
NMSimplex2 :: MinimizeMethod

-- | Minimization with derivatives.
minimizeD :: MinimizeMethodD -> Double -> Int -> Double -> Double -> ([Double] -> Double) -> ([Double] -> [Double]) -> [Double] -> ([Double], Matrix Double)

-- | Minimization with derivatives (vector version)
minimizeVD :: MinimizeMethodD -> Double -> Int -> Double -> Double -> (Vector Double -> Double) -> (Vector Double -> Vector Double) -> Vector Double -> (Vector Double, Matrix Double)
data MinimizeMethodD
ConjugateFR :: MinimizeMethodD
ConjugatePR :: MinimizeMethodD
VectorBFGS :: MinimizeMethodD
VectorBFGS2 :: MinimizeMethodD
SteepestDescent :: MinimizeMethodD

-- | Onedimensional minimization.
uniMinimize :: UniMinimizeMethod -> Double -> Int -> (Double -> Double) -> Double -> Double -> Double -> (Double, Matrix Double)
data UniMinimizeMethod
GoldenSection :: UniMinimizeMethod
BrentMini :: UniMinimizeMethod
QuadGolden :: UniMinimizeMethod

-- | <i>Deprecated: use minimize NMSimplex2 eps maxit sizes f xi</i>
minimizeNMSimplex :: ([Double] -> Double) -> [Double] -> [Double] -> Double -> Int -> ([Double], Matrix Double)

-- | <i>Deprecated: use minimizeD ConjugateFR eps maxit step tol f g xi</i>
minimizeConjugateGradient :: Double -> Double -> Double -> Int -> ([Double] -> Double) -> ([Double] -> [Double]) -> [Double] -> ([Double], Matrix Double)

-- | <i>Deprecated: use minimizeD VectorBFGS2 eps maxit step tol f g xi</i>
minimizeVectorBFGS2 :: Double -> Double -> Double -> Int -> ([Double] -> Double) -> ([Double] -> [Double]) -> [Double] -> ([Double], Matrix Double)
instance GHC.Internal.Enum.Bounded Numeric.GSL.Minimization.MinimizeMethod
instance GHC.Internal.Enum.Bounded Numeric.GSL.Minimization.MinimizeMethodD
instance GHC.Internal.Enum.Bounded Numeric.GSL.Minimization.UniMinimizeMethod
instance GHC.Internal.Enum.Enum Numeric.GSL.Minimization.MinimizeMethod
instance GHC.Internal.Enum.Enum Numeric.GSL.Minimization.MinimizeMethodD
instance GHC.Internal.Enum.Enum Numeric.GSL.Minimization.UniMinimizeMethod
instance GHC.Classes.Eq Numeric.GSL.Minimization.MinimizeMethod
instance GHC.Classes.Eq Numeric.GSL.Minimization.MinimizeMethodD
instance GHC.Classes.Eq Numeric.GSL.Minimization.UniMinimizeMethod
instance GHC.Internal.Show.Show Numeric.GSL.Minimization.MinimizeMethod
instance GHC.Internal.Show.Show Numeric.GSL.Minimization.MinimizeMethodD
instance GHC.Internal.Show.Show Numeric.GSL.Minimization.UniMinimizeMethod


-- | Solution of ordinary differential equation (ODE) initial value
--   problems.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Ordinary-Differential-Equations.html</a>
--   
--   A simple example:
--   
--   <pre>
--   import Numeric.GSL.ODE
--   import Numeric.LinearAlgebra
--   import Graphics.Plot(mplot)
--   
--   xdot t [x,v] = [v, -0.95*x - 0.1*v]
--   
--   ts = linspace 100 (0,20 :: Double)
--   
--   sol = odeSolve xdot [10,0] ts
--   
--   main = mplot (ts : toColumns sol)
--   </pre>
module Numeric.GSL.ODE

-- | A version of <a>odeSolveV</a> with reasonable default parameters and
--   system of equations defined using lists.
odeSolve :: (Double -> [Double] -> [Double]) -> [Double] -> Vector Double -> Matrix Double

-- | A version of <a>odeSolveVWith</a> with reasonable default step
--   control.
odeSolveV :: ODEMethod -> Double -> Double -> Double -> (Double -> Vector Double -> Vector Double) -> Vector Double -> Vector Double -> Matrix Double

-- | Evolution of the system with adaptive step-size control.
odeSolveVWith :: ODEMethod -> StepControl -> Double -> (Double -> Vector Double -> Vector Double) -> Vector Double -> Vector Double -> Matrix Double

-- | Stepping functions
data ODEMethod

-- | Embedded Runge-Kutta (2, 3) method.
RK2 :: ODEMethod

-- | 4th order (classical) Runge-Kutta. The error estimate is obtained by
--   halving the step-size. For more efficient estimate of the error, use
--   the embedded methods.
RK4 :: ODEMethod

-- | Embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good
--   general-purpose integrator.
RKf45 :: ODEMethod

-- | Embedded Runge-Kutta Cash-Karp (4, 5) method.
RKck :: ODEMethod

-- | Embedded Runge-Kutta Prince-Dormand (8,9) method.
RK8pd :: ODEMethod

-- | Implicit 2nd order Runge-Kutta at Gaussian points.
RK2imp :: Jacobian -> ODEMethod

-- | Implicit 4th order Runge-Kutta at Gaussian points.
RK4imp :: Jacobian -> ODEMethod

-- | Implicit Bulirsch-Stoer method of Bader and Deuflhard. The method is
--   generally suitable for stiff problems.
BSimp :: Jacobian -> ODEMethod

-- | Implicit Gaussian first order Runge-Kutta. Also known as implicit
--   Euler or backward Euler method. Error estimation is carried out by the
--   step doubling method.
RK1imp :: Jacobian -> ODEMethod

-- | A variable-coefficient linear multistep Adams method in Nordsieck
--   form. This stepper uses explicit Adams-Bashforth (predictor) and
--   implicit Adams-Moulton (corrector) methods in P(EC)^m functional
--   iteration mode. Method order varies dynamically between 1 and 12.
MSAdams :: ODEMethod

-- | A variable-coefficient linear multistep backward differentiation
--   formula (BDF) method in Nordsieck form. This stepper uses the explicit
--   BDF formula as predictor and implicit BDF formula as corrector. A
--   modified Newton iteration method is used to solve the system of
--   non-linear equations. Method order varies dynamically between 1 and 5.
--   The method is generally suitable for stiff problems.
MSBDF :: Jacobian -> ODEMethod
type Jacobian = Double -> Vector Double -> Matrix Double

-- | Adaptive step-size control functions
data StepControl

-- | abs. and rel. tolerance for x(t)
X :: Double -> Double -> StepControl

-- | abs. and rel. tolerance for x'(t)
X' :: Double -> Double -> StepControl

-- | include both via rel. tolerance scaling factors a_x, a_x'
XX' :: Double -> Double -> Double -> Double -> StepControl

-- | scale abs. tolerance of x(t) components
ScXX' :: Double -> Double -> Double -> Double -> Vector Double -> StepControl


-- | Polynomials.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/General-Polynomial-Equations.html#General-Polynomial-Equations</a>
module Numeric.GSL.Polynomials

-- | Solution of general polynomial equations, using
--   <i>gsl_poly_complex_solve</i>.
--   
--   For example, the three solutions of x^3 + 8 = 0
--   
--   <pre>
--   &gt;&gt;&gt; polySolve [8,0,0,1]
--   [(-2.0) :+ 0.0,1.0 :+ 1.7320508075688776,1.0 :+ (-1.7320508075688776)]
--   </pre>
--   
--   The example in the GSL manual: To find the roots of x^5 -1 = 0:
--   
--   <pre>
--   &gt;&gt;&gt; polySolve [-1, 0, 0, 0, 0, 1]
--   [(-0.8090169943749472) :+ 0.5877852522924731,
--   (-0.8090169943749472) :+ (-0.5877852522924731),
--   0.30901699437494756 :+ 0.9510565162951535,
--   0.30901699437494756 :+ (-0.9510565162951535),
--   1.0000000000000002 :+ 0.0]
--   </pre>
polySolve :: [Double] -> [Complex Double]


-- | Multidimensional root finding.
--   
--   
--   <a>http://www.gnu.org/software/gsl/manual/html_node/Multidimensional-Root_002dFinding.html</a>
--   
--   The example in the GSL manual:
--   
--   <pre>
--   &gt;&gt;&gt; let rosenbrock a b [x,y] = [ a*(1-x), b*(y-x^2) ]
--   
--   &gt;&gt;&gt; let (sol,path) = root Hybrids 1E-7 30 (rosenbrock 1 10) [-10,-5]
--   
--   &gt;&gt;&gt; sol
--   [1.0,1.0]
--   
--   &gt;&gt;&gt; disp 3 path
--   11x5
--    1.000  -10.000  -5.000  11.000  -1050.000
--    2.000   -3.976  24.827   4.976     90.203
--    3.000   -3.976  24.827   4.976     90.203
--    4.000   -3.976  24.827   4.976     90.203
--    5.000   -1.274  -5.680   2.274    -73.018
--    6.000   -1.274  -5.680   2.274    -73.018
--    7.000    0.249   0.298   0.751      2.359
--    8.000    0.249   0.298   0.751      2.359
--    9.000    1.000   0.878  -0.000     -1.218
--   10.000    1.000   0.989  -0.000     -0.108
--   11.000    1.000   1.000   0.000      0.000
--   </pre>
module Numeric.GSL.Root
uniRoot :: UniRootMethod -> Double -> Int -> (Double -> Double) -> Double -> Double -> (Double, Matrix Double)
data UniRootMethod
Bisection :: UniRootMethod
FalsePos :: UniRootMethod
Brent :: UniRootMethod
uniRootJ :: UniRootMethodJ -> Double -> Int -> (Double -> Double) -> (Double -> Double) -> Double -> (Double, Matrix Double)
data UniRootMethodJ
UNewton :: UniRootMethodJ
Secant :: UniRootMethodJ
Steffenson :: UniRootMethodJ

-- | Nonlinear multidimensional root finding using algorithms that do not
--   require any derivative information to be supplied by the user. Any
--   derivatives needed are approximated by finite differences.
root :: RootMethod -> Double -> Int -> ([Double] -> [Double]) -> [Double] -> ([Double], Matrix Double)
data RootMethod
Hybrids :: RootMethod
Hybrid :: RootMethod
DNewton :: RootMethod
Broyden :: RootMethod

-- | Nonlinear multidimensional root finding using both the function and
--   its derivatives.
rootJ :: RootMethodJ -> Double -> Int -> ([Double] -> [Double]) -> ([Double] -> [[Double]]) -> [Double] -> ([Double], Matrix Double)
data RootMethodJ
HybridsJ :: RootMethodJ
HybridJ :: RootMethodJ
Newton :: RootMethodJ
GNewton :: RootMethodJ
instance GHC.Internal.Enum.Bounded Numeric.GSL.Root.RootMethod
instance GHC.Internal.Enum.Bounded Numeric.GSL.Root.RootMethodJ
instance GHC.Internal.Enum.Bounded Numeric.GSL.Root.UniRootMethod
instance GHC.Internal.Enum.Bounded Numeric.GSL.Root.UniRootMethodJ
instance GHC.Internal.Enum.Enum Numeric.GSL.Root.RootMethod
instance GHC.Internal.Enum.Enum Numeric.GSL.Root.RootMethodJ
instance GHC.Internal.Enum.Enum Numeric.GSL.Root.UniRootMethod
instance GHC.Internal.Enum.Enum Numeric.GSL.Root.UniRootMethodJ
instance GHC.Classes.Eq Numeric.GSL.Root.RootMethod
instance GHC.Classes.Eq Numeric.GSL.Root.RootMethodJ
instance GHC.Classes.Eq Numeric.GSL.Root.UniRootMethod
instance GHC.Classes.Eq Numeric.GSL.Root.UniRootMethodJ
instance GHC.Internal.Show.Show Numeric.GSL.Root.RootMethod
instance GHC.Internal.Show.Show Numeric.GSL.Root.RootMethodJ
instance GHC.Internal.Show.Show Numeric.GSL.Root.UniRootMethod
instance GHC.Internal.Show.Show Numeric.GSL.Root.UniRootMethodJ


-- | This module reexports all available GSL functions.
--   
--   The GSL special functions are in the separate package hmatrix-special.
module Numeric.GSL

-- | This action removes the GSL default error handler (which aborts the
--   program), so that GSL errors can be handled by Haskell (using
--   Control.Exception) and ghci doesn't abort.
setErrorHandlerOff :: IO ()


-- | Simulated annealing routines.
--   
--   
--   <a>https://www.gnu.org/software/gsl/manual/html_node/Simulated-Annealing.html#Simulated-Annealing</a>
--   
--   Here is a translation of the simple example given in <a>the GSL
--   manual</a>:
--   
--   <pre>
--   import Numeric.GSL.SimulatedAnnealing
--   import Numeric.LinearAlgebra.HMatrix
--   
--   main = print $ simanSolve 0 1 exampleParams 15.5 exampleE exampleM exampleS (Just show)
--   
--   exampleParams = SimulatedAnnealingParams 200 1000 1.0 1.0 0.008 1.003 2.0e-6
--   
--   exampleE x = exp (-(x - 1)**2) * sin (8 * x)
--   
--   exampleM x y = abs $ x - y
--   
--   exampleS rands stepSize current = (rands ! 0) * 2 * stepSize - stepSize + current
--   </pre>
--   
--   The manual states:
--   
--   <pre>
--   The first example, in one dimensional Cartesian space, sets up an
--   energy function which is a damped sine wave; this has many local
--   minima, but only one global minimum, somewhere between 1.0 and
--   1.5. The initial guess given is 15.5, which is several local minima
--   away from the global minimum.
--   </pre>
--   
--   This global minimum is around 1.36.
module Numeric.GSL.SimulatedAnnealing

-- | Calling
--   
--   <pre>
--   simanSolve seed nrand params x0 e m step print
--   </pre>
--   
--   performs a simulated annealing search through a given space. So that
--   any configuration type may be used, the space is specified by
--   providing the functions <tt>e</tt> (the energy functional) and
--   <tt>m</tt> (the metric definition). <tt>x0</tt> is the initial
--   configuration of the system. The simulated annealing steps are
--   generated using the user-provided function <tt>step</tt>, which should
--   randomly construct a new system configuration.
--   
--   If <a>Nothing</a> is passed instead of a printing function, no
--   incremental output will be generated. Otherwise, the GSL-formatted
--   output, including the configuration description the user function
--   generates, will be printed to stdout.
--   
--   Each time the step function is called, it is supplied with a random
--   vector containing <tt>nrand</tt> <a>Double</a> values, uniformly
--   distributed in <tt>[0, 1)</tt>. It should use these values to generate
--   its new configuration.
simanSolve :: Int -> Int -> SimulatedAnnealingParams -> a -> (a -> Double) -> (a -> a -> Double) -> (Vector Double -> Double -> a -> a) -> Maybe (a -> String) -> a

-- | <a>SimulatedAnnealingParams</a> is a translation of the
--   <tt>gsl_siman_params_t</tt> structure documented in <a>the GSL
--   manual</a>, which controls the simulated annealing algorithm.
--   
--   The annealing process is parameterized by the Boltzmann distribution
--   and the <i>cooling schedule</i>. For more details, see <a>the relevant
--   section of the manual</a>.
data SimulatedAnnealingParams
SimulatedAnnealingParams :: CInt -> CInt -> Double -> Double -> Double -> Double -> Double -> SimulatedAnnealingParams

-- | The number of points to try for each step.
[n_tries] :: SimulatedAnnealingParams -> CInt

-- | The number of iterations at each temperature
[iters_fixed_T] :: SimulatedAnnealingParams -> CInt

-- | The maximum step size in the random walk
[step_size] :: SimulatedAnnealingParams -> Double

-- | Boltzmann distribution parameter
[boltzmann_k] :: SimulatedAnnealingParams -> Double

-- | Initial temperature
[cooling_t_initial] :: SimulatedAnnealingParams -> Double

-- | Cooling rate parameter
[cooling_mu_t] :: SimulatedAnnealingParams -> Double

-- | Final temperature
[cooling_t_min] :: SimulatedAnnealingParams -> Double
instance GHC.Classes.Eq Numeric.GSL.SimulatedAnnealing.SimulatedAnnealingParams
instance GHC.Internal.Read.Read Numeric.GSL.SimulatedAnnealing.SimulatedAnnealingParams
instance GHC.Internal.Show.Show Numeric.GSL.SimulatedAnnealing.SimulatedAnnealingParams
instance GHC.Internal.Foreign.Storable.Storable Numeric.GSL.SimulatedAnnealing.SimulatedAnnealingParams
