]> gitweb.michael.orlitzky.com - numerical-analysis.git/blob - src/Roots/Simple.hs
src/Roots/Simple.hs: fix monomorphism restriction warning.
[numerical-analysis.git] / src / Roots / Simple.hs
1 {-# LANGUAGE RebindableSyntax #-}
2 {-# LANGUAGE ScopedTypeVariables #-}
3
4 -- | The Roots.Simple module contains root-finding algorithms. That
5 -- is, procedures to (numerically) find solutions to the equation,
6 --
7 -- > f(x) = 0
8 --
9 -- where f is assumed to be continuous on the interval of interest.
10 --
11
12 module Roots.Simple (
13 bisect,
14 fixed_point,
15 fixed_point_error_ratios,
16 fixed_point_iteration_count,
17 has_root,
18 newtons_method,
19 secant_method,
20 trisect )
21 where
22
23 import Data.List (find)
24 import NumericPrelude hiding ( abs )
25 import Algebra.Absolute ( abs )
26 import qualified Algebra.Additive as Additive ( C )
27 import qualified Algebra.Algebraic as Algebraic ( C )
28 import qualified Algebra.Field as Field ( C )
29 import qualified Algebra.RealField as RealField ( C )
30 import qualified Algebra.RealRing as RealRing ( C )
31
32 import Normed ( Normed(..) )
33 import qualified Roots.Fast as F (
34 bisect,
35 fixed_point_iterations,
36 fixed_point_with_iterations,
37 has_root,
38 trisect )
39
40
41 -- | Does the (continuous) function @f@ have a root on the interval
42 -- [a,b]? If f(a) <] 0 and f(b) ]> 0, we know that there's a root in
43 -- [a,b] by the intermediate value theorem. Likewise when f(a) >= 0
44 -- and f(b) <= 0.
45 --
46 -- Examples:
47 --
48 -- >>> let f x = x**3
49 -- >>> has_root f (-1) 1 Nothing
50 -- True
51 --
52 -- This fails if we don't specify an @epsilon@, because cos(-2) ==
53 -- cos(2) doesn't imply that there's a root on [-2,2].
54 --
55 -- >>> has_root cos (-2) 2 Nothing
56 -- False
57 -- >>> has_root cos (-2) 2 (Just 0.001)
58 -- True
59 --
60 has_root :: (RealField.C a, RealRing.C b)
61 => (a -> b) -- ^ The function @f@
62 -> a -- ^ The \"left\" endpoint, @a@
63 -> a -- ^ The \"right\" endpoint, @b@
64 -> Maybe a -- ^ The size of the smallest subinterval
65 -- we'll examine, @epsilon@
66 -> Bool
67 has_root f a b epsilon =
68 F.has_root f a b epsilon Nothing Nothing
69
70
71 -- | We are given a function @f@ and an interval [a,b]. The bisection
72 -- method finds a root by splitting [a,b] in half repeatedly.
73 --
74 -- If one is found within some prescribed tolerance @epsilon@, it is
75 -- returned. Otherwise, the interval [a,b] is split into two
76 -- subintervals [a,c] and [c,b] of equal length which are then both
77 -- checked via the same process.
78 --
79 -- Returns 'Just' the value x for which f(x) == 0 if one is found,
80 -- or Nothing if one of the preconditions is violated.
81 --
82 -- Examples:
83 --
84 -- >>> let actual = 1.5707963267948966
85 -- >>> let Just root = bisect cos 1 2 0.001
86 -- >>> root
87 -- 1.5712890625
88 -- >>> abs (root - actual) < 0.001
89 -- True
90 --
91 -- >>> bisect sin (-1) 1 0.001
92 -- Just 0.0
93 --
94 bisect :: (RealField.C a, RealRing.C b)
95 => (a -> b) -- ^ The function @f@ whose root we seek
96 -> a -- ^ The \"left\" endpoint of the interval, @a@
97 -> a -- ^ The \"right\" endpoint of the interval, @b@
98 -> a -- ^ The tolerance, @epsilon@
99 -> Maybe a
100 bisect f a b epsilon =
101 F.bisect f a b epsilon Nothing Nothing
102
103
104 -- | We are given a function @f@ and an interval [a,b]. The trisection
105 -- method finds a root by splitting [a,b] into three
106 -- subintervals repeatedly.
107 --
108 -- If one is found within some prescribed tolerance @epsilon@, it is
109 -- returned. Otherwise, the interval [a,b] is split into two
110 -- subintervals [a,c] and [c,b] of equal length which are then both
111 -- checked via the same process.
112 --
113 -- Returns 'Just' the value x for which f(x) == 0 if one is found,
114 -- or Nothing if one of the preconditions is violated.
115 --
116 -- Examples:
117 --
118 -- >>> let actual = 1.5707963267948966
119 -- >>> let Just root = trisect cos 1 2 0.001
120 -- >>> root
121 -- 1.5713305898491083
122 -- >>> abs (root - actual) < 0.001
123 -- True
124 --
125 -- >>> trisect sin (-1) 1 0.001
126 -- Just 0.0
127 --
128 trisect :: (RealField.C a, RealRing.C b)
129 => (a -> b) -- ^ The function @f@ whose root we seek
130 -> a -- ^ The \"left\" endpoint of the interval, @a@
131 -> a -- ^ The \"right\" endpoint of the interval, @b@
132 -> a -- ^ The tolerance, @epsilon@
133 -> Maybe a
134 trisect f a b epsilon =
135 F.trisect f a b epsilon Nothing Nothing
136
137
138 -- | Find a fixed point of the function @f@ with the search starting
139 -- at x0. We delegate to the version that returns the number of
140 -- iterations and simply discard the number of iterations.
141 --
142 fixed_point :: (Normed a, Additive.C a, Algebraic.C b, RealField.C b)
143 => (a -> a) -- ^ The function @f@ to iterate.
144 -> b -- ^ The tolerance, @epsilon@.
145 -> a -- ^ The initial value @x0@.
146 -> a -- ^ The fixed point.
147 fixed_point f epsilon x0 =
148 snd $ F.fixed_point_with_iterations f epsilon x0
149
150
151 -- | Return the number of iterations required to find a fixed point of
152 -- the function @f@ with the search starting at x0 and tolerance
153 -- @epsilon@. We delegate to the version that returns the number of
154 -- iterations and simply discard the fixed point.
155 fixed_point_iteration_count :: (Normed a,
156 Additive.C a,
157 RealField.C b,
158 Algebraic.C b)
159 => (a -> a) -- ^ The function @f@ to iterate.
160 -> b -- ^ The tolerance, @epsilon@.
161 -> a -- ^ The initial value @x0@.
162 -> Int -- ^ The fixed point.
163 fixed_point_iteration_count f epsilon x0 =
164 fst $ F.fixed_point_with_iterations f epsilon x0
165
166
167 -- | Returns a list of ratios,
168 --
169 -- ||x^{*} - x_{n+1}|| / ||x^{*} - x_{n}||^{p}
170 --
171 -- of fixed point iterations for the function @f@ with initial guess
172 -- @x0@ and @p@ some positive power.
173 --
174 -- This is used to determine the rate of convergence.
175 --
176 fixed_point_error_ratios :: forall a b. (Normed a,
177 Additive.C a,
178 RealField.C b,
179 Algebraic.C b)
180 => (a -> a) -- ^ The function @f@ to iterate.
181 -> a -- ^ The initial value @x0@.
182 -> a -- ^ The true solution, @x_star@.
183 -> Integer -- ^ The power @p@.
184 -> [b] -- ^ The resulting sequence of x_{n}.
185 fixed_point_error_ratios f x0 x_star p =
186 zipWith (/) en_plus_one en_exp
187 where
188 xn = F.fixed_point_iterations f x0
189 en = map (\x -> norm (x_star - x)) xn :: [b]
190 en_plus_one = tail en
191 en_exp = map (^p) en
192
193
194
195 -- | The sequence x_{n} of values obtained by applying Newton's method
196 -- on the function @f@ and initial guess @x0@.
197 --
198 -- Examples:
199 --
200 -- Atkinson, p. 60.
201 --
202 -- >>> let f x = x^6 - x - 1
203 -- >>> let f' x = 6*x^5 - 1
204 -- >>> tail $ take 4 $ newton_iterations f f' 2
205 -- [1.6806282722513088,1.4307389882390624,1.2549709561094362]
206 --
207 newton_iterations :: (Field.C a)
208 => (a -> a) -- ^ The function @f@ whose root we seek
209 -> (a -> a) -- ^ The derivative of @f@
210 -> a -- ^ Initial guess, x-naught
211 -> [a]
212 newton_iterations f f' =
213 iterate next
214 where
215 next xn =
216 xn - ( (f xn) / (f' xn) )
217
218
219 -- | Use Newton's method to find a root of @f@ near the initial guess
220 -- @x0@. If your guess is bad, this will recurse forever!
221 --
222 -- Examples:
223 --
224 -- Atkinson, p. 60.
225 --
226 -- >>> let f x = x^6 - x - 1
227 -- >>> let f' x = 6*x^5 - 1
228 -- >>> let Just root = newtons_method f f' (1/1000000) 2
229 -- >>> root
230 -- 1.1347241385002211
231 -- >>> abs (f root) < 1/100000
232 -- True
233 --
234 -- >>> import Data.Number.BigFloat
235 -- >>> let eps = 1/(10^20) :: BigFloat Prec50
236 -- >>> let Just root = newtons_method f f' eps 2
237 -- >>> root
238 -- 1.13472413840151949260544605450647284028100785303643e0
239 -- >>> abs (f root) < eps
240 -- True
241 --
242 newtons_method :: (RealField.C a)
243 => (a -> a) -- ^ The function @f@ whose root we seek
244 -> (a -> a) -- ^ The derivative of @f@
245 -> a -- ^ The tolerance epsilon
246 -> a -- ^ Initial guess, x-naught
247 -> Maybe a
248 newtons_method f f' epsilon x0 =
249 find (\x -> abs (f x) < epsilon) x_n
250 where
251 x_n = newton_iterations f f' x0
252
253
254 -- | Takes a function @f@ of two arguments and repeatedly applies @f@
255 -- to the previous two values. Returns a list containing all
256 -- generated values, f(x0, x1), f(x1, x2), f(x2, x3)...
257 --
258 -- Examples:
259 --
260 -- >>> let fibs = iterate2 (+) 0 1
261 -- >>> take 15 fibs
262 -- [0,1,1,2,3,5,8,13,21,34,55,89,144,233,377]
263 --
264 iterate2 :: (a -> a -> a) -- ^ The function @f@
265 -> a -- ^ The initial value @x0@
266 -> a -- ^ The second value, @x1@
267 -> [a] -- ^ The result list, [x0, x1, ...]
268 iterate2 f x0 x1 =
269 x0 : x1 : (go x0 x1)
270 where
271 go prev2 prev1 =
272 let next = f prev2 prev1 in
273 next : go prev1 next
274
275
276 -- | The sequence x_{n} of values obtained by applying the secant
277 -- method on the function @f@ and initial guesses @x0@, @x1@.
278 --
279 -- The recursion more or less implements a two-parameter 'iterate',
280 -- although one list is passed to the next iteration (as opposed to
281 -- one function argument, with iterate). At each step, we peel the
282 -- first two elements off the list and then compute/append elements
283 -- three, four... onto the end of the list.
284 --
285 -- Examples:
286 --
287 -- Atkinson, p. 67.
288 --
289 -- >>> let f x = x^6 - x - 1
290 -- >>> take 4 $ secant_iterations f 2 1
291 -- [2.0,1.0,1.0161290322580645,1.190577768676638]
292 --
293 secant_iterations :: (Field.C a)
294 => (a -> a) -- ^ The function @f@ whose root we seek
295 -> a -- ^ Initial guess, x-naught
296 -> a -- ^ Second initial guess, x-one
297 -> [a]
298 secant_iterations f =
299 iterate2 g
300 where
301 g prev2 prev1 =
302 let x_change = prev1 - prev2
303 y_change = (f prev1) - (f prev2)
304 in
305 (prev1 - (f prev1 * (x_change / y_change)))
306
307
308 -- | Use the secant method to find a root of @f@ near the initial guesses
309 -- @x0@ and @x1@. If your guesses are bad, this will recurse forever!
310 --
311 -- Examples:
312 --
313 -- Atkinson, p. 67.
314 --
315 -- >>> let f x = x^6 - x - 1
316 -- >>> let Just root = secant_method f (1/10^9) 2 1
317 -- >>> root
318 -- 1.1347241384015196
319 -- >>> abs (f root) < (1/10^9)
320 -- True
321 --
322 secant_method :: (RealField.C a)
323 => (a -> a) -- ^ The function @f@ whose root we seek
324 -> a -- ^ The tolerance epsilon
325 -> a -- ^ Initial guess, x-naught
326 -> a -- ^ Second initial guess, x-one
327 -> Maybe a
328 secant_method f epsilon x0 x1
329 = find (\x -> abs (f x) < epsilon) x_n
330 where
331 x_n = secant_iterations f x0 x1