Change name of Maxima function 'chinese' to 'solve_congruences'.
[maxima/cygwin.git] / src / specfn.lisp
blobf6b85ae3deaff0a521ff82865313f9641d80b946
1 ;;; -*- Mode: Lisp; Package: Maxima; Syntax: Common-Lisp; Base: 10 -*- ;;;;
2 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3 ;;; The data in this file contains enhancements. ;;;;;
4 ;;; ;;;;;
5 ;;; Copyright (c) 1984,1987 by William Schelter,University of Texas ;;;;;
6 ;;; All rights reserved ;;;;;
7 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
8 ;;; (c) Copyright 1980 Massachusetts Institute of Technology ;;;
9 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
11 (in-package :maxima)
13 (macsyma-module specfn)
15 ;;*********************************************************************
16 ;;**************** ******************
17 ;;**************** Macsyma Special Function Routines ******************
18 ;;**************** ******************
19 ;;*********************************************************************
21 (load-macsyma-macros rzmac)
22 (load-macsyma-macros mhayat)
24 (defmacro mnumericalp (arg)
25 `(or (floatp ,arg) (and (or $numer $float) (integerp ,arg))))
27 ;; subtitle polylogarithm routines
29 (declare-top (special tlist))
31 (defun lisimp (expr vestigial z)
32 (declare (ignore vestigial))
33 (let ((s (simpcheck (car (subfunsubs expr)) z))
34 ($zerobern t)
35 (a))
36 (subargcheck expr 1 1 '$li)
37 (setq a (simpcheck (car (subfunargs expr)) z))
38 (or (cond ((zerop1 a) a)
39 ((and (mnump s) (ratgreaterp s 1) (eql a 1))
40 ;; li[s](a) = zeta(s) if s > 1 and if a = 1. We
41 ;; simplify this only if s is a rational number and a
42 ;; is the integer 1.
43 (ftake '%zeta s))
44 ((not (integerp s)) ())
45 ((= s 1)
46 (if (onep1 a)
47 (simp-domain-error
48 (intl:gettext "li: li[~:M](~:M) is undefined.") s a)
49 (neg (take '(%log) (sub 1 a)))))
50 ((= s 0) (div a (sub 1 a)))
51 ((< s 0) (lisimp-negative-integer s a))
52 ((and (integerp a) (> s 1)
53 (cond ((= a 1) (take '(%zeta) s))
54 ((= a -1)
55 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
56 (mul (add -1 (inv (expt 2 (- s 1))))
57 (take '(%zeta) s))))))
58 ((= s 2) (li2simp a))
59 ((= s 3) (li3simp a))
60 ((or (complex-float-numerical-eval-p a)
61 (complex-bigfloat-numerical-eval-p a))
62 (cond ((bigfloat:= 1 (bigfloat:to a))
63 ;; li[s](1) -> zeta(s)
64 (let ((result ($zeta s)))
65 (if (floatp a)
66 ($float result)
67 ($bfloat result))))
68 ((bigfloat:= -1 (bigfloat:to a))
69 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
70 (let ((result (mul (add -1 (inv (expt 2 (- s 1))))
71 (take '(%zeta) s))))
72 (if (floatp a)
73 ($float result)
74 ($bfloat result))))
75 ((integerp s)
76 (to (bigfloat::li-s-simp s (bigfloat:to a)))))))
77 (eqtest (subfunmakes '$li (ncons s) (ncons a))
78 expr))))
80 ;; Expand the Polylogarithm li[s](z) for a negative integer parameter s.
81 (defun lisimp-negative-integer (s z)
82 (let ((n (- s)))
83 (mul (inv (power (sub 1 z) (+ n 1)))
84 (let ((index1 (gensumindex))
85 ($simpsum t))
86 (dosum
87 (mul (power z index1)
88 (let ((index2 (gensumindex)))
89 (dosum
90 (mul (power -1 (add index2 1))
91 (take '(%binomial) (+ n 1) (sub index2 1))
92 (power (add 1 (sub index1 index2)) n))
93 index2 1 index1 t)))
94 index1 1 n t)))))
96 (defun li2simp (arg)
97 (cond ((mnumericalp arg)
98 ;; When arg is a float or rational, use the original li2numer
99 ;; using Spences function.
100 (li2numer (float arg)))
101 ((complex-float-numerical-eval-p arg)
102 ;; For complex args that should should result in float
103 ;; answers, use bigfloat::li2numer.
104 (to (bigfloat::li2numer (bigfloat:to ($rectform ($float arg))))))
105 ((or (bigfloat-numerical-eval-p arg)
106 (complex-bigfloat-numerical-eval-p arg))
107 (to (bigfloat::li2numer (bigfloat:to ($rectform ($bfloat arg))))))
108 ((alike1 arg '((rat) 1 2))
109 (add (div (take '(%zeta) 2) 2)
110 (mul '((rat simp) -1 2)
111 (power (take '(%log) 2) 2))))))
113 (defun li3simp (arg)
114 (cond ((or (float-numerical-eval-p arg)
115 (complex-float-numerical-eval-p arg))
116 (to (bigfloat::li3numer (bigfloat:to ($rectform ($float arg))))))
117 ((or (bigfloat-numerical-eval-p arg)
118 (complex-bigfloat-numerical-eval-p arg))
119 (to (bigfloat::li3numer (bigfloat:to ($rectform ($bfloat arg))))))
120 ((alike1 arg '((rat) 1 2))
121 (add (mul '((rat simp) 7 8) (take '(%zeta) 3))
122 (mul (div (take '(%zeta) 2) -2) (take '(%log) 2))
123 (mul '((rat simp) 1 6) (power (take '(%log) 2) 3))))))
125 ;; exponent in first term of taylor expansion of $li is one
126 (defun li-ord (subl)
127 (declare (ignore subl))
128 (ncons (rcone)))
130 ;; taylor expansion of $li is its definition:
131 ;; x + x^2/2^s + x^3/3^s + ...
132 (defun exp$li-fun (pw subl l) ; l is a irrelevant here
133 (setq subl (car subl)) ; subl is subscript of li
134 (prog ((e 0) ; e is exponent of current term
135 npw) ; npw is exponent of last term needed
136 (declare (fixnum e))
137 (setq npw (/ (float (car pw)) (float (cdr pw))))
138 (setq
139 l (cons '((0 . 1) 0 . 1)
140 nil))
141 a (setq e (1+ e))
142 (if (> e npw) (return l)
143 (rplacd (last l)
144 `(((,e . 1)
145 . ,(prep1 (m^ e (m- subl)))))))
146 (go a)))
149 ;; computes first pw terms of asymptotic expansion of $li[s](z)
151 ;; pw should be < (1/2)*s or gamma term is undefined
153 ;; Wood, D.C. (June 1992). The Computation of Polylogarithms. Technical Report 15-92
154 ;; University of Kent Computing Laboratory.
155 ;; http://www.cs.kent.ac.uk/pubs/1992/110
156 ;; equation 11.1
157 (defun li-asymptotic-expansion (pw s z)
158 (m+l (loop for k from 0 to pw collect
159 (if (and ($integerp k)
160 (mgqp 0 (m+ s 1 (m* -2 k) )))
161 0 ;; gamma in denominator below is infinite, this term is 0
162 (m* (m^ -1 k)
163 (m- 1 (m^ 2 (m- 1 (m* 2 k))))
164 (m^ (m* 2 '$%pi) (m* 2 k))
165 (m// ($bern (m* 2 k))
166 `((mfactorial) ,(m* 2 k)))
167 (m// (m^ `((%log) ,(m- z)) (m- s (m* 2 k)))
168 ($gamma (m+ s 1 (m* -2 k)))))))))
170 ;; Numerical evaluation for Chebyschev expansions of the first kind
172 (defun cheby (x chebarr)
173 (let ((bn+2 0.0) (bn+1 0.0))
174 (do ((i (floor (aref chebarr 0)) (1- i)))
175 ((< i 1) (- bn+1 (* bn+2 x)))
176 (setq bn+2
177 (prog1 bn+1 (setq bn+1 (+ (aref chebarr i)
178 (- (* 2.0 x bn+1) bn+2))))))))
180 (defun cheby-prime (x chebarr)
181 (- (cheby x chebarr)
182 (* (aref chebarr 1) 0.5)))
184 ;; These should really be calculated with minimax rational approximations.
185 ;; Someone has done LI[2] already, and this should be updated; I haven't
186 ;; seen any results for LI[3] yet.
188 (defun li2numer (y)
189 ;; Spence's function can be used to compute li[2] for 0 <= x <= 1.
190 ;; To compute the rest, we need the following identities:
192 ;; li[2](x) = -li[2](1/x)-log(-x)^2/2-%pi^2/6
193 ;; li[2](x) = li[2](1/(1-x)) + log(1-x)*log((1-x)/x^2)/2 - %pi^2/6
195 ;; The first tells us how to compute li[2] for x > 1. The result is complex.
196 ;; For x < 0, the second can be used, and the result is real.
198 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog2/17/01/01/)
199 (labels ((li2 (x)
200 (cond ((< x 0)
201 (+ (li2 (/ (- 1 x)))
202 (* 0.5 (log (- 1 x)) (log (/ (- 1 x) (* x x))))
203 (- (/ (cl:expt (float pi) 2) 6))))
204 ((< x 1)
205 (slatec:dspenc x))
206 ((= x 1)
207 (/ (cl:expt (float pi) 2) 6))
209 ;; li[2](x) = -li[2](1/x)-log(-x)^2/2-%pi^2/6
210 (- (+ (li2 (/ x))
211 (/ (cl:expt (cl:log (- x)) 2) 2)
212 (/ (cl:expt (float pi) 2) 6)))))))
213 (complexify (li2 y))))
216 (defvar *li2* (make-array 15. :initial-contents '(14.0 1.93506430 .166073033 2.48793229e-2
217 4.68636196e-3 1.0016275e-3 2.32002196e-4
218 5.68178227e-5 1.44963006e-5 3.81632946e-6
219 1.02990426e-6 2.83575385e-7 7.9387055e-8
220 2.2536705e-8 6.474338e-9)
221 :element-type 'flonum))
224 (defvar *li3* (make-array 15. :initial-contents '(14.0 1.95841721 8.51881315e-2 8.55985222e-3
225 1.21177214e-3 2.07227685e-4 3.99695869e-5
226 8.38064066e-6 1.86848945e-6 4.36660867e-7
227 1.05917334e-7 2.6478920e-8 6.787e-9
228 1.776536e-9 4.73417e-10)
229 :element-type 'flonum))
231 (defvar *s12* (make-array 18. :initial-contents '(17.0 1.90361778 .431311318 .100022507
232 2.44241560e-2 6.22512464e-3 1.64078831e-3
233 4.44079203e-4 1.22774942e-4 3.45398128e-5
234 9.85869565e-6 2.84856995e-6 8.31708473e-7
235 2.45039499e-7 7.2764962e-8 2.1758023e-8 6.546158e-9
236 1.980328e-9)
237 :element-type 'flonum))
239 (defun chebyli2 (x)
240 (* x (cheby-prime (/ (1+ (* x 4)) 3) *li2*)))
242 (defun chebyli3 (x)
243 (* x (cheby-prime (/ (1+ (* 4 x)) 3) *li3*)))
245 (defun chebys12 (x)
246 (* (/ (expt x 2) 4)
247 (cheby-prime (/ (1+ (* 4 x)) 3) *s12*)))
249 ;; subtitle polygamma routines
251 ;; gross efficiency hack, exp is a function of *k*, *k* should be mbind'ed
253 (defun msum (exp lo hi)
254 (if (< hi lo)
256 (let ((sum 0))
257 (do ((*k* lo (1+ *k*)))
258 ((> *k* hi) sum)
259 (declare (special *k*))
260 (setq sum (add2 sum (meval exp)))))))
263 (defun pole-err (exp)
264 (cond (errorsw (throw 'errorsw t))
265 (t (merror (intl:gettext "Pole encountered in: ~M") exp))))
268 (defprop $psi psisimp specsimp)
270 ;; Integral of psi function psi[n](x)
271 (putprop '$psi
272 `((n x)
274 ,(lambda (n x)
275 (cond
276 ((and ($integerp n) (>= n 0))
277 (cond
278 ((= n 0) `((%log_gamma) ,x))
279 (t `((mqapply) (($psi array) ((mplus) -1 ,n)) ,x))))
280 (t nil))))
281 'integral)
283 (defun psisimp (expr a z)
284 (let ((s (simpcheck (car (subfunsubs expr)) z)))
285 (subargcheck expr 1 1 '$psi)
286 (setq a (simpcheck (car (subfunargs expr)) z))
287 (and (setq z (integer-representation-p a))
288 (< z 1)
289 (pole-err expr))
290 (eqtest (psisimp1 s a) expr)))
292 ;; This gets pretty hairy now.
294 (defun psisimp1 (s a)
295 (let ((*k*))
296 (declare (special *k*))
298 (and (integerp s) (>= s 0) (mnumericalp a)
299 (let (($float2bf t)) ($float (mfuncall '$bfpsi s a 18))))
300 (and (integerp s) (>= s 0) ($bfloatp a)
301 (mfuncall '$bfpsi s a $fpprec))
302 (and (not $numer) (not $float) (integerp s) (> s -1)
303 (cond
304 ((integerp a)
305 (and (not (> a $maxpsiposint)) ; integer values
306 (m*t (expt -1 s) (factorial s)
307 (m- (msum (inv (m^t '*k* (1+ s))) 1 (1- a))
308 (cond ((zerop s) '$%gamma)
309 (($zeta (1+ s))))))))
310 ((or (not (ratnump a)) (ratgreaterp a $maxpsiposint)) ())
311 ((ratgreaterp a 0)
312 (cond
313 ((ratgreaterp a 1)
314 (let* ((int ($entier a)) ; reduction to fractional values
315 (frac (m-t a int)))
316 (m+t
317 (psisimp1 s frac)
318 (if (> int $maxpsiposint)
319 (subfunmakes '$psi (ncons s) (ncons int))
320 (m*t (expt -1 s) (factorial s)
321 (msum (m^t (m+t (m-t a int) '*k*)
322 (1- (- s)))
323 0 (1- int)))))))
324 ((= s 0)
325 (let ((p (cadr a)) (q (caddr a)))
326 (cond
327 ((or (> p $maxpsifracnum)
328 (> q $maxpsifracdenom) (bignump p) (bignump q)) ())
329 ((and (= p 1)
330 (cond ((= q 2)
331 (m+ (m* -2 '((%log) 2)) (m- '$%gamma)))
332 ((= q 3)
333 (m+ (m* '((rat simp) -1 2)
334 (m^t 3 '((rat simp) -1 2)) '$%pi)
335 (m* '((rat simp) -3 2) '((%log) 3))
336 (m- '$%gamma)))
337 ((= q 4)
338 (m+ (m* '((rat simp) -1 2) '$%pi)
339 (m* -3 '((%log) 2)) (m- '$%gamma)))
340 ((= q 6)
341 (m- (m+ (m* '((rat simp) 3 2) '((%log) 3))
342 (m* 2 '((%log) 2))
343 (m* '((rat simp) 1 2) '$%pi
344 (m^t 3 '((rat simp) 1 2)))
345 '$%gamma))))))
346 ((and (= p 2) (= q 3))
347 (m+ (m* '((rat simp) 1 2)
348 (m^t 3 '((rat simp) -1 2)) '$%pi)
349 (m* '((rat simp) -3 2) '((%log) 3))
350 (m- '$%gamma)))
351 ((and (= p 3) (= q 4))
352 (m+ (m* '((rat simp) 1 2) '$%pi)
353 (m* -3 '((%log) 2)) (m- '$%gamma)))
354 ((and (= p 5) (= q 6))
355 (m- (m* '((rat simp) 1 2) '$%pi
356 (m^t 3 '((rat simp) 1 2)))
357 (m+ (m* '((rat simp) 3 2) '((%log) 3))
358 (m* 2 '((%log) 2))
359 '$%gamma)))
360 ;; Gauss's Formula
361 ((let ((f (m* `((%cos) ,(m* 2 a '$%pi '*k*))
362 `((%log) ,(m-t 2 (m* 2 `((%cos)
363 ,(m//t (m* 2 '$%pi '*k*)
364 q))))))))
365 (m+t (msum f 1 (1- (truncate q 2)))
366 (let ((*k* (truncate q 2)))
367 (declare (special *k*))
368 (m*t (meval f)
369 (cond ((oddp q) 1)
370 ('((rat simp) 1 2)))))
371 (m-t (m+ (m* '$%pi '((rat simp) 1 2)
372 `((%cot) ((mtimes simp) ,a $%pi)))
373 `((%log) ,q)
374 '$%gamma))))))))
375 ((alike1 a '((rat) 1 2))
376 (m*t (expt -1 (1+ s)) (factorial s)
377 (1- (expt 2 (1+ s))) (simplify ($zeta (1+ s)))))
378 ((and (ratgreaterp a '((rat) 1 2))
379 (ratgreaterp 1 a))
380 (m*t
381 (expt -1 s)
382 (m+t (psisimp1 s (m- 1 a))
383 (let ((dif (m* '$%pi
384 ($diff `((%cot) ,(m* '$%pi '$z)) '$z s)))
385 ($z (m-t a)))
386 (declare (special $z))
387 (meval dif)))))))
388 ((ratgreaterp a $maxpsinegint) ;;; Reflection Formula
389 (m*t
390 (expt -1 s)
391 (m+t (m+t (psisimp1 s (m- a))
392 (let ((dif (m* '$%pi
393 ($diff `((%cot) ,(m* '$%pi '$z)) '$z s)))
394 ($z (m-t a)))
395 (declare (special $z))
396 (meval dif)))
397 (m*t (factorial s) (m^t (m-t a) (1- (- s)))))))))
398 (subfunmakes '$psi (ncons s) (ncons a)))))
401 ;; subtitle polygamma tayloring routines
403 ;; These routines are specially coded to be as fast as possible given the
404 ;; current $TAYLOR; too bad they have to be so ugly.
406 (declare-top (special var))
408 (defun expgam-fun (pw temp)
409 (setq temp (get-datum (get-key-var (car var))))
410 (let-pw temp pw
411 (pstimes
412 (let-pw temp (e1+ pw)
413 (psexpt-fn (getexp-fun '(($psi) -1) var (e1+ pw))))
414 (make-ps var (ncons pw) '(((-1 . 1) 1 . 1))))))
416 (defun expplygam-funs (pw subl l) ; l is a irrelevant here
417 (setq subl (car subl))
418 (if (or (not (integerp subl)) (< subl -1))
419 (tay-err "Unable to expand at a subscript in")
420 (prog ((e 0) (sf-sign 0) npw sf-last)
421 (declare (fixnum e) (fixnum sf-sign))
422 (setq npw (/ (float (car pw)) (float (cdr pw))))
423 (setq
424 l (cond ((= subl -1)
425 `(((1 . 1) . ,(prep1 '((mtimes) -1 $%gamma)))))
426 ((= subl 0)
427 (cons '((-1 . 1) -1 . 1)
428 (if (> 0.0 npw) ()
429 `(((0 . 1)
430 . ,(prep1 '((mtimes) -1 $%gamma)))))))
431 (t (setq sf-last (factorial subl))
432 `(((,(- (1+ subl)) . 1)
433 ,(* (expt -1 (1+ subl))
434 (factorial subl)) . 1))))
435 e (if (< subl 1) (- subl) -1)
436 sf-sign (if (< subl 1) -1 (expt -1 subl)))
437 a (setq e (1+ e) sf-sign (- sf-sign))
438 (if (> e npw) (return l)
439 (rplacd (last l)
440 `(((,e . 1)
441 . ,(rctimes (rcplygam e sf-sign subl sf-last)
442 (prep1 ($zeta (+ (1+ subl) e))))))))
443 (go a))))
445 (defun rcplygam (k sf-sign subl sf-last)
446 (declare (fixnum k) )
447 (cond ((= subl -1) (cons sf-sign k))
448 ((= subl 0) (cons sf-sign 1))
449 (t (prog1
450 (cons (* sf-sign sf-last) 1)
451 (setq sf-last
452 (quot (* sf-last (+ subl (1+ k)))
453 (1+ k)))))))
455 (defun plygam-ord (subl)
456 (if (equal (car subl) -1) (ncons (rcone))
457 `((,(m- (m1+ (car subl))) . 1))))
459 (defun plygam-pole (a c func)
460 (if (rcmintegerp c)
461 (let ((ps (get-lexp (m- a (rcdisrep c)) () t)))
462 (rplacd (cddr ps) (cons `((0 . 1) . ,c) (cdddr ps)))
463 (if (atom func) (gam-const a ps func)
464 (plygam-const a ps func)))
465 (prep1 (simplifya
466 (if (atom func) `((%gamma) ,(rcdisrep c))
467 `((mqapply) ,func ,(rcdisrep c)))
468 () ))))
470 (defun gam-const (a arg func)
471 (let ((const (ps-lc* arg)) (arg-c))
472 (cond ((not (rcintegerp const))
473 (taylor2 (diff-expand `((%gamma) ,a) tlist)))
475 (setq const (car const))
476 (if (pscoefp arg) (setq arg-c (get-lexp (m+t a (- const)) (rcone) (signp le const))))
477 (if (and arg-c (not (psp arg-c)))
478 (taylor2 (simplify `((%gamma) ,const)))
479 (let ((datum (get-datum (get-key-var (gvar (or arg-c arg)))))
480 (ord (if arg-c (le (terms arg-c)) (le (n-term (terms arg))))))
481 (setq func (current-trunc datum))
482 (if (> const 0)
483 (pstimes (let-pw datum (e- func ord) (expand (m+t a (- const)) '%gamma))
484 (let-pw datum (e+ func ord)
485 (tsprsum (m+t a (m-t '%%taylor-index%%))
486 `(%%taylor-index%% 1 ,const) '%product)))
487 (pstimes (expand (m+t a (- const)) '%gamma)
488 (let-pw datum (e+ func ord)
489 (psexpt (tsprsum (m+t a '%%taylor-index%%)
490 `(%%taylor-index%% 0 ,(- (1+ const))) '%product)
491 (rcmone)))))))))))
493 (defun plygam-const (a arg func)
494 (let ((const (ps-lc* arg)) (sub (cadr func)))
495 (cond
496 ((or (not (integerp sub)) (< sub -1))
497 (tay-err "Unable to expand at a subscript in"))
498 ((not (rcintegerp const))
499 (taylor2 (diff-expand `((mqapply) ,func ,a) tlist)))
500 (t (setq const (car const))
501 (psplus
502 (expand (m+t a (- const)) func)
503 (if (> const 0)
504 (pstimes
505 (cons (* (expt -1 sub) (factorial sub)) 1)
506 (tsprsum `((mexpt) ,(m+t a (m-t '%%taylor-index%%)) ,(- (1+ sub)))
507 `(%%taylor-index%% 1 ,const) '%sum))
508 (pstimes
509 (cons (* (expt -1 (1+ sub)) (factorial sub)) 1)
510 (tsprsum `((mexpt) ,(m+t a '%%taylor-index%%) ,(- (1+ sub)))
511 `(%%taylor-index%% 0 ,(- (1+ const))) '%sum))))))))
513 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
514 ;;; Lambert W function
515 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
517 ;; References
519 ;; Corless, R. M., Gonnet, D. E. G., Jeffrey, D. J., Knuth, D. E. (1996).
520 ;; "On the Lambert W function". Advances in Computational Mathematics 5:
521 ;; pp 329-359
523 ;; http://www.apmaths.uwo.ca/~djeffrey/Offprints/W-adv-cm.pdf.
524 ;; or http://www.apmaths.uwo.ca/~rcorless/frames/PAPERS/LambertW/
526 ;; D. J. Jeffrey, D. E. G. Hare, R. M. Corless
527 ;; Unwinding the branches of the Lambert W function
528 ;; The Mathematical Scientist, 21, pp 1-7, (1996)
529 ;; http://www.apmaths.uwo.ca/~djeffrey/Offprints/wbranch.pdf
531 ;; Winitzki, S. Uniform Approximations for Transcendental Functions.
532 ;; In Part 1 of Computational Science and its Applications - ICCSA 2003,
533 ;; Lecture Notes in Computer Science, Vol. 2667, Springer-Verlag,
534 ;; Berlin, 2003, 780-789. DOI 10.1007/3-540-44839-X_82
535 ;; http://homepages.physik.uni-muenchen.de/~Winitzki/papers/
537 ;; Darko Verebic,
538 ;; Having Fun with Lambert W(x) Function
539 ;; arXiv:1003.1628v1, March 2010, http://arxiv.org/abs/1003.1628
541 ;; See also http://en.wikipedia.org/wiki/Lambert's_W_function
543 (defmfun $lambert_w (z)
544 (simplify (list '(%lambert_w) (resimplify z))))
546 ;;; Set properties to give full support to the parser and display
547 (defprop $lambert_w %lambert_w alias)
548 (defprop $lambert_w %lambert_w verb)
549 (defprop %lambert_w $lambert_w reversealias)
550 (defprop %lambert_w $lambert_w noun)
552 ;;; lambert_w is a simplifying function
553 (defprop %lambert_w simp-lambertw operators)
555 ;;; Derivative of lambert_w
556 (defprop %lambert_w
557 ((x)
558 ((mtimes)
559 ((mexpt) $%e ((mtimes ) -1 ((%lambert_w) x)))
560 ((mexpt) ((mplus) 1 ((%lambert_w) x)) -1)))
561 grad)
563 ;;; Integral of lambert_w
564 ;;; integrate(W(x),x) := x*(W(x)^2-W(x)+1)/W(x)
565 (defprop %lambert_w
566 ((x)
567 ((mtimes)
569 ((mplus)
570 ((mexpt) ((%lambert_w) x) 2)
571 ((mtimes) -1 ((%lambert_w) x))
573 ((mexpt) ((%lambert_w) x) -1)))
574 integral)
576 (defun simp-lambertw (x yy z)
577 (declare (ignore yy))
578 (oneargcheck x)
579 (setq x (simpcheck (cadr x) z))
580 (cond ((equal x 0) 0)
581 ((equal x 0.0) 0.0)
582 ((zerop1 x) ($bfloat 0)) ;bfloat case
583 ((alike1 x '$%e)
584 ;; W(%e) = 1
586 ((alike1 x '((mtimes simp) ((rat simp) -1 2) ((%log simp) 2)))
587 ;; W(-log(2)/2) = -log(2)
588 '((mtimes simp) -1 ((%log simp) 2)))
589 ((alike1 x '((mtimes simp) -1 ((mexpt simp) $%e -1)))
590 ;; W(-1/e) = -1
592 ((alike1 x '((mtimes) ((rat) -1 2) $%pi))
593 ;; W(-%pi/2) = %i*%pi/2
594 '((mtimes simp) ((rat simp) 1 2) $%i $%pi))
595 ;; numerical evaluation
596 ((complex-float-numerical-eval-p x)
597 ;; x may be an integer. eg "lambert_w(1),numer;"
598 (if (integerp x)
599 (to (bigfloat::lambert-w-k 0 (bigfloat:to ($float x))))
600 (to (bigfloat::lambert-w-k 0 (bigfloat:to x)))))
601 ((complex-bigfloat-numerical-eval-p x)
602 (to (bigfloat::lambert-w-k 0 (bigfloat:to x))))
603 (t (list '(%lambert_w simp) x))))
605 ;; An approximation of the k-branch of generalized Lambert W function
606 ;; k integer
607 ;; z real or complex lisp float
608 ;; Used as initial guess for Halley's iteration.
609 ;; When W(z) is real, ensure that guess is real.
610 (defun init-lambert-w-k (k z)
611 (let ( ; parameters for k = +/- 1 near branch pont z=-1/%e
612 (branch-eps 0.2e0)
613 (branch-point (/ -1 %e-val))) ; branch pont z=-1/%e
614 (cond
615 ; For principal branch k=0, use expression by Winitzki
616 ((= k 0) (init-lambert-w-0 z))
617 ; For k=1 branch, near branch point z=-1/%e with im(z) < 0
618 ((and (= k 1)
619 (< (imagpart z) 0)
620 (< (abs (- branch-point z)) branch-eps))
621 (bigfloat::lambert-branch-approx z))
622 ; For k=-1 branch, z real with -1/%e < z < 0
623 ; W(z) is real in this range
624 ((and (= k -1) (realp z) (> z branch-point) (< z 0))
625 (init-lambert-w-minus1 z))
626 ; For k=-1 branch, near branch point z=-1/%e with im(z) >= 0
627 ((and (= k -1)
628 (>= (imagpart z) 0)
629 (< (abs (- branch-point z)) branch-eps))
630 (bigfloat::lambert-branch-approx z))
631 ; Default to asymptotic expansion Corless et al (4.20)
632 ; W_k(z) = log(z) + 2.pi.i.k - log(log(z)+2.pi.i.k)
633 (t (let ((two-pi-i-k (complex 0.0e0 (* 2 pi k))))
634 (+ (log z)
635 two-pi-i-k
636 (* -1 (log (+ (log z) two-pi-i-k )))))))))
638 ;; Complex value of the principal branch of Lambert's W function in
639 ;; the entire complex plane with relative error less than 1%, given
640 ;; standard branch cuts for sqrt(z) and log(z).
641 ;; Winitzki (2003)
642 (defun init-lambert-w-0 (z)
643 (let ((a 2.344e0) (b 0.8842e0) (c 0.9294e0) (d 0.5106e0) (e -1.213e0)
644 (y (sqrt (+ (* 2 %e-val z ) 2)) ) ) ; y=sqrt(2*%e*z+2)
645 ; w = (2*log(1+B*y)-log(1+C*log(1+D*y))+E)/(1+1/(2*log(1+B*y)+2*A)
647 (+ (* 2 (log (+ 1 (* b y))))
648 (* -1 (log (+ 1 (* c (log (+ 1 (* d y)))))))
650 (+ 1
651 (/ 1 (+ (* 2 (log (+ 1 (* b y)))) (* 2 a)))))))
653 ;; Approximate k=-1 branch of Lambert's W function over -1/e < z < 0.
654 ;; W(z) is real, so we ensure the starting guess for Halley iteration
655 ;; is also real.
656 ;; Verebic (2010)
657 (defun init-lambert-w-minus1 (z)
658 (cond
659 ((not (realp z))
660 (merror "z not real in init-lambert-w-minus1"))
661 ((or (< z (/ -1 %e-val)) (plusp z))
662 (merror "z outside range of approximation in init-lambert-w-minus1"))
663 ;; In the region where W(z) is real
664 ;; -1/e < z < C, use power series about branch point -1/e ~ -0.36787
665 ;; C = -0.3 seems a reasonable crossover point
666 ((< z -0.3)
667 (bigfloat::lambert-branch-approx z))
668 ;; otherwise C <= z < 0, use iteration W(z) ~ ln(-z)-ln(-W(z))
669 ;; nine iterations are sufficient over -0.3 <= z < 0
670 (t (let* ((ln-z (log (- z))) (maxiter 9) (w ln-z))
671 (dotimes (k maxiter w)
672 (setq w (- ln-z (log (- w)))))))))
674 (in-package #:bigfloat)
676 ;; Approximate Lambert W(k,z) for k=1 and k=-1 near branch point z=-1/%e
677 ;; using power series in y=-sqrt(2*%e*z+2)
678 ;; for im(z) < 0, approximates k=1 branch
679 ;; for im(z) >= 0, approximates k=-1 branch
681 ;; Corless et al (1996) (4.22)
682 ;; Verebic (2010)
684 ;; z is a real or complex bigfloat:
685 (defun lambert-branch-approx (z)
686 (let ((y (- (sqrt (+ (* 2 (%e z) z ) 2)))) ; y=-sqrt(2*%e*z+2)
687 (b0 -1) (b1 1) (b2 -1/3) (b3 11/72))
688 (+ b0 (* y (+ b1 (* y (+ b2 (* b3 y))))))))
690 ;; Algorithm based in part on Corless et al (1996).
692 ;; It is Halley's iteration applied to w*exp(w).
695 ;; w[j] exp(w[j]) - z
696 ;; w[j+1] = w[j] - -------------------------------------------------
697 ;; (w[j]+2)(w[j] exp(w[j]) -z)
698 ;; exp(w[j])(w[j]+1) - ---------------------------
699 ;; 2 w[j] + 2
701 ;; The algorithm has cubic convergence. Once convergence begins, the
702 ;; number of digits correct at step k is roughly 3 times the number
703 ;; which were correct at step k-1.
705 ;; Convergence can stall using convergence test abs(w[j+1]-w[j]) < prec,
706 ;; as happens for generalized_lambert_w(-1,z) near branch point z = -1/%e
707 ;; Therefore also stop iterating if abs(w[j]*exp(w[j]) - z) << abs(z)
708 (defun lambert-w-k (k z &key (maxiter 50))
709 (let ((w (init-lambert-w-k k z)) we w1e delta (prec (* 4 (epsilon z))))
710 (dotimes (i maxiter (maxima::merror "lambert-w-k did not converge"))
711 (setq we (* w (exp w)))
712 (when (<= (abs (- z we)) (* 4 (epsilon z) (abs z))) (return))
713 (setq w1e (* (1+ w) (exp w)))
714 (setq delta (/ (- we z)
715 (- w1e (/ (* (+ w 2) (- we z)) (+ 2 (* 2 w))))))
716 (decf w delta)
717 (when (<= (abs (/ delta w)) prec) (return)))
718 ;; Check iteration converged to correct branch
719 (check-lambert-w-k k w z)
722 (defmethod init-lambert-w-k ((k integer) (z number))
723 (maxima::init-lambert-w-k k z))
725 (defmethod init-lambert-w-k ((k integer) (z bigfloat))
726 (bfloat-init-lambert-w-k k z))
728 (defmethod init-lambert-w-k ((k integer) (z complex-bigfloat))
729 (bfloat-init-lambert-w-k k z))
731 (defun bfloat-init-lambert-w-k (k z)
732 "Approximate generalized_lambert_w(k,z) for bigfloat: z as initial guess"
733 (let ((branch-point -0.36787944117144)) ; branch point -1/%e
734 (cond
735 ;; if k=-1, z very close to -1/%e and imag(z)>=0, use power series
736 ((and (= k -1)
737 (or (zerop (imagpart z))
738 (plusp (imagpart z)))
739 (< (abs (- z branch-point)) 1e-10))
740 (lambert-branch-approx z))
741 ;; if k=1, z very close to -1/%e and imag(z)<0, use power series
742 ((and (= k 1)
743 (minusp (imagpart z))
744 (< (abs (- z branch-point)) 1e-10))
745 (lambert-branch-approx z))
746 ;; Initialize using float value if z is representable as a float
747 ((< (abs z) 1.0e100)
748 (if (complexp z)
749 (bigfloat (lambert-w-k k (cl:complex (float (realpart z) 1.0)
750 (float (imagpart z) 1.0))))
751 (bigfloat (lambert-w-k k (float z 1.0)))))
752 ;; For large z, use Corless et al (4.20)
753 ;; W_k(z) ~ log(z) + 2.pi.i.k - log(log(z)+2.pi.i.k)
755 (let ((log-z (log z)))
756 (if (= k 0)
757 (- log-z (log log-z))
758 (let* ((i (make-instance 'complex-bigfloat :imag (intofp 1)))
759 (two-pi-i-k (* 2 (%pi z) i k)))
760 (- (+ log-z two-pi-i-k)
761 (log (+ log-z two-pi-i-k))))))))))
763 ;; Check Lambert W iteration converged to the correct branch
764 ;; W_k(z) + ln W_k(z) = ln z, for k = -1 and z in [-1/e,0)
765 ;; = ln z + 2 pi i k, otherwise
766 ;; See Jeffrey, Hare, Corless (1996), eq (12)
767 ;; k integer
768 ;; z, w bigfloat: numbers
769 (defun check-lambert-w-k (k w z)
770 (let ((tolerance 1.0e-6))
772 (cond
773 ;; k=-1 branch with z and w real.
774 ((and (= k -1) (realp z) (minusp z) (>= z (/ -1 (%e z))))
775 (if (and (realp w)
776 (<= w -1)
777 (< (abs (+ w (log w) (- (log z)))) tolerance))
779 nil))
781 ; i k = (W_k(z) + ln W_k(z) - ln(z)) / 2 pi
782 (let (ik)
783 (setq ik (/ (+ w (log w) (- (log z))) (* 2 (%pi z))))
784 (if (and (< (realpart ik) tolerance)
785 (< (abs (- k (imagpart ik))) tolerance))
787 nil))))
789 (maxima::merror "Lambert W iteration converged to wrong branch"))))
791 (in-package :maxima)
793 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
794 ;;; Generalized Lambert W function
795 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
797 (defmfun $generalized_lambert_w (k z)
798 (simplify (list '(%generalized_lambert_w) (resimplify k) (resimplify z))))
800 ;;; Set properties to give full support to the parser and display
801 (defprop $generalized_lambert_w %generalized_lambert_w alias)
802 (defprop $generalized_lambert_w %generalized_lambert_w verb)
803 (defprop %generalized_lambert_w $generalized_lambert_w reversealias)
804 (defprop %generalized_lambert_w $generalized_lambert_w noun)
806 ;;; lambert_w is a simplifying function
807 (defprop %generalized_lambert_w simp-generalized-lambertw operators)
809 ;;; Derivative of lambert_w
810 (defprop %generalized_lambert_w
811 ((k x)
813 ((mtimes)
814 ((mexpt) $%e ((mtimes ) -1 ((%generalized_lambert_w) k x)))
815 ((mexpt) ((mplus) 1 ((%generalized_lambert_w) k x)) -1)))
816 grad)
818 ;;; Integral of lambert_w
819 ;;; integrate(W(k,x),x) := x*(W(k,x)^2-W(k,x)+1)/W(k,x)
820 (defprop %generalized_lambert_w
821 ((k x)
823 ((mtimes)
825 ((mplus)
826 ((mexpt) ((%generalized_lambert_w) k x) 2)
827 ((mtimes) -1 ((%generalized_lambert_w) k x))
829 ((mexpt) ((%generalized_lambert_w) k x) -1)))
830 integral)
832 (defun simp-generalized-lambertw (expr ignored z)
833 (declare (ignore ignored))
834 (twoargcheck expr)
835 (let ((k (simpcheck (cadr expr) z))
836 (x (simpcheck (caddr expr) z)))
837 (cond
838 ;; Numerical evaluation for real or complex x
839 ((and (integerp k) (complex-float-numerical-eval-p x))
840 ;; x may be an integer. eg "generalized_lambert_w(0,1),numer;"
841 (if (integerp x)
842 (to (bigfloat::lambert-w-k k (bigfloat:to ($float x))))
843 (to (bigfloat::lambert-w-k k (bigfloat:to x)))))
844 ;; Numerical evaluation for real or complex bigfloat x
845 ((and (integerp k) (complex-bigfloat-numerical-eval-p x))
846 (to (bigfloat::lambert-w-k k (bigfloat:to x))))
847 (t (list '(%generalized_lambert_w simp) k x)))))
849 (in-package "BIGFLOAT")
851 (defvar *debug-li-eval* nil)
853 (defun li-using-powers-of-log (n x)
854 ;; When x is on or near the unit circle the other
855 ;; approaches don't work. Use the expansion in powers of
856 ;; log(z) (from cephes cpolylog)
858 ;; li[s](z) = sum(Z(s-j)*(log(z))^j/j!, j = 0, inf)
860 ;; where Z(j) = zeta(j) for j != 1. For j = 1:
862 ;; Z(1) = -log(-log(z)) + sum(1/k, k, 1, s - 1)
865 ;; This is similar to
866 ;; http://functions.wolfram.com/10.08.06.0024.01, but that
867 ;; identity is clearly undefined if v is a positive
868 ;; integer because zeta(1) is undefined.
870 ;; Thus,
872 ;; li[3](z) = Z(3) + Z(2)*log(z) + Z(1)*log(z)^2/2!
873 ;; + Z(0)*log(z)^3/3! + sum(Z(-k)*log(z)^(k+4)/(k+4)!,k,1,inf);
875 ;; But Z(-k) = zeta(-k) is 0 if k is even. So
877 ;; li[3](z) = Z(3) + Z(2)*log(z) + Z(1)*log(z)^2/2!
878 ;; + Z(0)*log(z)^3/3! + sum(Z(-(2*k+1))*log(z)^(2*k+4)/(2*k+4)!,k,0,inf);
879 (flet ((zfun (j)
880 (cond ((= j 1)
881 (let ((sum (- (log (- (log x))))))
882 (+ sum
883 (loop for k from 1 below n
884 sum (/ k)))))
886 (to (maxima::$zeta (maxima::to (float j (realpart x)))))))))
887 (let* ((eps (epsilon x))
888 (logx (log x))
889 (sum 0))
890 (do* ((k 0 (1+ k))
891 (top 1 (* top logx))
892 (bot 1 (* bot k))
893 (term (* (/ top bot) (to (zfun (- n k))))
894 (* (/ top bot) (to (zfun (- n k))))))
895 ((and (> k 4)
896 (oddp (- n k)) ;; even terms are all zero
897 (<= (abs term) (* (abs sum) eps))))
898 ;;(format t "~3d: ~A / ~A = ~A~%" k top bot term)
899 (incf sum term))
900 sum)))
903 (defun li3numer (x)
904 ;; If |x| < series-threshold, the series is used.
905 (let ((series-threshold 0.8))
906 (cond ((zerop x)
907 0.0)
908 ((= x 1)
909 (maxima::$zeta (maxima::to (float 3 x))))
910 ((= x -1)
911 ;; li[3](-1) = -(1-2^(1-3))*li[3](1)
912 ;; = -3/4*zeta(3)
914 ;; From the formula
916 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
918 ;; (See http://functions.wolfram.com/10.08.03.0003.01)
919 (* -3/4 (to (maxima::$zeta (maxima::to (float 3 x))))))
920 ((> (abs x) 1)
921 ;; For z not in the interval (0, 1) and for integral n, we
922 ;; have the identity:
924 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
925 ;; + 2 * sum(li[2*r](-1)/(n-2*r)!*log(-z)^(n-2*r), r, 1, floor(n/2))
927 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog/17/02/01/01/0008/)
929 ;; In particular for n = 3:
931 ;; li[3](z) = li[3](1/z) - log(-z)/6*(log(-z)^2+%pi^2)
932 (let* ((lg (log (- x)))
933 (dpi (%pi x))
934 (result (- (li3numer (/ x))
935 (* (/ lg 6)
936 (+ (* lg lg) (* dpi dpi))))))
937 result))
938 ((> (abs x) series-threshold)
939 (let ((result (li-using-powers-of-log 3 x)))
940 ;; For real x, li-using-power-of-log can return a complex
941 ;; number with a tiny imaginary part. Get rid of that
942 ;; when x is real.
943 (if (realp x)
944 (realpart result)
945 result)))
946 ;; Don't use the identity below because the identity seems
947 ;; to be incorrect. For example, for x = -0.862 it returns
948 ;; a complex value with an imaginary part that is not close
949 ;; to zero as expected.
950 #+nil
951 ((> (abs x) series-threshold)
952 ;; The series converges too slowly so use the identity:
954 ;; li[3](-x/(1-x)) + li[3](1-x) + li[3](x)
955 ;; = li[3](1) + %pi^2/6*log(1-x) - 1/2*log(x)*(log(1-x))^2 + 1/6*(log(1-x))^3
957 ;; Or
959 ;; li[3](x) = li[3](1) + %pi^2/6*log(1-x) - 1/2*log(x)*(log(1-x))^2 + 1/6*(log(1-x))^3
960 ;; - li[3](-x/(1-x)) - li[3](1-x)
962 ;; (See http://functions.wolfram.com/10.08.17.0048.01)
963 (let* ((dpi (%pi x))
964 (u (log x))
965 (s (/ (* u u u) 6))
966 (xc (- 1 x)))
967 (decf s (* 0.5 u u (log xc)))
968 (incf s (/ (* dpi dpi u) 6))
969 (decf s (li3numer (- (/ xc x))))
970 (decf s (li3numer xc))
971 (incf s (li3numer (float 1 x)))))
973 ;; Sum the power series. threshold determines when the
974 ;; summation has converted.
975 (let* ((threshold (epsilon x))
976 (p (* x x x))
977 (term (/ p 27)))
978 (incf term (* 0.125 x x))
979 (incf term x)
980 (do* ((k 4 (1+ k))
981 (p1 (* p x) (* p1 x))
982 (h (/ p1 (* k k k)) (/ p1 (* k k k)))
983 (s h (+ s h)))
984 ((<= (abs (/ h s)) threshold)
985 (+ s term))))))))
987 (defun li2numer (z)
988 ;; The series threshold to above sqrt(1/2) because li[2](%i) needs
989 ;; the value of li[2](1/2-%i/2), and the magnitude of the argument
990 ;; is sqrt(1/2) = 0.707. If the threshold is below this, we get
991 ;; into an infinite recursion oscillating between the two args.
992 (let ((series-threshold .75))
993 (cond ((zerop z)
995 ((= z 1)
996 ;; %pi^2/6. This follows from the series.
997 (/ (expt (%pi z) 2) 6))
998 ((= z -1)
999 ;; -%pi^2/12. From the formula
1001 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
1003 ;; (See http://functions.wolfram.com/10.08.03.0003.01)
1004 (/ (expt (%pi z) 2) -12))
1005 ((> (abs z) 1)
1006 ;; Use
1007 ;; li[2](z) = -li[2](1/z) - 1/2*log(-z)^2 - %pi^2/6,
1009 ;; valid for all z not in the intervale (0, 1).
1011 ;; (See http://functions.wolfram.com/10.08.17.0013.01)
1012 (- (+ (li2numer (/ z))
1013 (* 0.5 (expt (log (- z)) 2))
1014 (/ (expt (%pi z) 2) 6))))
1015 ;; this converges faster when close to unit circle
1016 ((> (abs z) series-threshold)
1017 (li-using-powers-of-log 2 z))
1018 ;; no longer in use
1019 ;;(> (abs z) series-threshold)
1020 ;; For 0.5 <= |z|, where the series would not converge very quickly, use
1022 ;; li[2](z) = li[2](1/(1-z)) + 1/2*log(1-z)^2 - log(-z)*log(1-z) - %pi^2/6
1024 ;; (See http://functions.wolfram.com/10.08.17.0016.01)
1025 ;;(let* ((1-z (- 1 z))
1026 ;; (ln (log 1-z)))
1027 ;; (+ (li2numer (/ 1-z))
1028 ;; (* 0.5 ln ln)
1029 ;; (- (* (log (- z))
1030 ;; ln))
1031 ;; (- (/ (expt (%pi z) 2) 6)))))
1033 ;; Series evaluation:
1035 ;; li[2](z) = sum(z^k/k^2, k, 1, inf);
1036 (let ((eps (epsilon z)))
1037 (do* ((k 0 (1+ k))
1038 (term z (* term (/ (* z k k)
1039 (expt (1+ k) 2))))
1040 (sum z (+ term sum)))
1041 ((<= (abs (/ term sum)) eps)
1042 sum)))))))
1044 (defun polylog-power-series (s z)
1045 ;; Series evaluation:
1047 ;; li[s](z) = sum(z^k/k^s, k, 1, inf);
1048 (let ((eps (epsilon z)))
1049 (do* ((k 1 (1+ k))
1050 (term z (* term z (expt (/ (- k 1) k) s)))
1051 (sum z (+ term sum)))
1052 ((<= (abs (/ term sum)) eps)
1053 ;; Return the value and the number of terms used, for
1054 ;; debugging and for helping in determining the series
1055 ;; threshold.
1056 (values sum k)))))
1058 (defun polylog-log-series (s z)
1059 ;; When x is on or near the unit circle the other
1060 ;; approaches don't work. Use the expansion in powers of
1061 ;; log(z) (from cephes cpolylog)
1063 ;; li[s](z) = sum(Z(s-j)*(log(z))^j/j!, j = 0, inf)
1065 ;; where Z(j) = zeta(j) for j != 1. For j = 1:
1067 ;; Z(1) = -log(-log(z)) + sum(1/k, k, 1, s - 1)
1068 (flet ((zfun (j)
1069 ;; Compute Z(j)
1070 (cond ((= j 1)
1071 (let ((sum (- (log (- (log z))))))
1072 (+ sum
1073 (loop for k from 1 below s
1074 sum (/ k)))))
1076 (to (maxima::$zeta (maxima::to (float j (realpart z)))))))))
1077 (let* ((eps (epsilon z))
1078 (logx (log z))
1079 (logx^2 (* logx logx))
1080 (top logx)
1081 (bot 1)
1082 (sum (zfun s)))
1083 ;; Compute sum(Z(s-j)*log(z)^j/j!, j = 1, s)
1084 (do* ((k 1 (1+ k))
1085 (zf (zfun (- s k)) (zfun (- s k)))
1086 (term (* (/ top bot) zf)
1087 (* (/ top bot) zf)))
1088 ((> k s))
1089 (when *debug-li-eval*
1090 (format t "~3d: ~A / ~A * ~A => ~A~%" k top bot zf term))
1091 (incf sum term)
1092 (setf bot (* bot (1+ k)))
1093 (setf top (* top logx)))
1095 (when *debug-li-eval*
1096 (format t "s = ~A, sum = ~S top, bot = ~S ~S~%"
1097 s sum top bot))
1098 ;; Compute the sum for j = s+1 and up. But since
1099 ;; zeta(-k) is 0 for k even, we only every other term.
1100 (do* ((k (+ s 1) (+ k 2))
1101 (zf (zfun (- s k)) (zfun (- s k)))
1102 (term (* (/ top bot) zf)
1103 (* (/ top bot) zf)))
1104 ((<= (abs term) (* (abs sum) eps))
1105 ;; Return the result and the number of terms used for
1106 ;; helping in determining the series threshold and the
1107 ;; log-series threshold.
1109 ;; Note that if z is real and less than 0, li[s](z) is
1110 ;; real. The series can return a tiny complex value in
1111 ;; this case, so we want to clear that out before
1112 ;; returning the answer.
1113 (values (if (and (realp z) (minusp z))
1114 (realpart sum)
1115 sum)
1117 (when *debug-li-eval*
1118 (format t "~3d: ~A / ~A = ~A~%" k top bot term))
1119 (incf sum term)
1120 (setf bot (* bot (+ k 1) (+ k 2)))
1121 (setf top (* top logx^2))))))
1123 (defun polylog-inversion-formula (s z)
1124 ;; For z not in the interval (0, 1) and for integral n, we
1125 ;; have the identity:
1127 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1128 ;; + 2 * sum(li[2*r](-1)/(n-2*r)!*log(-z)^(n-2*r), r, 1, floor(n/2))
1130 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog/17/02/01/01/0008/)
1132 ;; Or
1134 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1135 ;; + 2 * sum(li[2*m-2*r](-1)/(n-2*m+2*r)!*log(-z)^(n-2*m+2*r), r, 0, m - 1)
1137 ;; where m = floor(n/2). Thus, n-2*m = 0 if n is even and 1 if n is odd.
1139 ;; For n = 2*m, we have
1141 ;; li[2*m](z) = -log(-z)^(2*m)/(2*m)! - li[2*m](1/z)
1142 ;; + 2 * sum(li[2*r](-1)/(2*m-2*r)!*log(-z)^(2*m-2*r), r, 1, m)
1143 ;; = -log(-z)^(2*m)/(2*m)! - li[2*m](1/z)
1144 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,m-1);
1146 ;; For n = 2*m+1, we have
1148 ;; li[2*m+1](z) = -log(-z)^(2*m+1)/(2*m+1)! + li[2*m+1](1/z)
1149 ;; + 2 * sum(li[2*r](-1)/(2*m-2*r + 1)!*log(-z)^(2*m-2*r + 1), r, 1, m)
1150 ;; = -log(-z)^(2*m+1)/(2*m+1)! + li[2*m+1](1/z)
1151 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,m-1);
1152 ;; Thus,
1154 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1155 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,floor(n/2)-1);
1156 (let* ((lgz (log (- z)))
1157 (lgz^2 (* lgz lgz))
1158 (half-s (floor s 2))
1159 (neg-1 (float -1 (realpart z)))
1160 (sum 0))
1161 (if (evenp s)
1162 (do* ((r 0 (1+ r))
1163 (top (if (oddp s) lgz 1) (* top lgz^2))
1164 (bot 1 (* bot (+ r r -1) (+ r r)))
1165 (term (* (li-s-simp (* 2 (- half-s r)) neg-1)
1166 (/ top bot))
1167 (* (li-s-simp (* 2 (- half-s r)) neg-1)
1168 (/ top bot))))
1169 ((>= r half-s))
1170 (incf sum term)
1171 (when *debug-li-eval*
1172 (format t "r = ~4d: ~A / ~A, ~A; ~A~%" r top bot term sum)))
1173 (do* ((r 0 (1+ r))
1174 (top (if (oddp s) lgz 1) (* top lgz^2))
1175 (bot 1 (* bot (+ r r) (+ r r 1)))
1176 (term (* (li-s-simp (* 2 (- half-s r)) neg-1)
1177 (/ top bot))
1178 (* (li-s-simp (* 2 (- half-s r)) neg-1)
1179 (/ top bot))))
1180 ((>= r half-s))
1181 (incf sum term)
1182 (when *debug-li-eval*
1183 (format t "r = ~4d: ~A / ~A, ~A; ~A~%" r top bot term sum))))
1184 (+ (+ sum sum)
1185 (- (/ (expt lgz s)
1186 (maxima::take '(maxima::mfactorial) s)))
1187 (* (expt -1 (- s 1))
1188 (li-s-simp s (/ z))))))
1190 (defun li-s-simp (s z)
1191 (let ((series-threshold 0.5)
1192 (log-series-threshold 2))
1193 (cond ((zerop z)
1194 (maxima::to (to 0.0)))
1195 ((= z 1)
1196 (maxima::$zeta (maxima::to (float s z))))
1197 ((= z -1)
1198 (- (* (- 1 (expt 2 (- 1 s)))
1199 (to (li-s-simp s (- z))))))
1200 ((<= (abs z) series-threshold)
1201 (values (polylog-power-series s z)))
1202 ((<= (abs z) log-series-threshold)
1203 (values (polylog-log-series s z)))
1204 ((> (abs z) 1.5)
1205 (polylog-inversion-formula s z)))))