1 ;;; -*- Mode: Lisp; Package: Maxima; Syntax: Common-Lisp; Base: 10 -*- ;;;;
2 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3 ;;; The data in this file contains enhancements. ;;;;;
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 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
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
))
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
44 ((not (integerp s
)) ())
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
))
55 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
56 (mul (add -
1 (inv (expt 2 (- s
1))))
57 (take '(%zeta
) s
))))))
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
)))
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))))
76 (to (bigfloat::li-s-simp s
(bigfloat:to a
)))))))
77 (eqtest (subfunmakes '$li
(ncons s
) (ncons a
))
80 ;; Expand the Polylogarithm li[s](z) for a negative integer parameter s.
81 (defun lisimp-negative-integer (s z
)
83 (mul (inv (power (sub 1 z
) (+ n
1)))
84 (let ((index1 (gensumindex))
88 (let ((index2 (gensumindex)))
90 (mul (power -
1 (add index2
1))
91 (take '(%binomial
) (+ n
1) (sub index2
1))
92 (power (add 1 (sub index1 index2
)) n
))
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))))))
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
127 (declare (ignore subl
))
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
137 (setq npw
(/ (float (car pw
)) (float (cdr pw
))))
139 l
(cons '((0 .
1) 0 .
1)
142 (if (> e npw
) (return l
)
145 .
,(prep1 (m^ e
(m- subl
)))))))
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
157 (defun li-asymptotic-expansion (pw s z
)
158 (m+l
(loop for k from
0 to pw collect
160 (m- 1 (m^
2 (m- 1 (m* 2 k
))))
161 (m^
(m* 2 '$%pi
) (m* 2 k
))
162 (m// ($bern
(m* 2 k
))
163 `((mfactorial) ,(m* 2 k
)))
164 (m// (m^
`((%log
) ,(m- z
)) (m- 2 (m* 2 k
)))
165 ($gamma
(m+ s
1 (m* -
2 k
))))))))
167 ;; Numerical evaluation for Chebyschev expansions of the first kind
169 (defun cheby (x chebarr
)
170 (let ((bn+2 0.0) (bn+1 0.0))
171 (do ((i (floor (aref chebarr
0)) (1- i
)))
172 ((< i
1) (- bn
+1 (* bn
+2 x
)))
174 (prog1 bn
+1 (setq bn
+1 (+ (aref chebarr i
)
175 (- (* 2.0 x bn
+1) bn
+2))))))))
177 (defun cheby-prime (x chebarr
)
179 (* (aref chebarr
1) 0.5)))
181 ;; These should really be calculated with minimax rational approximations.
182 ;; Someone has done LI[2] already, and this should be updated; I haven't
183 ;; seen any results for LI[3] yet.
186 ;; Spence's function can be used to compute li[2] for 0 <= x <= 1.
187 ;; To compute the rest, we need the following identities:
189 ;; li[2](x) = -li[2](1/x)-log(-x)^2/2-%pi^2/6
190 ;; li[2](x) = li[2](1/(1-x)) + log(1-x)*log((1-x)/x^2)/2 - %pi^2/6
192 ;; The first tells us how to compute li[2] for x > 1. The result is complex.
193 ;; For x < 0, the second can be used, and the result is real.
195 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog2/17/01/01/)
199 (* 0.5 (log (- 1 x
)) (log (/ (- 1 x
) (* x x
))))
200 (- (/ (cl:expt
(float pi
) 2) 6))))
204 (/ (cl:expt
(float pi
) 2) 6))
206 ;; li[2](x) = -li[2](1/x)-log(-x)^2/2-%pi^2/6
208 (/ (cl:expt
(cl:log
(- x
)) 2) 2)
209 (/ (cl:expt
(float pi
) 2) 6)))))))
210 (complexify (li2 y
))))
213 (defvar *li2
* (make-array 15.
:initial-contents
'(14.0
1.93506430 .166073033 2.48793229e-2
214 4.68636196e-3 1.0016275e-3 2.32002196e-4
215 5.68178227e-5 1.44963006e-5 3.81632946e-6
216 1.02990426e-6 2.83575385e-7 7.9387055e-8
217 2.2536705e-8 6.474338e-9)
218 :element-type
'flonum
))
221 (defvar *li3
* (make-array 15.
:initial-contents
'(14.0
1.95841721 8.51881315e-2 8.55985222e-3
222 1.21177214e-3 2.07227685e-4 3.99695869e-5
223 8.38064066e-6 1.86848945e-6 4.36660867e-7
224 1.05917334e-7 2.6478920e-8 6.787e-9
225 1.776536e-9 4.73417e-10)
226 :element-type
'flonum
))
228 (defvar *s12
* (make-array 18.
:initial-contents
'(17.0
1.90361778 .431311318 .100022507
229 2.44241560e-2 6.22512464e-3 1.64078831e-3
230 4.44079203e-4 1.22774942e-4 3.45398128e-5
231 9.85869565e-6 2.84856995e-6 8.31708473e-7
232 2.45039499e-7 7.2764962e-8 2.1758023e-8 6.546158e-9
234 :element-type
'flonum
))
237 (* x
(cheby-prime (/ (1+ (* x
4)) 3) *li2
*)))
240 (* x
(cheby-prime (/ (1+ (* 4 x
)) 3) *li3
*)))
244 (cheby-prime (/ (1+ (* 4 x
)) 3) *s12
*)))
246 ;; subtitle polygamma routines
248 ;; gross efficiency hack, exp is a function of *k*, *k* should be mbind'ed
250 (defun msum (exp lo hi
)
254 (do ((*k
* lo
(1+ *k
*)))
256 (declare (special *k
*))
257 (setq sum
(add2 sum
(meval exp
)))))))
260 (defun pole-err (exp)
261 (cond (errorsw (throw 'errorsw t
))
262 (t (merror (intl:gettext
"Pole encountered in: ~M") exp
))))
265 (declare-top (special $maxpsiposint $maxpsinegint $maxpsifracnum $maxpsifracdenom
))
267 (defprop $psi psisimp specsimp
)
269 ;; Integral of psi function psi[n](x)
275 ((and ($integerp n
) (>= n
0))
277 ((= n
0) `((%log_gamma
) ,x
))
278 (t `((mqapply) (($psi array
) ((mplus) -
1 ,n
)) ,x
))))
282 (mapcar #'(lambda (var val
)
283 (and (not (boundp var
)) (setf (symbol-value var
) val
)))
284 '($maxpsiposint $maxpsinegint $maxpsifracnum $maxpsifracdenom
)
287 (defun psisimp (expr a z
)
288 (let ((s (simpcheck (car (subfunsubs expr
)) z
)))
289 (subargcheck expr
1 1 '$psi
)
290 (setq a
(simpcheck (car (subfunargs expr
)) z
))
291 (and (setq z
(integer-representation-p a
))
294 (eqtest (psisimp1 s a
) expr
)))
296 ;; This gets pretty hairy now.
298 (defun psisimp1 (s a
)
300 (declare (special *k
*))
302 (and (integerp s
) (>= s
0) (mnumericalp a
)
303 (let (($float2bf t
)) ($float
(mfuncall '$bfpsi s a
18))))
304 (and (integerp s
) (>= s
0) ($bfloatp a
)
305 (mfuncall '$bfpsi s a $fpprec
))
306 (and (not $numer
) (not $float
) (integerp s
) (> s -
1)
309 (and (not (> a $maxpsiposint
)) ; integer values
310 (m*t
(expt -
1 s
) (factorial s
)
311 (m- (msum (inv (m^t
'*k
* (1+ s
))) 1 (1- a
))
312 (cond ((zerop s
) '$%gamma
)
313 (($zeta
(1+ s
))))))))
314 ((or (not (ratnump a
)) (ratgreaterp a $maxpsiposint
)) ())
318 (let* ((int ($entier a
)) ; reduction to fractional values
322 (if (> int $maxpsiposint
)
323 (subfunmakes '$psi
(ncons s
) (ncons int
))
324 (m*t
(expt -
1 s
) (factorial s
)
325 (msum (m^t
(m+t
(m-t a int
) '*k
*)
329 (let ((p (cadr a
)) (q (caddr a
)))
331 ((or (> p $maxpsifracnum
)
332 (> q $maxpsifracdenom
) (bignump p
) (bignump q
)) ())
335 (m+ (m* -
2 '((%log
) 2)) (m- '$%gamma
)))
337 (m+ (m* '((rat simp
) -
1 2)
338 (m^t
3 '((rat simp
) -
1 2)) '$%pi
)
339 (m* '((rat simp
) -
3 2) '((%log
) 3))
342 (m+ (m* '((rat simp
) -
1 2) '$%pi
)
343 (m* -
3 '((%log
) 2)) (m- '$%gamma
)))
345 (m- (m+ (m* '((rat simp
) 3 2) '((%log
) 3))
347 (m* '((rat simp
) 1 2) '$%pi
348 (m^t
3 '((rat simp
) 1 2)))
350 ((and (= p
2) (= q
3))
351 (m+ (m* '((rat simp
) 1 2)
352 (m^t
3 '((rat simp
) -
1 2)) '$%pi
)
353 (m* '((rat simp
) -
3 2) '((%log
) 3))
355 ((and (= p
3) (= q
4))
356 (m+ (m* '((rat simp
) 1 2) '$%pi
)
357 (m* -
3 '((%log
) 2)) (m- '$%gamma
)))
358 ((and (= p
5) (= q
6))
359 (m- (m* '((rat simp
) 1 2) '$%pi
360 (m^t
3 '((rat simp
) 1 2)))
361 (m+ (m* '((rat simp
) 3 2) '((%log
) 3))
365 ((let ((f (m* `((%cos
) ,(m* 2 a
'$%pi
'*k
*))
366 `((%log
) ,(m-t 2 (m* 2 `((%cos
)
367 ,(m//t
(m* 2 '$%pi
'*k
*)
369 (m+t
(msum f
1 (1- (truncate q
2)))
370 (let ((*k
* (truncate q
2)))
371 (declare (special *k
*))
374 ('((rat simp
) 1 2)))))
375 (m-t (m+ (m* '$%pi
'((rat simp
) 1 2)
376 `((%cot
) ((mtimes simp
) ,a $%pi
)))
379 ((alike1 a
'((rat) 1 2))
380 (m*t
(expt -
1 (1+ s
)) (factorial s
)
381 (1- (expt 2 (1+ s
))) (simplify ($zeta
(1+ s
)))))
382 ((and (ratgreaterp a
'((rat) 1 2))
386 (m+t
(psisimp1 s
(m- 1 a
))
388 ($diff
`((%cot
) ,(m* '$%pi
'$z
)) '$z s
)))
390 (declare (special $z
))
392 ((ratgreaterp a $maxpsinegint
) ;;; Reflection Formula
395 (m+t
(m+t
(psisimp1 s
(m- a
))
397 ($diff
`((%cot
) ,(m* '$%pi
'$z
)) '$z s
)))
399 (declare (special $z
))
401 (m*t
(factorial s
) (m^t
(m-t a
) (1- (- s
)))))))))
402 (subfunmakes '$psi
(ncons s
) (ncons a
)))))
405 ;; subtitle polygamma tayloring routines
407 ;; These routines are specially coded to be as fast as possible given the
408 ;; current $TAYLOR; too bad they have to be so ugly.
410 (declare-top (special var subl
*last
* sign last-exp
))
412 (defun expgam-fun (pw temp
)
413 (setq temp
(get-datum (get-key-var (car var
))))
416 (let-pw temp
(e1+ pw
)
417 (psexpt-fn (getexp-fun '(($psi
) -
1) var
(e1+ pw
))))
418 (make-ps var
(ncons pw
) '(((-1 .
1) 1 .
1))))))
420 (defun expplygam-funs (pw subl l
) ; l is a irrelevant here
421 (setq subl
(car subl
))
422 (if (or (not (integerp subl
)) (< subl -
1))
423 (tay-err "Unable to expand at a subscript in")
424 (prog ((e 0) (sign 0) npw
)
425 (declare (fixnum e
) (fixnum sign
))
426 (setq npw
(/ (float (car pw
)) (float (cdr pw
))))
429 `(((1 .
1) .
,(prep1 '((mtimes) -
1 $%gamma
)))))
431 (cons '((-1 .
1) -
1 .
1)
434 .
,(prep1 '((mtimes) -
1 $%gamma
)))))))
435 (t (setq *last
* (factorial subl
))
436 `(((,(- (1+ subl
)) .
1)
437 ,(* (expt -
1 (1+ subl
))
438 (factorial subl
)) .
1))))
439 e
(if (< subl
1) (- subl
) -
1)
440 sign
(if (< subl
1) -
1 (expt -
1 subl
)))
441 a
(setq e
(1+ e
) sign
(- sign
))
442 (if (> e npw
) (return l
)
445 .
,(rctimes (rcplygam e
)
446 (prep1 ($zeta
(+ (1+ subl
) e
))))))))
450 (declare (fixnum k
) )
451 (cond ((= subl -
1) (cons sign k
))
452 ((= subl
0) (cons sign
1))
454 (cons (* sign
*last
*) 1)
456 (quot (* *last
* (+ subl
(1+ k
)))
459 (defun plygam-ord (subl)
460 (if (equal (car subl
) -
1) (ncons (rcone))
461 `((,(m- (m1+ (car subl
))) .
1))))
463 (defun plygam-pole (a c func
)
465 (let ((ps (get-lexp (m- a
(rcdisrep c
)) () t
)))
466 (rplacd (cddr ps
) (cons `((0 .
1) .
,c
) (cdddr ps
)))
467 (if (atom func
) (gam-const a ps func
)
468 (plygam-const a ps func
)))
470 (if (atom func
) `((%gamma
) ,(rcdisrep c
))
471 `((mqapply) ,func
,(rcdisrep c
)))
474 (defun gam-const (a arg func
)
475 (let ((const (ps-lc* arg
)) (arg-c))
476 (cond ((not (rcintegerp const
))
477 (taylor2 (diff-expand `((%gamma
) ,a
) tlist
)))
479 (setq const
(car const
))
480 (if (pscoefp arg
) (setq arg-c
(get-lexp (m+t a
(- const
)) (rcone) (signp le const
))))
481 (if (and arg-c
(not (psp arg-c
)))
482 (taylor2 (simplify `((%gamma
) ,const
)))
483 (let ((datum (get-datum (get-key-var (gvar (or arg-c arg
)))))
484 (ord (if arg-c
(le (terms arg-c
)) (le (n-term (terms arg
))))))
485 (setq func
(current-trunc datum
))
487 (pstimes (let-pw datum
(e- func ord
) (expand (m+t a
(- const
)) '%gamma
))
488 (let-pw datum
(e+ func ord
)
489 (tsprsum (m+t a
(m-t '%%taylor-index%%
))
490 `(%%taylor-index%%
1 ,const
) '%product
)))
491 (pstimes (expand (m+t a
(- const
)) '%gamma
)
492 (let-pw datum
(e+ func ord
)
493 (psexpt (tsprsum (m+t a
'%%taylor-index%%
)
494 `(%%taylor-index%%
0 ,(- (1+ const
))) '%product
)
497 (defun plygam-const (a arg func
)
498 (let ((const (ps-lc* arg
)) (sub (cadr func
)))
500 ((or (not (integerp sub
)) (< sub -
1))
501 (tay-err "Unable to expand at a subscript in"))
502 ((not (rcintegerp const
))
503 (taylor2 (diff-expand `((mqapply) ,func
,a
) tlist
)))
504 (t (setq const
(car const
))
506 (expand (m+t a
(- const
)) func
)
509 (cons (* (expt -
1 sub
) (factorial sub
)) 1)
510 (tsprsum `((mexpt) ,(m+t a
(m-t '%%taylor-index%%
)) ,(- (1+ sub
)))
511 `(%%taylor-index%%
1 ,const
) '%sum
))
513 (cons (* (expt -
1 (1+ sub
)) (factorial sub
)) 1)
514 (tsprsum `((mexpt) ,(m+t a
'%%taylor-index%%
) ,(- (1+ sub
)))
515 `(%%taylor-index%%
0 ,(- (1+ const
))) '%sum
))))))))
517 (declare-top (unspecial var subl
*last
* sign last-exp
))
519 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
520 ;;; Lambert W function
521 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
525 ;; Corless, R. M., Gonnet, D. E. G., Jeffrey, D. J., Knuth, D. E. (1996).
526 ;; "On the Lambert W function". Advances in Computational Mathematics 5:
529 ;; http://www.apmaths.uwo.ca/~djeffrey/Offprints/W-adv-cm.pdf.
530 ;; or http://www.apmaths.uwo.ca/~rcorless/frames/PAPERS/LambertW/
532 ;; D. J. Jeffrey, D. E. G. Hare, R. M. Corless
533 ;; Unwinding the branches of the Lambert W function
534 ;; The Mathematical Scientist, 21, pp 1-7, (1996)
535 ;; http://www.apmaths.uwo.ca/~djeffrey/Offprints/wbranch.pdf
537 ;; Winitzki, S. Uniform Approximations for Transcendental Functions.
538 ;; In Part 1 of Computational Science and its Applications - ICCSA 2003,
539 ;; Lecture Notes in Computer Science, Vol. 2667, Springer-Verlag,
540 ;; Berlin, 2003, 780-789. DOI 10.1007/3-540-44839-X_82
541 ;; http://homepages.physik.uni-muenchen.de/~Winitzki/papers/
544 ;; Having Fun with Lambert W(x) Function
545 ;; arXiv:1003.1628v1, March 2010, http://arxiv.org/abs/1003.1628
547 ;; See also http://en.wikipedia.org/wiki/Lambert's_W_function
549 (defmfun $lambert_w
(z)
550 (simplify (list '(%lambert_w
) (resimplify z
))))
552 ;;; Set properties to give full support to the parser and display
553 (defprop $lambert_w %lambert_w alias
)
554 (defprop $lambert_w %lambert_w verb
)
555 (defprop %lambert_w $lambert_w reversealias
)
556 (defprop %lambert_w $lambert_w noun
)
558 ;;; lambert_w is a simplifying function
559 (defprop %lambert_w simp-lambertw operators
)
561 ;;; Derivative of lambert_w
565 ((mexpt) $%e
((mtimes ) -
1 ((%lambert_w
) x
)))
566 ((mexpt) ((mplus) 1 ((%lambert_w
) x
)) -
1)))
569 ;;; Integral of lambert_w
570 ;;; integrate(W(x),x) := x*(W(x)^2-W(x)+1)/W(x)
576 ((mexpt) ((%lambert_w
) x
) 2)
577 ((mtimes) -
1 ((%lambert_w
) x
))
579 ((mexpt) ((%lambert_w
) x
) -
1)))
582 (defun simp-lambertw (x yy z
)
583 (declare (ignore yy
))
585 (setq x
(simpcheck (cadr x
) z
))
586 (cond ((equal x
0) 0)
588 ((zerop1 x
) ($bfloat
0)) ;bfloat case
592 ((alike1 x
'((mtimes simp
) ((rat simp
) -
1 2) ((%log simp
) 2)))
593 ;; W(-log(2)/2) = -log(2)
594 '((mtimes simp
) -
1 ((%log simp
) 2)))
595 ((alike1 x
'((mtimes simp
) -
1 ((mexpt simp
) $%e -
1)))
598 ((alike1 x
'((mtimes) ((rat) -
1 2) $%pi
))
599 ;; W(-%pi/2) = %i*%pi/2
600 '((mtimes simp
) ((rat simp
) 1 2) $%i $%pi
))
601 ;; numerical evaluation
602 ((complex-float-numerical-eval-p x
)
603 ;; x may be an integer. eg "lambert_w(1),numer;"
605 (to (bigfloat::lambert-w-k
0 (bigfloat:to
($float x
))))
606 (to (bigfloat::lambert-w-k
0 (bigfloat:to x
)))))
607 ((complex-bigfloat-numerical-eval-p x
)
608 (to (bigfloat::lambert-w-k
0 (bigfloat:to x
))))
609 (t (list '(%lambert_w simp
) x
))))
611 ;; An approximation of the k-branch of generalized Lambert W function
613 ;; z real or complex lisp float
614 ;; Used as initial guess for Halley's iteration.
615 ;; When W(z) is real, ensure that guess is real.
616 (defun init-lambert-w-k (k z
)
617 (let ( ; parameters for k = +/- 1 near branch pont z=-1/%e
619 (branch-point (/ -
1 %e-val
))) ; branch pont z=-1/%e
621 ; For principal branch k=0, use expression by Winitzki
622 ((= k
0) (init-lambert-w-0 z
))
623 ; For k=1 branch, near branch point z=-1/%e with im(z) < 0
626 (< (abs (- branch-point z
)) branch-eps
))
627 (bigfloat::lambert-branch-approx z
))
628 ; For k=-1 branch, z real with -1/%e < z < 0
629 ; W(z) is real in this range
630 ((and (= k -
1) (realp z
) (> z branch-point
) (< z
0))
631 (init-lambert-w-minus1 z
))
632 ; For k=-1 branch, near branch point z=-1/%e with im(z) >= 0
635 (< (abs (- branch-point z
)) branch-eps
))
636 (bigfloat::lambert-branch-approx z
))
637 ; Default to asymptotic expansion Corless et al (4.20)
638 ; W_k(z) = log(z) + 2.pi.i.k - log(log(z)+2.pi.i.k)
639 (t (let ((two-pi-i-k (complex 0.0e0
(* 2 pi k
))))
642 (* -
1 (log (+ (log z
) two-pi-i-k
)))))))))
644 ;; Complex value of the principal branch of Lambert's W function in
645 ;; the entire complex plane with relative error less than 1%, given
646 ;; standard branch cuts for sqrt(z) and log(z).
648 (defun init-lambert-w-0 (z)
649 (let ((a 2.344e0
) (b 0.8842e0
) (c 0.9294e0
) (d 0.5106e0
) (e -
1.213e0
)
650 (y (sqrt (+ (* 2 %e-val z
) 2)) ) ) ; y=sqrt(2*%e*z+2)
651 ; w = (2*log(1+B*y)-log(1+C*log(1+D*y))+E)/(1+1/(2*log(1+B*y)+2*A)
653 (+ (* 2 (log (+ 1 (* b y
))))
654 (* -
1 (log (+ 1 (* c
(log (+ 1 (* d y
)))))))
657 (/ 1 (+ (* 2 (log (+ 1 (* b y
)))) (* 2 a
)))))))
659 ;; Approximate k=-1 branch of Lambert's W function over -1/e < z < 0.
660 ;; W(z) is real, so we ensure the starting guess for Halley iteration
663 (defun init-lambert-w-minus1 (z)
666 (merror "z not real in init-lambert-w-minus1"))
667 ((or (< z
(/ -
1 %e-val
)) (plusp z
))
668 (merror "z outside range of approximation in init-lambert-w-minus1"))
669 ;; In the region where W(z) is real
670 ;; -1/e < z < C, use power series about branch point -1/e ~ -0.36787
671 ;; C = -0.3 seems a reasonable crossover point
673 (bigfloat::lambert-branch-approx z
))
674 ;; otherwise C <= z < 0, use iteration W(z) ~ ln(-z)-ln(-W(z))
675 ;; nine iterations are sufficient over -0.3 <= z < 0
676 (t (let* ((ln-z (log (- z
))) (maxiter 9) (w ln-z
))
677 (dotimes (k maxiter w
)
678 (setq w
(- ln-z
(log (- w
)))))))))
680 (in-package #-gcl
#:bigfloat
#+gcl
"BIGFLOAT")
682 ;; Approximate Lambert W(k,z) for k=1 and k=-1 near branch point z=-1/%e
683 ;; using power series in y=-sqrt(2*%e*z+2)
684 ;; for im(z) < 0, approximates k=1 branch
685 ;; for im(z) >= 0, approximates k=-1 branch
687 ;; Corless et al (1996) (4.22)
690 ;; z is a real or complex bigfloat:
691 (defun lambert-branch-approx (z)
692 (let ((y (- (sqrt (+ (* 2 (%e z
) z
) 2)))) ; y=-sqrt(2*%e*z+2)
693 (b0 -
1) (b1 1) (b2 -
1/3) (b3 11/72))
694 (+ b0
(* y
(+ b1
(* y
(+ b2
(* b3 y
))))))))
696 ;; Algorithm based in part on Corless et al (1996).
698 ;; It is Halley's iteration applied to w*exp(w).
701 ;; w[j] exp(w[j]) - z
702 ;; w[j+1] = w[j] - -------------------------------------------------
703 ;; (w[j]+2)(w[j] exp(w[j]) -z)
704 ;; exp(w[j])(w[j]+1) - ---------------------------
707 ;; The algorithm has cubic convergence. Once convergence begins, the
708 ;; number of digits correct at step k is roughly 3 times the number
709 ;; which were correct at step k-1.
711 ;; Convergence can stall using convergence test abs(w[j+1]-w[j]) < prec,
712 ;; as happens for generalized_lambert_w(-1,z) near branch point z = -1/%e
713 ;; Therefore also stop iterating if abs(w[j]*exp(w[j]) - z) << abs(z)
714 (defun lambert-w-k (k z
&key
(maxiter 50))
715 (let ((w (init-lambert-w-k k z
)) we w1e delta
(prec (* 4 (epsilon z
))))
716 (dotimes (i maxiter
(maxima::merror
"lambert-w-k did not converge"))
717 (setq we
(* w
(exp w
)))
718 (when (<= (abs (- z we
)) (* 4 (epsilon z
) (abs z
))) (return))
719 (setq w1e
(* (1+ w
) (exp w
)))
720 (setq delta
(/ (- we z
)
721 (- w1e
(/ (* (+ w
2) (- we z
)) (+ 2 (* 2 w
))))))
723 (when (<= (abs (/ delta w
)) prec
) (return)))
724 ;; Check iteration converged to correct branch
725 (check-lambert-w-k k w z
)
728 (defmethod init-lambert-w-k ((k integer
) (z number
))
729 (maxima::init-lambert-w-k k z
))
731 (defmethod init-lambert-w-k ((k integer
) (z bigfloat
))
732 (bfloat-init-lambert-w-k k z
))
734 (defmethod init-lambert-w-k ((k integer
) (z complex-bigfloat
))
735 (bfloat-init-lambert-w-k k z
))
737 (defun bfloat-init-lambert-w-k (k z
)
738 "Approximate generalized_lambert_w(k,z) for bigfloat: z as initial guess"
739 (let ((branch-point -
0.36787944117144)) ; branch point -1/%e
741 ;; if k=-1, z very close to -1/%e and imag(z)>=0, use power series
743 (or (zerop (imagpart z
))
744 (plusp (imagpart z
)))
745 (< (abs (- z branch-point
)) 1e-10))
746 (lambert-branch-approx z
))
747 ;; if k=1, z very close to -1/%e and imag(z)<0, use power series
749 (minusp (imagpart z
))
750 (< (abs (- z branch-point
)) 1e-10))
751 (lambert-branch-approx z
))
752 ;; Initialize using float value if z is representable as a float
755 (bigfloat (lambert-w-k k
(cl:complex
(float (realpart z
) 1.0)
756 (float (imagpart z
) 1.0))))
757 (bigfloat (lambert-w-k k
(float z
1.0)))))
758 ;; For large z, use Corless et al (4.20)
759 ;; W_k(z) ~ log(z) + 2.pi.i.k - log(log(z)+2.pi.i.k)
761 (let ((log-z (log z
)))
763 (- log-z
(log log-z
))
764 (let* ((i (make-instance 'complex-bigfloat
:imag
(intofp 1)))
765 (two-pi-i-k (* 2 (%pi z
) i k
)))
766 (- (+ log-z two-pi-i-k
)
767 (log (+ log-z two-pi-i-k
))))))))))
769 ;; Check Lambert W iteration converged to the correct branch
770 ;; W_k(z) + ln W_k(z) = ln z, for k = -1 and z in [-1/e,0)
771 ;; = ln z + 2 pi i k, otherwise
772 ;; See Jeffrey, Hare, Corless (1996), eq (12)
774 ;; z, w bigfloat: numbers
775 (defun check-lambert-w-k (k w z
)
776 (let ((tolerance #-gcl
1.0e-6
777 #+gcl
(cl:float
1/1000000)))
780 ;; k=-1 branch with z and w real.
781 ((and (= k -
1) (realp z
) (minusp z
) (>= z
(/ -
1 (%e z
))))
784 (< (abs (+ w
(log w
) (- (log z
)))) tolerance
))
788 ; i k = (W_k(z) + ln W_k(z) - ln(z)) / 2 pi
790 (setq ik
(/ (+ w
(log w
) (- (log z
))) (* 2 (%pi z
))))
791 (if (and (< (realpart ik
) tolerance
)
792 (< (abs (- k
(imagpart ik
))) tolerance
))
796 (maxima::merror
"Lambert W iteration converged to wrong branch"))))
800 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
801 ;;; Generalized Lambert W function
802 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
804 (defmfun $generalized_lambert_w
(k z
)
805 (simplify (list '(%generalized_lambert_w
) (resimplify k
) (resimplify z
))))
807 ;;; Set properties to give full support to the parser and display
808 (defprop $generalized_lambert_w %generalized_lambert_w alias
)
809 (defprop $generalized_lambert_w %generalized_lambert_w verb
)
810 (defprop %generalized_lambert_w $generalized_lambert_w reversealias
)
811 (defprop %generalized_lambert_w $generalized_lambert_w noun
)
813 ;;; lambert_w is a simplifying function
814 (defprop %generalized_lambert_w simp-generalized-lambertw operators
)
816 ;;; Derivative of lambert_w
817 (defprop %generalized_lambert_w
821 ((mexpt) $%e
((mtimes ) -
1 ((%generalized_lambert_w
) k x
)))
822 ((mexpt) ((mplus) 1 ((%generalized_lambert_w
) k x
)) -
1)))
825 ;;; Integral of lambert_w
826 ;;; integrate(W(k,x),x) := x*(W(k,x)^2-W(k,x)+1)/W(k,x)
827 (defprop %generalized_lambert_w
833 ((mexpt) ((%generalized_lambert_w
) k x
) 2)
834 ((mtimes) -
1 ((%generalized_lambert_w
) k x
))
836 ((mexpt) ((%generalized_lambert_w
) k x
) -
1)))
839 (defun simp-generalized-lambertw (expr ignored z
)
840 (declare (ignore ignored
))
842 (let ((k (simpcheck (cadr expr
) z
))
843 (x (simpcheck (caddr expr
) z
)))
845 ;; Numerical evaluation for real or complex x
846 ((and (integerp k
) (complex-float-numerical-eval-p x
))
847 ;; x may be an integer. eg "generalized_lambert_w(0,1),numer;"
849 (to (bigfloat::lambert-w-k k
(bigfloat:to
($float x
))))
850 (to (bigfloat::lambert-w-k k
(bigfloat:to x
)))))
851 ;; Numerical evaluation for real or complex bigfloat x
852 ((and (integerp k
) (complex-bigfloat-numerical-eval-p x
))
853 (to (bigfloat::lambert-w-k k
(bigfloat:to x
))))
854 (t (list '(%generalized_lambert_w simp
) k x
)))))
856 (in-package "BIGFLOAT")
858 (defvar *debug-li-eval
* nil
)
860 (defun li-using-powers-of-log (n x
)
861 ;; When x is on or near the unit circle the other
862 ;; approaches don't work. Use the expansion in powers of
863 ;; log(z) (from cephes cpolylog)
865 ;; li[s](z) = sum(Z(s-j)*(log(z))^j/j!, j = 0, inf)
867 ;; where Z(j) = zeta(j) for j != 1. For j = 1:
869 ;; Z(1) = -log(-log(z)) + sum(1/k, k, 1, s - 1)
872 ;; This is similar to
873 ;; http://functions.wolfram.com/10.08.06.0024.01, but that
874 ;; identity is clearly undefined if v is a positive
875 ;; integer because zeta(1) is undefined.
879 ;; li[3](z) = Z(3) + Z(2)*log(z) + Z(1)*log(z)^2/2!
880 ;; + Z(0)*log(z)^3/3! + sum(Z(-k)*log(z)^(k+4)/(k+4)!,k,1,inf);
882 ;; But Z(-k) = zeta(-k) is 0 if k is even. So
884 ;; li[3](z) = Z(3) + Z(2)*log(z) + Z(1)*log(z)^2/2!
885 ;; + Z(0)*log(z)^3/3! + sum(Z(-(2*k+1))*log(z)^(2*k+4)/(2*k+4)!,k,0,inf);
888 (let ((sum (- (log (- (log x
))))))
890 (loop for k from
1 below n
893 (to (maxima::$zeta
(maxima::to
(float j
(realpart x
)))))))))
894 (let* ((eps (epsilon x
))
900 (term (* (/ top bot
) (to (zfun (- n k
))))
901 (* (/ top bot
) (to (zfun (- n k
))))))
903 (oddp (- n k
)) ;; even terms are all zero
904 (<= (abs term
) (* (abs sum
) eps
))))
905 ;;(format t "~3d: ~A / ~A = ~A~%" k top bot term)
911 ;; If |x| < series-threshold, the series is used.
912 (let ((series-threshold 0.8))
916 (maxima::$zeta
(maxima::to
(float 3 x
))))
918 ;; li[3](-1) = -(1-2^(1-3))*li[3](1)
923 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
925 ;; (See http://functions.wolfram.com/10.08.03.0003.01)
926 (* -
3/4 (to (maxima::$zeta
(maxima::to
(float 3 x
))))))
928 ;; For z not in the interval (0, 1) and for integral n, we
929 ;; have the identity:
931 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
932 ;; + 2 * sum(li[2*r](-1)/(n-2*r)!*log(-z)^(n-2*r), r, 1, floor(n/2))
934 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog/17/02/01/01/0008/)
936 ;; In particular for n = 3:
938 ;; li[3](z) = li[3](1/z) - log(-z)/6*(log(-z)^2+%pi^2)
939 (let* ((lg (log (- x
)))
941 (result (- (li3numer (/ x
))
943 (+ (* lg lg
) (* dpi dpi
))))))
945 ((> (abs x
) series-threshold
)
946 (let ((result (li-using-powers-of-log 3 x
)))
947 ;; For real x, li-using-power-of-log can return a complex
948 ;; number with a tiny imaginary part. Get rid of that
953 ;; Don't use the identity below because the identity seems
954 ;; to be incorrect. For example, for x = -0.862 it returns
955 ;; a complex value with an imaginary part that is not close
956 ;; to zero as expected.
958 ((> (abs x
) series-threshold
)
959 ;; The series converges too slowly so use the identity:
961 ;; li[3](-x/(1-x)) + li[3](1-x) + li[3](x)
962 ;; = li[3](1) + %pi^2/6*log(1-x) - 1/2*log(x)*(log(1-x))^2 + 1/6*(log(1-x))^3
966 ;; 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
967 ;; - li[3](-x/(1-x)) - li[3](1-x)
969 ;; (See http://functions.wolfram.com/10.08.17.0048.01)
974 (decf s
(* 0.5 u u
(log xc
)))
975 (incf s
(/ (* dpi dpi u
) 6))
976 (decf s
(li3numer (- (/ xc x
))))
977 (decf s
(li3numer xc
))
978 (incf s
(li3numer (float 1 x
)))))
980 ;; Sum the power series. threshold determines when the
981 ;; summation has converted.
982 (let* ((threshold (epsilon x
))
985 (incf term
(* 0.125 x x
))
988 (p1 (* p x
) (* p1 x
))
989 (h (/ p1
(* k k k
)) (/ p1
(* k k k
)))
991 ((<= (abs (/ h s
)) threshold
)
995 ;; The series threshold to above sqrt(1/2) because li[2](%i) needs
996 ;; the value of li[2](1/2-%i/2), and the magnitude of the argument
997 ;; is sqrt(1/2) = 0.707. If the threshold is below this, we get
998 ;; into an infinite recursion oscillating between the two args.
999 (let ((series-threshold .75))
1003 ;; %pi^2/6. This follows from the series.
1004 (/ (expt (%pi z
) 2) 6))
1006 ;; -%pi^2/12. From the formula
1008 ;; li[s](-1) = (2^(1-s)-1)*zeta(s)
1010 ;; (See http://functions.wolfram.com/10.08.03.0003.01)
1011 (/ (expt (%pi z
) 2) -
12))
1014 ;; li[2](z) = -li[2](1/z) - 1/2*log(-z)^2 - %pi^2/6,
1016 ;; valid for all z not in the intervale (0, 1).
1018 ;; (See http://functions.wolfram.com/10.08.17.0013.01)
1019 (- (+ (li2numer (/ z
))
1020 (* 0.5 (expt (log (- z
)) 2))
1021 (/ (expt (%pi z
) 2) 6))))
1022 ;; this converges faster when close to unit circle
1023 ((> (abs z
) series-threshold
)
1024 (li-using-powers-of-log 2 z
))
1026 ;;(> (abs z) series-threshold)
1027 ;; For 0.5 <= |z|, where the series would not converge very quickly, use
1029 ;; li[2](z) = li[2](1/(1-z)) + 1/2*log(1-z)^2 - log(-z)*log(1-z) - %pi^2/6
1031 ;; (See http://functions.wolfram.com/10.08.17.0016.01)
1032 ;;(let* ((1-z (- 1 z))
1034 ;; (+ (li2numer (/ 1-z))
1036 ;; (- (* (log (- z))
1038 ;; (- (/ (expt (%pi z) 2) 6)))))
1040 ;; Series evaluation:
1042 ;; li[2](z) = sum(z^k/k^2, k, 1, inf);
1043 (let ((eps (epsilon z
)))
1045 (term z
(* term
(/ (* z k k
)
1047 (sum z
(+ term sum
)))
1048 ((<= (abs (/ term sum
)) eps
)
1051 (defun polylog-power-series (s z
)
1052 ;; Series evaluation:
1054 ;; li[s](z) = sum(z^k/k^s, k, 1, inf);
1055 (let ((eps (epsilon z
)))
1057 (term z
(* term z
(expt (/ (- k
1) k
) s
)))
1058 (sum z
(+ term sum
)))
1059 ((<= (abs (/ term sum
)) eps
)
1060 ;; Return the value and the number of terms used, for
1061 ;; debugging and for helping in determining the series
1065 (defun polylog-log-series (s z
)
1066 ;; When x is on or near the unit circle the other
1067 ;; approaches don't work. Use the expansion in powers of
1068 ;; log(z) (from cephes cpolylog)
1070 ;; li[s](z) = sum(Z(s-j)*(log(z))^j/j!, j = 0, inf)
1072 ;; where Z(j) = zeta(j) for j != 1. For j = 1:
1074 ;; Z(1) = -log(-log(z)) + sum(1/k, k, 1, s - 1)
1078 (let ((sum (- (log (- (log z
))))))
1080 (loop for k from
1 below s
1083 (to (maxima::$zeta
(maxima::to
(float j
(realpart z
)))))))))
1084 (let* ((eps (epsilon z
))
1086 (logx^
2 (* logx logx
))
1090 ;; Compute sum(Z(s-j)*log(z)^j/j!, j = 1, s)
1092 (zf (zfun (- s k
)) (zfun (- s k
)))
1093 (term (* (/ top bot
) zf
)
1094 (* (/ top bot
) zf
)))
1096 (when *debug-li-eval
*
1097 (format t
"~3d: ~A / ~A * ~A => ~A~%" k top bot zf term
))
1099 (setf bot
(* bot
(1+ k
)))
1100 (setf top
(* top logx
)))
1102 (when *debug-li-eval
*
1103 (format t
"s = ~A, sum = ~S top, bot = ~S ~S~%"
1105 ;; Compute the sum for j = s+1 and up. But since
1106 ;; zeta(-k) is 0 for k even, we only every other term.
1107 (do* ((k (+ s
1) (+ k
2))
1108 (zf (zfun (- s k
)) (zfun (- s k
)))
1109 (term (* (/ top bot
) zf
)
1110 (* (/ top bot
) zf
)))
1111 ((<= (abs term
) (* (abs sum
) eps
))
1112 ;; Return the result and the number of terms used for
1113 ;; helping in determining the series threshold and the
1114 ;; log-series threshold.
1116 ;; Note that if z is real and less than 0, li[s](z) is
1117 ;; real. The series can return a tiny complex value in
1118 ;; this case, so we want to clear that out before
1119 ;; returning the answer.
1120 (values (if (and (realp z
) (minusp z
))
1124 (when *debug-li-eval
*
1125 (format t
"~3d: ~A / ~A = ~A~%" k top bot term
))
1127 (setf bot
(* bot
(+ k
1) (+ k
2)))
1128 (setf top
(* top logx^
2))))))
1130 (defun polylog-inversion-formula (s z
)
1131 ;; For z not in the interval (0, 1) and for integral n, we
1132 ;; have the identity:
1134 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1135 ;; + 2 * sum(li[2*r](-1)/(n-2*r)!*log(-z)^(n-2*r), r, 1, floor(n/2))
1137 ;; (See http://functions.wolfram.com/ZetaFunctionsandPolylogarithms/PolyLog/17/02/01/01/0008/)
1141 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1142 ;; + 2 * sum(li[2*m-2*r](-1)/(n-2*m+2*r)!*log(-z)^(n-2*m+2*r), r, 0, m - 1)
1144 ;; where m = floor(n/2). Thus, n-2*m = 0 if n is even and 1 if n is odd.
1146 ;; For n = 2*m, we have
1148 ;; li[2*m](z) = -log(-z)^(2*m)/(2*m)! - li[2*m](1/z)
1149 ;; + 2 * sum(li[2*r](-1)/(2*m-2*r)!*log(-z)^(2*m-2*r), r, 1, m)
1150 ;; = -log(-z)^(2*m)/(2*m)! - li[2*m](1/z)
1151 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,m-1);
1153 ;; For n = 2*m+1, we have
1155 ;; li[2*m+1](z) = -log(-z)^(2*m+1)/(2*m+1)! + li[2*m+1](1/z)
1156 ;; + 2 * sum(li[2*r](-1)/(2*m-2*r + 1)!*log(-z)^(2*m-2*r + 1), r, 1, m)
1157 ;; = -log(-z)^(2*m+1)/(2*m+1)! + li[2*m+1](1/z)
1158 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,m-1);
1161 ;; li[n](z) = -log(-z)^n/n! + (-1)^(n-1)*li[n](1/z)
1162 ;; + 2 * sum((li[2*m-2*r](-1)*log(-z)^(2*r+1))/(2*r+1)!,r,0,floor(n/2)-1);
1163 (let* ((lgz (log (- z
)))
1165 (half-s (floor s
2))
1166 (neg-1 (float -
1 (realpart z
)))
1170 (top (if (oddp s
) lgz
1) (* top lgz^
2))
1171 (bot 1 (* bot
(+ r r -
1) (+ r r
)))
1172 (term (* (li-s-simp (* 2 (- half-s r
)) neg-1
)
1174 (* (li-s-simp (* 2 (- half-s r
)) neg-1
)
1178 (when *debug-li-eval
*
1179 (format t
"r = ~4d: ~A / ~A, ~A; ~A~%" r top bot term sum
)))
1181 (top (if (oddp s
) lgz
1) (* top lgz^
2))
1182 (bot 1 (* bot
(+ r r
) (+ r r
1)))
1183 (term (* (li-s-simp (* 2 (- half-s r
)) neg-1
)
1185 (* (li-s-simp (* 2 (- half-s r
)) neg-1
)
1189 (when *debug-li-eval
*
1190 (format t
"r = ~4d: ~A / ~A, ~A; ~A~%" r top bot term sum
))))
1193 (maxima::take
'(maxima::mfactorial
) s
)))
1194 (* (expt -
1 (- s
1))
1195 (li-s-simp s
(/ z
))))))
1197 (defun li-s-simp (s z
)
1198 (let ((series-threshold 0.5)
1199 (log-series-threshold 2))
1201 (maxima::to
(to 0.0)))
1203 (maxima::$zeta
(maxima::to
(float s z
))))
1205 (- (* (- 1 (expt 2 (- 1 s
)))
1206 (to (li-s-simp s
(- z
))))))
1207 ((<= (abs z
) series-threshold
)
1208 (values (polylog-power-series s z
)))
1209 ((<= (abs z
) log-series-threshold
)
1210 (values (polylog-log-series s z
)))
1212 (polylog-inversion-formula s z
)))))