Rename tab to datatab in numeric/interpol.mac
[maxima.git] / src / hyp.lisp
blob93bd7de6948b8751baa03de2639b0f18567c1404
1 ;;; -*- LISP -*-
2 ;;; ** (c) Copyright 1979 Massachusetts Institute of Technology **
4 (in-package :maxima)
6 ;;; References:
7 ;;;
8 ;;; Definite integration using the generalized hypergeometric functions
9 ;;; Avgoustis, Ioannis Dimitrios
10 ;;; Thesis. 1977. M.S.--Massachusetts Institute of Technology. Dept.
11 ;;; of Electrical Engineering and Computer Science
12 ;;; http://dspace.mit.edu/handle/1721.1/16269
13 ;;;
14 ;;; Avgoustis, I. D., Symbolic Laplace Transforms of Special Functions,
15 ;;; Proceedings of the 1977 MACSYMA Users' Conference, pp 21-41
17 (macsyma-module hyp)
19 (defvar *debug-hyp* nil)
21 (defmvar $prefer_whittaker nil)
23 ;; When T give result in terms of gamma_incomplete and not gamma_incomplete_lower
24 (defmvar $prefer_gamma_incomplete nil)
26 ;; When NIL do not automatically expand polynomials as a result
27 (defmvar $expand_polynomials t)
29 (eval-when
30 (:execute :compile-toplevel)
31 (defmacro simp (x) `(simplifya ,x ()))
33 ;; The macro MABS has been renamed to HYP-MABS in order to
34 ;; avoid conflict with the Maxima symbol MABS. The other
35 ;; M* macros defined here should probably be similarly renamed
36 ;; for consistency. jfa 03/27/2002
38 (defmacro hyp-mabs (x) `(simp `((mabs) ,,x)))
40 (defmacro msqrt (x) `(m^t ,x 1//2))
42 (defmacro mexpt (x) `(m^t '$%e ,x))
44 (defmacro mlog (x) `(simp `((%log) ,,x)))
46 (defmacro msin (x) `(simp `((%sin) ,,x)))
48 (defmacro mcos (x) `(simp `((%cos) ,,x)))
50 (defmacro masin (x) `(simp `((%asin) ,,x)))
52 (defmacro matan (x) `(simp `((%atan) ,,x)))
54 (defmacro mgamma (x) `(simp `((%gamma) ,,x)))
56 (defmacro mbinom (x y) `(simp `((%binomial) ,,x ,,y)))
58 (defmacro merf (x) `(simp `((%erf) ,,x)))
60 (defmacro =1//2 (x) `(alike1 ,x 1//2))
63 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
65 ;;; Functions moved from hypgeo.lisp to this place.
66 ;;; These functions are no longer used in hypgeo.lisp.
68 ;; Gamma function
69 (defun gm (expr)
70 (simplifya (list '(%gamma) expr) nil))
72 ;; sin(x)
73 (defun sin% (arg)
74 (list '(%sin) arg))
76 ;; cos(x)
77 (defun cos% (arg)
78 (list '(%cos) arg))
80 ;; Test if X is a number, either Lisp number or a maxima rational.
81 (defun nump (x)
82 (cond ((atom x)
83 (numberp x))
84 ((not (atom x))
85 (eq (caar (simplifya x nil)) 'rat))))
87 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
89 (defun hyp-integerp (x)
90 ;; In this file, maxima-integerp was used in many places. But it
91 ;; seems that this code expects maxima-integerp to return T when it
92 ;; is really an integer, not something that was declared an integer.
93 ;; But I'm not really sure if this is true everywhere, but it is
94 ;; true in some places.
96 ;; Thus, we replace all calls to maxima-integerp with hyp-integerp,
97 ;; which, for now, returns T only when the arg is an integer.
98 ;; Should we do something more?
99 (and (maxima-integerp x) (integerp x)))
101 ;; Main entry point for simplification of hypergeometric functions.
103 ;; F(a1,a2,a3,...;b1,b2,b3;z)
105 ;; L1 is a (maxima) list of an's, L2 is a (maxima) list of bn's.
106 (defmfun $hgfred (arg-l1 arg-l2 arg)
107 (flet ((arg-ok (a)
108 (and (listp a)
109 (eq (caar a) 'mlist))))
110 (unless (arg-ok arg-l1)
111 (merror (intl:gettext "hgfred: first argument must be a list; found: ~:M") arg-l1))
112 (unless (arg-ok arg-l2)
113 (merror (intl:gettext "hgfred: second argument must be a list; found: ~:M") arg-l2)))
115 ;; Do we really want $radexpand set to '$all? This is probably a
116 ;; bad idea in general, but we'll leave this in for now until we can
117 ;; verify find all of the code that does or does not need this and
118 ;; until we can verify all of the test cases are correct.
119 (let (;;($radexpand '$all)
120 (*par* arg)
121 (*checkcoefsignlist* nil))
122 (hgfsimp-exec (cdr arg-l1) (cdr arg-l2) arg)))
124 (defun hgfsimp-exec (arg-l1 arg-l2 arg)
125 (let* ((l1 (copy-tree arg-l1))
126 (l2 (copy-tree arg-l2))
127 ($exponentialize nil)
128 (res (hgfsimp l1 l2 arg)))
129 ;; I think hgfsimp returns FAIL and UNDEF for cases where it
130 ;; couldn't reduce the function.
131 (cond ((eq res 'fail)
132 (fpqform l1 l2 arg))
133 ((eq res 'undef)
134 '$und)
136 res))))
138 (defun hgfsimp (arg-l1 arg-l2 arg)
139 (prog (resimp listcmdiff)
140 (setq arg-l1 (macsimp arg-l1)
141 arg-l2 (macsimp arg-l2)
142 resimp (simpg arg-l1 arg-l2 arg))
143 (cond ((not (eq (and (consp resimp) (car resimp)) 'fail))
144 (return resimp))
145 ((and (not (zerop1 arg)) ; Do not call splitfpq for a zero argument
146 (setq listcmdiff
147 (intdiffl1l2 (cadr resimp) (caddr resimp))))
148 (return (splitpfq listcmdiff
149 (cadr resimp)
150 (caddr resimp)
151 arg)))
153 (return (dispatch-spec-simp (cadr resimp)
154 (caddr resimp)
155 arg))))))
157 (defun macsimp (expr)
158 (mapcar #'(lambda (index) ($expand index)) expr))
160 ;; Simplify the parameters. If L1 and L2 have common elements, remove
161 ;; them from both L1 and L2.
162 (defun simpg (arg-l1 arg-l2 arg)
163 (let ((il (zl-intersection arg-l1 arg-l2)))
164 (cond ((null il)
165 (simpg-exec arg-l1 arg-l2 arg))
167 (simpg-exec (del il arg-l1)
168 (del il arg-l2)
169 arg)))))
171 (defun del (a b)
172 (cond ((null a) b)
174 (del (cdr a) (delete (car a) b :count 1 :test #'equal)))))
176 ;; Handle the simple cases where the result is either a polynomial, or
177 ;; is undefined because we divide by zero.
178 (defun simpg-exec (arg-l1 arg-l2 arg)
179 (let (n)
180 (cond ((zerop-in-l arg-l1)
181 ;; A zero in the first index means the series terminates
182 ;; after the first term, so the result is always 1.
184 ((setq n (hyp-negp-in-l arg-l1))
185 ;; A negative integer in the first series means we have a
186 ;; polynomial.
187 (create-poly arg-l1 arg-l2 n arg))
188 ((and (or (zerop-in-l arg-l2)
189 (hyp-negp-in-l arg-l2))
190 (every #'mnump arg-l1)
191 (every #'mnump arg-l2))
192 ;; A zero or negative number in the second index means we
193 ;; eventually divide by zero, so we're undefined. But only
194 ;; do this if both indices contain numbers. See Bug
195 ;; 1858964 for discussion.
196 'undef)
198 ;; We failed so more complicated stuff needs to be done.
199 (append (list 'fail) (list arg-l1) (list arg-l2))))))
201 (defun intdiffl1l2 (arg-l1 arg-l2)
202 (cond ((null arg-l1)
203 nil)
205 (intdiff arg-l1 arg-l2))))
207 ;; For each element x in arg-l1 and y in arg-l2, compute d = x - y.
208 ;; Find the smallest such non-negative integer d and return (list x
209 ;; d).
210 (defun intdiff (arg-l1 arg-l2)
211 (let ((result nil))
212 ;; Compute all possible differences between elements in arg-l1 and
213 ;; arg-l2. Only save the ones where the difference is a positive
214 ;; integer
215 (dolist (x arg-l1)
216 (dolist (y arg-l2)
217 (let ((diff (sub x y)))
218 (when (nni diff)
219 (push (list x diff) result)))))
220 ;; Find the smallest one and return it.
221 (let ((min (first result)))
222 (dolist (x (rest result))
223 (when (< (second x) (second min))
224 (setf min x)))
225 min)))
227 ;; Create the appropriate polynomial for the hypergeometric function.
228 (defun create-poly (arg-l1 arg-l2 n arg)
229 (let ((len1 (length arg-l1))
230 (len2 (length arg-l2)))
231 ;; n is the smallest (in magnitude) negative integer in L1. To
232 ;; make everything come out right, we need to make sure this value
233 ;; is first in L1. This is ok, the definition of the
234 ;; hypergeometric function does not depend on the order of values
235 ;; in L1.
236 (setf arg-l1 (cons n (remove n arg-l1 :count 1)))
237 (cond ((and (equal len1 2)
238 (equal len2 1))
239 (2f1polys arg-l1 arg-l2 n arg))
240 ((and (equal len1 1)
241 (equal len2 1))
242 (1f1polys arg-l2 n arg))
243 ((and (equal len1 2)
244 (zerop len2))
245 (2f0polys arg-l1 n arg))
246 (t (create-any-poly arg-l1 arg-l2 (mul -1 n) arg)))))
248 (defun 1f1polys (arg-l2 n arg)
249 (let* ((c (car arg-l2))
250 (n (mul -1 n))
251 (fact1 (mul (power 2 n)
252 (take '(mfactorial) n)
253 (inv (power -1 n))))
254 ;; For all of the polynomials here, I think it's ok to
255 ;; replace sqrt(z^2) with z because when everything is
256 ;; expanded out they evaluate to exactly the same thing. So
257 ;; $radexpand $all is ok here.
258 (fact2 (let (($radexpand '$all))
259 (mul (power 2 '((rat simp) 1 2))
260 (power arg '((rat simp) 1 2))))))
261 (cond ((alike1 c '((rat simp) 1 2))
262 ;; A&S 22.5.56
263 ;; hermite(2*n,x) = (-1)^n*(2*n)!/n!*M(-n,1/2,x^2)
265 ;; So
266 ;; M(-n,1/2,x) = n!/(2*n)!*(-1)^n*hermite(2*n,sqrt(x))
268 ;; But hermite(m,x) = 2^(m/2)*He(sqrt(2)*sqrt(x)), so
270 ;; M(-n,1/2,x) = (-1)^n*n!*2^n/(2*n)!*He(2*n,sqrt(2)*sqrt(x))
271 (mul fact1
272 (inv (take '(mfactorial) (add n n)))
273 (hermpol (add n n) fact2)))
274 ((alike1 c '((rat simp) 3 2))
275 ;; A&S 22.5.57
276 ;; hermite(2*n+1,x) = (-1)^n*(2*n+1)!/n!*M(-n,3/2,x^2)*2*x
278 ;; So
279 ;; M(-n,3/2,x) = n!/(2*n+1)!*(-1)^n*hermite(2*n+1,sqrt(x))/2/sqrt(x)
281 ;; and in terms of He, we get
283 ;; M(-n,3/2,x) = (-1)^n*n!*2^(n-1/2)/(2*n+1)!/sqrt(x)*He(2*n+1,sqrt(2)*sqrt(x))
284 (mul fact1
285 (inv (power 2 '((rat simp) 1 2)))
286 (inv (take '(mfactorial) (add n n 1)))
287 (hermpol (add n n 1) fact2)
288 ;; Similarly, $radexpand here is ok to convert sqrt(z^2) to z.
289 (let (($radexpand '$all))
290 (inv (power arg '((rat simp) 1 2))))))
291 ((alike1 c (neg (add n n)))
292 ;; 1F1(-n; -2*n; z)
293 (mul (power -1 n)
294 (inv (take '(%binomial) (add n n) n))
295 (lagpol n (sub c 1) arg)))
297 ;; A&S 22.5.54:
299 ;; gen_laguerre(n,alpha,x) =
300 ;; binomial(n+alpha,n)*hgfred([-n],[alpha+1],x);
302 ;; Or hgfred([-n],[alpha],x) =
303 ;; gen_laguerre(n,alpha-1,x)/binomial(n+alpha-1,n)
305 ;; But 1/binomial(n+alpha-1,n) = n!*(alpha-1)!/(n+alpha-1)!
306 ;; = n! / (alpha * (alpha + 1) * ... * (alpha + n - 1)
307 ;; = n! / poch(alpha, n)
309 ;; See Bug 1858939.
311 ;; However, if c is not a number leave the result in terms
312 ;; of gamma functions. I (rtoy) think this makes for a
313 ;; simpler result, especially if n is rather large. If the
314 ;; user really wants the answer expanded out, makefact and
315 ;; minfactorial will do that.
316 (if (and (integerp n)
317 (numberp c))
318 (if (or (zerop c)
319 (and (minusp c) (> c (- n))))
320 (merror (intl:gettext "hgfred: 1F1(~M; ~M; ~M) not defined.")
321 (- n) c arg)
322 (mul (take '(mfactorial) n)
323 (inv (take '($pochhammer) c n))
324 (lagpol n (sub c 1) arg)))
325 (let (($gamma_expand t)) ; Expand Gamma function
326 (mul (take '(mfactorial) n)
327 (take '(%gamma) c)
328 (inv (take '(%gamma) (add c n)))
329 (lagpol n (sub c 1) arg))))))))
331 ;; Hermite polynomial. Note: The Hermite polynomial used here is the
332 ;; He polynomial, defined as (A&S 22.5.18, 22.5.19)
334 ;; He(n,x) = 2^(-n/2)*H(n,x/sqrt(2))
336 ;; or
338 ;; H(n,x) = 2^(n/2)*He(x*sqrt(2))
340 ;; We want to use H, as used in specfun, so we need to convert it.
341 (defun hermpol (n arg)
342 (let ((fact (inv (power 2 (div n 2))))
343 (x (mul arg (inv (power 2 '((rat simp) 1 2))))))
344 (mul fact
345 (if (and $expand_polynomials (integerp n))
346 (mfuncall '$hermite n x)
347 (list '($hermite simp) n x)))))
349 ;; Generalized Laguerre polynomial
350 (defun lagpol (n a arg)
351 (if (and (numberp a) (zerop a))
352 (if (and $expand_polynomials (integerp n))
353 (mfuncall '$laguerre n arg)
354 (list '($laguerre simp) n arg))
355 (if (and $expand_polynomials (integerp n))
356 (mfuncall '$gen_laguerre n a arg)
357 (list '($gen_laguerre simp) n a arg))))
359 (defun 2f0polys (arg-l1 n arg)
360 (let ((a (car arg-l1))
361 (b (cadr arg-l1)))
362 (when (alike1 (sub b a) '((rat simp) -1 2))
363 (rotatef a b))
364 (cond ((alike1 (sub b a) '((rat simp) 1 2))
365 ;; 2F0(-n,-n+1/2,z) or 2F0(-n-1/2,-n,z)
366 (interhermpol n a b arg))
368 ;; 2F0(a,b;z)
369 (let ((x (mul -1 (inv arg)))
370 (order (mul -1 n)))
371 (mul (take '(mfactorial) order)
372 (inv (power x order))
373 (inv (power -1 order))
374 (lagpol order (mul -1 (add b order)) x)))))))
376 ;; Compute 2F0(-n,-n+1/2;z) and 2F0(-n-1/2,-n;z) in terms of Hermite
377 ;; polynomials.
379 ;; Ok. I couldn't find any references giving expressions for this, so
380 ;; here's a quick derivation.
382 ;; 2F0(-n,-n+1/2;z) = sum(pochhammer(-n,k)*pochhammer(-n+1/2,k)*z^k/k!, k, 0, n)
384 ;; It's easy to show pochhammer(-n,k) = (-1)^k*n!/(n-k)!
385 ;; Also, it's straightforward but tedious to show that
386 ;; pochhammer(-n+1/2,k) = (-1)^k*(2*n)!*(n-k)!/2^(2*k)/n!/(2*n-2*k)!
388 ;; Thus,
389 ;; 2F0 = (2*n)!*sum(z^k/2^(2*k)/k!/(2*n-2*k)!)
391 ;; Compare this to the expression for He(2*n,x) (A&S 22.3.11):
393 ;; He(2*n,x) = (2*n)! * x^(2*n) * sum((-1)^k*x^(-2*k)/2^k/k!/(2*n-2*k)!)
395 ;; Hence,
397 ;; 2F0(-n,-n+1/2;z) = y^n * He(2*n,y)
399 ;; where y = sqrt(-2/x)
401 ;; For 2F0(-n-1/2,-n;z) = sum(pochhammer(-n,k)*pochhammer(-n-1/2,k)*z^k/k!)
402 ;; we find that
404 ;; pochhammer(-n-1/2,k) = pochhammer(-(n+1)+1/2,k)
405 ;; =
407 ;; So 2F0 = (2*n+1)!*sum(z^k/z^(2*k)/k!/(2*n+1-2*k)!)
409 ;; and finally
411 ;; 2F0(-n-1/2,-n;z) = y^(2*n+1) * He(2*n+1,y)
413 ;; with y as above.
414 (defun interhermpol (n a b x)
415 (let ((arg (power (div 2 (mul -1 x)) '((rat simp) 1 2)))
416 (order (cond ((alike1 a n)
417 (mul -2 n))
418 ((alike1 b n)
419 (sub 1 (add n n))))))
420 ;; 2F0(-n,-n+1/2;z) = y^(-2*n)*He(2*n,y)
421 ;; 2F0(-n-1/2,-n;z) = y^(-(2*n+1))*He(2*n+1,y)
423 ;; where y = sqrt(-2/var);
424 (mul (power arg (mul -1 order))
425 (hermpol order arg))))
427 ;; F(n,b;c;z), where n is a negative integer (number or symbolic).
428 ;; The order of the arguments must be checked by the calling routine.
429 (defun 2f1polys (arg-l1 arg-l2 n arg)
430 (prog (l v lgf)
431 ;; Since F(a,b;c;z) = F(b,a;c;z), make sure L1 has the negative
432 ;; integer first, so we have F(-n,d;c;z)
433 ;; Remark: 2f1polys is only called from create-poly. create-poly calls
434 ;; 2f1polys with the correct order of arg-l1. This test is not necessary.
435 ; (cond ((not (alike1 (car arg-l1) n))
436 ; (setq arg-l1 (reverse arg-l1))))
438 (cond ((mnump *par*)
439 ;; The argument of the hypergeometric function is a number.
440 ;; Avoid the following check which does not work for this case.
441 (setq v (div (add (cadr arg-l1) n) 2)))
443 ;; Check if (b+n)/2 is free of the argument.
444 ;; At this point of the code there is no check of the return value
445 ;; of vfvp. When nil we have no solution and the result is wrong.
446 (setq l (vfvp (div (add (cadr arg-l1) n) 2) arg))
447 (setq v (cdr (assoc 'v l :test #'equal)))))
449 (cond ((and (or (not (integerp n))
450 (not $expand_polynomials))
451 ;; Assuming we have F(-n,b;c;z), then v is (b+n)/2.
452 ;; See if it can be a Legendre function.
453 ;; We call legpol-core because we know that arg-l1 has the
454 ;; correct order. This avoids questions about the parameter b
455 ;; from legpol-core, because legpol calls legpol-core with
456 ;; both order of arguments.
457 (setq lgf (legpol-core (car arg-l1)
458 (cadr arg-l1)
459 (car arg-l2)
460 arg)))
461 (return lgf))
462 ((and (or (not (integerp n))
463 (not $expand_polynomials))
464 (alike1 (sub (car arg-l2) v) '((rat simp) 1 2)))
465 ;; A&S 15.4.5:
466 ;; F(-n, n+2*a; a+1/2; x) = n!*gegen(n, a, 1-2*x)/pochhammer(2*a,n)
468 ;; So v = a, and (car arg-l2) = a + 1/2.
469 (return (mul
470 (cond ((zerop1 v) 1)
471 (t (mul (take '(mfactorial) (mul -1 n))
472 (inv (take '($pochhammer)
473 (mul 2 v)
474 (mul -1 n))))))
475 (gegenpol (mul -1 n)
477 (sub 1 (mul 2 *par*))))))
479 ;; A&S 15.4.6 says
480 ;; F(-n, n + a + 1 + b; a + 1; x)
481 ;; = n!*jacobi_p(n,a,b,1-2*x)/pochhammer(a+1,n)
482 (return (mul (take '(mfactorial) (mul -1 n))
483 (inv (take '($pochhammer) (car arg-l2) (mul -1 n)))
484 (jacobpol (mul -1 n)
485 (add (car arg-l2) -1)
486 (sub (mul 2 v) (car arg-l2))
487 (sub 1 (mul 2 *par*)))))))))
489 ;; Jacobi polynomial
490 (defun jacobpol (n a b x)
491 (if (and $expand_polynomials (integerp n))
492 (mfuncall '$jacobi_p n a b x)
493 (list '($jacobi_p simp) n a b x)))
495 ;; Gegenbauer (Ultraspherical) polynomial. We use ultraspherical to
496 ;; match specfun.
497 (defun gegenpol (n v x)
498 (cond ((equal v 0) (tchebypol n x))
500 (if (and $expand_polynomials (integerp n))
501 (mfuncall '$ultraspherical n v x)
502 `(($ultraspherical simp) ,n ,v ,x)))))
504 ;; Legendre polynomial
505 (defun legenpol (n x)
506 (if (and $expand_polynomials (integerp n))
507 (mfuncall '$legendre_p n x)
508 `(($legendre_p simp) ,n ,x)))
510 ;; Chebyshev polynomial
511 (defun tchebypol (n x)
512 (if (and $expand_polynomials (integerp n))
513 (mfuncall '$chebyshev_t n x)
514 `(($chebyshev_t simp) ,n ,x)))
516 ;; Expand the hypergeometric function as a polynomial. No real checks
517 ;; are made to ensure the hypergeometric function reduces to a
518 ;; polynomial.
519 (defmfun $hgfpoly (arg-l1 arg-l2 arg)
520 (let ((*par* arg)
521 (n (hyp-negp-in-l (cdr arg-l1))))
522 (create-any-poly (cdr arg-l1) (cdr arg-l2) (- n) arg)))
524 (defun create-any-poly (arg-l1 arg-l2 n arg)
525 (prog (result exp prodnum proden)
526 (setq result 1 prodnum 1 proden 1 exp 1)
527 loop
528 (cond ((zerop n) (return result)))
529 (setq prodnum (mul prodnum (mull arg-l1))
530 proden (mul proden (mull arg-l2)))
531 (setq result
532 (add result
533 (mul prodnum
534 (power arg exp)
535 (inv proden)
536 (inv (factorial exp)))))
537 (setq n (sub n 1)
538 exp (add exp 1)
539 arg-l1 (incr1 arg-l1)
540 arg-l2 (incr1 arg-l2))
541 (go loop)))
543 ;; Compute the product of the elements of the list L.
544 (defun mull (l)
545 (reduce #'mul l :initial-value 1))
547 ;; Add 1 to each element of the list L
548 (defun incr1 (l)
549 (mapcar #'(lambda (x) (add x 1)) l))
551 ;; Figure out the orders of generalized hypergeometric function we
552 ;; have and call the right simplifier.
553 (defun dispatch-spec-simp (arg-l1 arg-l2 arg)
554 (let ((len1 (length arg-l1))
555 (len2 (length arg-l2)))
556 (cond ((and (< len1 2)
557 (< len2 2))
558 ;; pFq where p and q < 2.
559 (simp2>f<2 arg-l1 arg-l2 len1 len2 arg))
560 ((and (equal len1 2)
561 (equal len2 1))
562 ;; 2F1
563 (simp2f1 arg-l1 arg-l2 arg))
564 ((and (equal len1 2)
565 (equal len2 0))
566 ;; 2F0(a,b; ; z)
567 (cond ((and (maxima-integerp (car arg-l1))
568 (member ($sign (car arg-l1)) '($neg $nz)))
569 ;; 2F0(-n,b; ; z), n a positive integer
570 (2f0polys arg-l1 (car arg-l1) arg))
571 ((and (maxima-integerp (cadr arg-l1))
572 (member ($sign (cadr arg-l1)) '($neg $nz)))
573 ;; 2F0(a,-n; ; z), n a positive integer
574 (2f0polys (reverse arg-l1) (cadr arg-l1) arg))
576 (fpqform arg-l1 arg-l2 arg))))
577 ((and (= len1 1)
578 (= len2 2))
579 ;; Some 1F2 forms
580 (simp1f2 arg-l1 arg-l2 arg))
581 ((and (= len1 2)
582 (= len2 2))
583 (simp2f2 arg-l1 arg-l2 arg))
585 ;; We don't have simplifiers for any other hypergeometric
586 ;; function.
587 (fpqform arg-l1 arg-l2 arg)))))
589 ;; Handle the cases where the number of indices is less than 2.
590 (defun simp2>f<2 (arg-l1 arg-l2 len1 len2 arg)
591 (cond ((and (zerop len1) (zerop len2))
592 ;; hgfred([],[],z) = e^z
593 (power '$%e arg))
594 ((and (zerop len1) (equal len2 1))
595 (cond
596 ((zerop1 arg)
597 ;; hgfred([],[b],0) = 1
598 (add arg 1))
600 ;; hgfred([],[b],z)
602 ;; The hypergeometric series is then
604 ;; 1+sum(z^k/k!/[b*(b+1)*...(b+k-1)], k, 1, inf)
606 ;; = 1+sum(z^k/k!*gamma(b)/gamma(b+k), k, 1, inf)
607 ;; = sum(z^k/k!*gamma(b)/gamma(b+k), k, 0, inf)
608 ;; = gamma(b)*sum(z^k/k!/gamma(b+k), k, 0, inf)
610 ;; Note that bessel_i(b,z) has the series
612 ;; (z/2)^(b)*sum((z^2/4)^k/k!/gamma(b+k+1), k, 0, inf)
614 ;; bessel_i(b-1,2*sqrt(z))
615 ;; = (sqrt(z))^(b-1)*sum(z^k/k!/gamma(b+k),k,0,inf)
616 ;; = z^((b-1)/2)*hgfred([],[b],z)/gamma(b)
618 ;; So this hypergeometric series is a Bessel I function:
620 ;; hgfred([],[b],z) = bessel_i(b-1,2*sqrt(z))*z^((1-b)/2)*gamma(b)
621 (bestrig (car arg-l2) arg))))
622 ((zerop len2)
623 ;; hgfred([a],[],z) = 1 + sum(binomial(a+k,k)*z^k) = 1/(1-z)^a
624 (power (sub 1 arg) (mul -1 (car arg-l1))))
626 ;; The general case of 1F1, the confluent hypergeomtric function.
627 (confl arg-l1 arg-l2 arg))))
629 ;; Computes
631 ;; bessel_i(a-1,2*sqrt(x))*gamma(a)*x^((1-a)/2)
633 ;; if x > 0
635 ;; or
637 ;; bessel_j(a-1,2*sqrt(x))*gamma(a)*x^((1-a)/2)
639 ;; if x < 0.
641 ;; If a is half of an odd integer and small enough, the Bessel
642 ;; functions are expanded in terms of trig or hyperbolic functions.
644 (defun bestrig (b x)
645 ;; I think it's ok to have $radexpand $all here so that sqrt(z^2) is
646 ;; converted to z.
647 (let (($radexpand '$all))
648 (if (mminusp x)
649 ;; gamma(b)*(-x)^((1-b)/2)*bessel_j(b-1,2*sqrt(-x))
650 (sratsimp (mul (power (neg x) (div (sub 1 b) 2))
651 (take '(%gamma) b)
652 (take '(%bessel_j)
653 (sub b 1)
654 (mul 2 (power (neg x) '((rat simp) 1 2))))))
655 ;; gamma(b)*(x)^((1-b)/2)*bessel_i(b-1,2*sqrt(x))
656 (sratsimp (mul (power x (div (sub 1 b) 2))
657 (take '(%gamma) b)
658 (take '(%bessel_i)
659 (sub b 1)
660 (mul 2 (power x '((rat simp) 1 2)))))))))
662 ;; Kummer's transformation. A&S 13.1.27
664 ;; M(a,b,z) = e^z*M(b-a,b,-z)
665 (defun kummer (arg-l1 arg-l2 arg)
666 (mul (list '(mexpt) '$%e arg)
667 (confl (list (sub (car arg-l2) (car arg-l1)))
668 arg-l2 (mul -1 arg))))
670 ;; Return non-NIL if any element of the list L is zero.
672 (defun zerop-in-l (l)
673 (some #'(lambda (x)
674 (and (numberp x) (zerop x)))
677 ;; If the list L contains a negative integer, return the most positive
678 ;; of the negative integers. Otherwise return NIL.
679 (defun hyp-negp-in-l (l)
680 (let ((max-neg nil))
681 (dolist (x l)
682 (when (and (integerp x) (minusp x))
683 (if max-neg
684 (setf max-neg (max max-neg x))
685 (setf max-neg x))))
686 max-neg))
688 ;; Compute the intersection of L1 and L2, possibly destructively
689 ;; modifying L2. Preserves duplications in L1.
690 (defun zl-intersection (arg-l1 arg-l2)
691 (cond ((null arg-l1) nil)
692 ((member (car arg-l1) arg-l2 :test #'equal)
693 (cons (car arg-l1)
694 (zl-intersection (cdr arg-l1)
695 (delete (car arg-l1) arg-l2 :count 1 :test #'equal))))
696 (t (zl-intersection (cdr arg-l1) arg-l2))))
698 ;; Whittaker M function. A&S 13.1.32 defines it as
700 ;; %m[k,u](z) = exp(-z/2)*z^(u+1/2)*M(1/2+u-k,1+2*u,z)
702 ;; where M is the confluent hypergeometric function.
703 (defun whitfun (k m arg)
704 (list '(mqapply) (list '($%m array) k m) arg))
706 (defun simp1f2 (arg-l1 arg-l2 arg)
707 "Simplify 1F2([a], [b,c], arg). ARG-L1 is the list [a], and ARG-L2 is
708 the list [b, c]. The dependent variable is the (special variable)
709 VAR."
710 (let ((a (car arg-l1)))
711 (destructuring-bind (b c)
712 arg-l2
713 (cond
714 ((and (alike1 a 1//2)
715 (alike1 b 3//2)
716 (alike1 c 3//2))
717 ;; Sine integral, expintegral_si(z) = z*%f[1,2]([1/2],[3/2,3/2], -z^2/4)
718 ;; (See http://functions.wolfram.com/06.37.26.0001.01)
719 ;; Subst z = 2*sqrt(-y) into the above to get
721 ;; expintegral_si(2*sqrt(-y)) = 2*%f[1,2]([1/2],[3/2,3/2], y)*sqrt(-y)
723 ;; Hence %f[1,2]([1/2],[3/2,3/2], y) = expintegral_si(2*sqrt(-y))/2/sqrt(-y)
724 (div (ftake '%expintegral_si (mul 2 (power (neg arg) 1//2)))
725 (mul 2 (power (neg arg) 1//2))))
727 (fpqform arg-l1 arg-l2 arg))))))
729 (defun simp2f2 (arg-l1 arg-l2 arg)
730 (destructuring-bind (a1 a2)
731 arg-l1
732 (destructuring-bind (b1 b2)
733 arg-l2
734 (cond
735 ((and (= a1 1)
736 (= a2 1)
737 (= b1 2)
738 (= b2 2))
739 ;; Cosine integral
740 ;; expintegral_ci(z) = -z^2/4*%f[2,2]([1,1],[2,2], -z^2/4) + log(z) + %gamma
741 ;; (See http://functions.wolfram.com/06.38.26.0001.01)
743 ;; Subst z = 2*sqrt(-y) into the above to get
745 ;; expintegral_ci(2*sqrt(-y)) =
746 ;; %f[2,2]([1,1],[2,2],y)*y + log(2*sqrt(-y)) + %gamma
748 ;; Hence
750 ;; %f[2,2]([1,1],[2,2],y) =
751 ;; (expintegral_ci(2*sqrt(-y)) - log(2*sqrt(-y)) - %gamma)/y
753 ;; But see also http://functions.wolfram.com/06.35.26.0001.01
754 ;; which shows that expintegral_ei can be written in terms of
755 ;; %f[2,2]([1,1],[2,2],z) too. Not sure how to choose
756 ;; between the two representations.
757 (div (sub (ftake '%expintegral_ci (mul 2 (power (neg arg) 1//2)))
758 (add (ftake '%log (mul 2 (power (neg arg) 1//2)))
759 '$%gamma))
760 arg))
762 (fpqform arg-l1 arg-l2 arg))))))
764 (defvar $trace2f1 nil
765 "Enables simple tracing of simp2f1 so you can see how it decides
766 what approach to use to simplify hypergeometric functions")
768 (defun simp2f1 (arg-l1 arg-l2 arg)
769 (prog (a b c lgf)
770 (setq a (car arg-l1) b (cadr arg-l1) c (car arg-l2))
772 (cond ((zerop1 arg)
773 ;; F(a,b; c; 0) = 1
774 (return (add arg 1))))
776 (when $trace2f1
777 (format t "Tracing SIMP2F1~%")
778 (format t " Test a or b negative integer ...~%"))
780 ;; Look if a or b is a symbolic negative integer. The routine
781 ;; 2f1polys handles this case.
782 (cond ((and (maxima-integerp a) (member ($sign a) '($neg $nz)))
783 (return (2f1polys arg-l1 arg-l2 a arg))))
784 (cond ((and (maxima-integerp b) (member ($sign b) '($neg $nz)))
785 (return (2f1polys (list b a) arg-l2 b arg))))
787 (when $trace2f1
788 (format t " Test F(1,1,2)...~%"))
790 (cond ((and (alike1 a 1)
791 (alike1 b 1)
792 (alike1 c 2))
793 ;; F(1,1;2;z) = -log(1-z)/z, A&S 15.1.3
794 (when $trace2f1
795 (format t " Yes~%"))
796 (return (mul (inv (mul -1 arg))
797 (take '(%log) (add 1 (mul -1 arg)))))))
799 (when $trace2f1
800 (format t " Test c = 1/2 or c = 3/2...~%"))
802 (cond ((or (alike1 c '((rat simp) 3 2))
803 (alike1 c '((rat simp) 1 2)))
804 ;; F(a,b; 3/2; z) or F(a,b;1/2;z)
805 (cond ((setq lgf (trig-log (list a b) (list c) arg))
806 (when $trace2f1
807 (format t " Yes: trig-log~%"))
808 (return lgf)))))
810 (when $trace2f1
811 (format t " Test |a-b|=1/2...~%"))
813 (cond ((or (alike1 (sub a b) '((rat simp) 1 2))
814 (alike1 (sub b a) '((rat simp) 1 2)))
815 ;; F(a,b;c;z) where |a-b|=1/2
816 (cond ((setq lgf (hyp-cos a b c arg))
817 (when $trace2f1
818 (format t " Yes: hyp-cos~%"))
819 (return lgf)))))
821 (when $trace2f1
822 (format t " Test a,b are integers, c is a numerical integer...~%"))
824 (cond ((and (maxima-integerp a)
825 (maxima-integerp b)
826 (hyp-integerp c))
827 ;; F(a,b;c;z) when a, and b are integers (or are declared
828 ;; to be integers) and c is a integral number.
829 (setf lgf (simpr2f1 (list a b) (list c) arg))
830 (unless (symbolp lgf) ; Should be more specific! (DK 01/2010)
831 (when $trace2f1
832 (format t " Yes: simpr2f1~%"))
833 (return lgf))))
835 (when $trace2f1
836 (format t " Test a+b and c+1/2 are numerical integers...~%"))
838 (cond ((and (hyp-integerp (add c (inv 2)))
839 (hyp-integerp (add a b)))
840 ;; F(a,b;c;z) where a+b is an integer and c+1/2 is an
841 ;; integer.
842 (when $trace2f1
843 (format t " Yes: step4~%"))
844 (return (step4 a b c arg))))
846 (when $trace2f1
847 (format t " Test a-b+1/2 is a numerical integer...~%"))
849 (cond ((hyp-integerp (add (sub a b) (inv 2)))
850 ;; F(a,b;c,z) where a-b+1/2 is an integer
851 (cond ((setq lgf (step7 a b c arg))
852 (unless (atom lgf)
853 (when $trace2f1
854 (format t " Yes: step7~%"))
855 (return lgf))))))
857 #+nil
858 (when (and (hyp-integerp (add c 1//2))
859 (or (and (hyp-integerp (add a 1//2))
860 (hyp-integerp b))
861 (and (hyp-integerp (add b 1//2))
862 (hyp-integerp a))))
863 (when $trace2f1
864 (format t " Test for atanh: a+1/2, b, and c+1/2 are integers~%"))
865 (return (hyp-atanh a b c arg)))
867 (when (hyp-integerp (add c 1//2))
868 (when $trace2f1
869 (format t " Test for atanh: c+1/2 is an integer~%"))
870 (cond ((and (hyp-integerp (add a 1//2))
871 (hyp-integerp b))
872 (when $trace2f1
873 (format t " atanh with integers a+1/2 and b ~%"))
874 (return (hyp-atanh a b c arg)))
875 ((and (hyp-integerp (add b 1//2))
876 (hyp-integerp a))
877 (when $trace2f1
878 (format t " atanh with integers a and b+1/2 ~%"))
879 (return (hyp-atanh b a c arg)))))
881 (when $trace2f1
882 (format t " Test for Legendre function...~%"))
884 (cond ((setq lgf (legfun a b c arg))
885 (unless (atom lgf)
886 ;; LEGFUN returned something interesting, so we're done.
887 (when $trace2f1
888 (format t " Yes: case 1~%"))
889 (return lgf))
890 ;; LEGFUN didn't return anything, so try it with the args
891 ;; reversed, since F(a,b;c;z) is F(b,a;c;z).
892 (setf lgf (legfun b a c arg))
893 (when lgf
894 (when $trace2f1
895 (format t " Yes: case 2~%"))
896 (return lgf))))
898 (when $trace2f1
899 (format t "'simp2f1-will-continue-in~%"))
900 (return (fpqform arg-l1 arg-l2 arg))))
902 ;; As best as I (rtoy) can tell, step7 is meant to handle an extension
903 ;; of hyp-cos, which handles |a-b|=1/2 and either a+b-1/2 = c or
904 ;; c=a+b+1/2.
906 ;; Based on the code, step7 wants a = s + m/n and c = 2*s+k/l. For
907 ;; hyp-cos, we want c-2*a to be a integer. Or k/l-2*m/n is an
908 ;; integer.
909 #+(or)
910 (progn
911 (defun step7 (a b c)
912 (prog (l m n k mn kl sym sym1 r)
913 ;; Write a = sym + mn, c = sym1 + kl
914 (setq l (s+c a)
915 sym (cdras 'f l)
916 mn (cdras 'c l)
917 l (s+c c)
918 syrm1 (cdras 'f l))
919 ;; We only handle the case where 2*sym = sym1.
920 (cond ((not (equal (mul sym 2) sym1))
921 (return nil)))
922 (setq kl (cdras 'c l))
923 ;; a-b+1/2 is an integer.
924 (setq l (s+c b)
925 r (sub (add (inv 2) (cdras 'c l)) mn)
926 m ($num mn)
927 n ($denom mn)
928 k ($num kl)
929 l ($denom kl))
930 ;; We have a = m*s+m/n, c = 2*m*s+k/l.
931 (cond ((equal (* 2 l) n)
932 (cond ((hyp-integerp (/ (- k m) n))
933 (return (hyp-algv k l m n a b c))))))
934 (cond ((hyp-integerp (/ k (* 2 l)))
935 (cond ((hyp-integerp (/ m n))
936 (return (hyp-algv k l m n a b c)))
937 (t (return nil))))
938 ((hyp-integerp (/ m n))
939 (return nil))
940 ((hyp-integerp (/ (- (* k n) (* 2 l m)) (* 2 l n)))
941 (return (hyp-algv k l m n a b c))))
942 (return nil)))
944 (defun getxy (k l m n)
945 (prog (x y)
946 (setq y 0)
947 loop
948 (cond ((hyp-integerp (setq x (truncate (+ y
949 (truncate k l)
950 (* -2 (quot m n)))
951 2)))
952 (return (list x y))))
953 (incf y 2)
954 (go loop)))
956 (defun hyp-algv (k l m n a b c)
957 (prog (x y xy a-b w)
958 (setq a-b (sub a b))
959 (setq xy (getxy k l m n)
960 x (car xy)
961 y (cadr xy))
962 (cond ((< x 0)(go out)))
963 (cond ((< x y)
964 (cond ((< (add a-b x (inv 2)) 0)
965 (return (f88 x y a b c fun)))
966 (t (return (f87 x y a c fun)))))
968 (cond ((< (add a-b x (inv 2)) 0)
969 (return (f90 x y a c fun)))
970 (t (return (f89 x y a c fun))))))
972 (setq w (* x -1))
973 (cond ((< (- (add a-b (inv 2)) w) 0)
974 (return (f92 x y a c fun)))
975 (t (return (f91 x y a c fun))))))
977 (defun f87 (x y a c fun )
978 (mul
979 (inv (mul (factf c y)
980 (factf (sub (add c y) (add a x)) (- x y))
981 (factf (sub (add c y) (add a x (inv 2)))
982 (sub (add a x (inv 2)) (add a (inv 2))))))
983 (power 'ell (sub 1 c))
984 (power (sub 1 'ell)(sub (add y c) (add a (inv 2))))
985 ($diff (mul (power 'ell (add a x))
986 (power (sub 1 'ell)(mul -1 a))
987 ($diff (mul (power 'ell (sub (add (inv 2) x) y))
988 ($diff (mul (power 'ell (sub (add c y) 1))
989 (power (sub 1 'ell)
990 (sub (add (inv 2)
991 (mul 2 a)
992 (* 2 x))
993 (add c y)))
994 fun)
995 'ell x))
996 'ell (- x y)))
997 'ell y)))
999 (defun f88 (x y a b c fun )
1000 (mul
1001 (inv (mul (factf c y)
1002 (factf (sub (add c y) (add a x)) (- x y))
1003 (factf (add a (inv 2) x)
1004 (sub b (sub x (sub a (inv 2)))))))
1005 (power 'ell (sub 1 c))
1006 (power (sub 1 'ell)(sub (add y c) (add a (inv 2))))
1007 ($diff (mul (power 'ell (add a x))
1008 (power (sub 1 'ell)(mul -1 a))
1009 ($diff (mul (power 'ell (sub c (sub x (sub (inv 2) (mul a 2))))))
1010 (power (sub 1 'ell) (sub (add a x b)(sub c y)))
1011 ($diff (mul (power 'ell (sub b 1 ))
1013 fun)
1014 'ell (sub b (sub a (sub x (inv 2)))))
1015 'ell (- x y)))
1016 'ell y)))
1019 ;; A new version of step7.
1020 (defun step7 (a b c arg)
1021 ;; To get here, we know that a-b+1/2 is an integer. To make further
1022 ;; progress, we want a+b-1/2-c to be an integer too.
1024 ;; Let a-b+1/2 = p and a+b+1/2-c = q where p and q are (numerical)
1025 ;; integers.
1027 ;; A&S 15.2.3 and 15.2.5 allows us to increase or decrease a. Then
1028 ;; F(a,b;c;z) can be written in terms of F(a',b;c;z) where a' = a-p
1029 ;; and a'-b = a-p-b = 1/2. Also, a'+b+1/2-c = a-p+b+1/2-c = q-p =
1030 ;; r, which is also an integer.
1032 ;; A&S 15.2.4 and 15.2.6 allows us to increase or decrease c. Thus,
1033 ;; we can write F(a',b;c;z) in terms of F(a',b;c',z) where c' =
1034 ;; c+r. Now a'-b=1/2 and a'+b+1/2-c' = a-p+b+1/2-c-r =
1035 ;; a+b+1/2-c-p-r = q-p-(q-p)=0.
1037 ;; Thus F(a',b;c';z) is exactly the form we want for hyp-cos. In
1038 ;; fact, it's A&S 15.1.14: F(a,a+1/2,;1+2a;z) =
1039 ;; 2^(2*a)*(1+sqrt(1-z))^(-2*a).
1040 (let ((q (sub (add a b (inv 2))
1041 c)))
1042 (unless (hyp-integerp q)
1043 ;; Wrong form, so return NIL
1044 (return-from step7 nil))
1045 ;; Since F(a,b;c;z) = F(b,a;c;z), we need to figure out which form
1046 ;; to use. The criterion will be the fewest number of derivatives
1047 ;; needed.
1048 (let* ((p1 (add (sub a b) (inv 2)))
1049 (r1 (sub q p1))
1050 (p2 (add (sub b a) (inv 2)))
1051 (r2 (sub q p2)))
1052 (when $trace2f1
1053 (format t "step 7:~%")
1054 (format t " q, p1, r1 = ~A ~A ~A~%" q p1 r1)
1055 (format t " p2, r2 = ~A ~A~%" p2 r2))
1056 (cond ((<= (+ (abs p1) (abs r1))
1057 (+ (abs p2) (abs r2)))
1058 (step7-core a b c arg))
1060 (step7-core b a c arg))))))
1062 (defun step7-core (a b c arg)
1063 (let* ((p (add (sub a b) (inv 2)))
1064 (q (sub (add a b (inv 2))
1066 (r (sub q p))
1067 (a-prime (sub a p))
1068 (c-prime (add 1 (mul 2 a-prime))))
1069 ;; Ok, p and q are integers. We can compute something. There are
1070 ;; four cases to handle depending on the sign of p and r.
1072 ;; We need to differentiate some hypergeometric forms, so use 'ell
1073 ;; as the variable.
1074 (let ((fun (hyp-cos a-prime (add a-prime 1//2) c-prime 'ell)))
1075 ;; fun is F(a',a'+1/2;2*a'+1;z), or NIL
1076 (when fun
1077 (when $trace2f1
1078 (format t "step7-core~%")
1079 (format t " a,b,c = ~A ~A ~A~%" a b c)
1080 (format t " p,q,r = ~A ~A ~A~%" p q r)
1081 (format t " a', c' = ~A ~A~%" a-prime c-prime)
1082 (format t " F(a',a'+1/2; 1+2*a';z) =~%")
1083 (maxima-display fun))
1084 ;; Compute the result, and substitute the actual argument into
1085 ;; result.
1086 (maxima-substitute arg 'ell
1087 (cond ((>= p 0)
1088 (cond ((>= r 0)
1089 (step-7-pp a-prime b c-prime p r 'ell fun))
1091 (step-7-pm a-prime b c-prime p r 'ell fun))))
1093 (cond ((>= r 0)
1094 (step-7-mp a-prime b c-prime p r 'ell fun))
1096 (step-7-mm a-prime b c-prime p r 'ell fun))))))))))
1098 ;; F(a,b;c;z) in terms of F(a',b;c';z)
1100 ;; F(a'+p,b;c'-r;z) where p >= 0, r >= 0.
1101 (defun step-7-pp (a b c p r z fun)
1102 ;; Apply A&S 15.2.4 and 15.2.3
1103 (let ((res (as-15.2.4 a b c r z fun)))
1104 (as-15.2.3 a b (sub c r) p z res)))
1106 ;; p >= 0, r < 0
1108 ;; Let r' = -r
1109 ;; F(a'+p,b;c'-r;z) = F(a'+p,b;c'+r';z)
1110 (defun step-7-pm (a b c p r z fun)
1111 ;; Apply A&S 15.2.6 and 15.2.3
1112 (let ((res (as-15.2.6 a b c (- r) z fun)))
1113 (as-15.2.3 a b (sub c r) p z res)))
1115 ;; p < 0, r >= 0
1117 ;; Let p' = -p
1118 ;; F(a'+p,b;c'-r;z) = F(a'-p',b;c'-r;z)
1119 (defun step-7-mp (a b c p r z fun)
1120 ;; Apply A&S 15.2.4 and 15.2.5
1121 (let ((res (as-15.2.4 a b c r z fun)))
1122 (as-15.2.5 a b (sub c r) (- p) z res)))
1124 ;; p < 0 r < 0
1126 ;; Let p' = - p, r' = -r
1128 ;; F(a'+p,b;c'-r;z) = F(a'-p',b;c'+r';z)
1129 (defun step-7-mm (a b c p r z fun)
1130 ;; Apply A&S 15.2.6 and A&S 15.2.5
1131 (let ((res (as-15.2.6 a b c (- r) z fun)))
1132 (as-15.2.5 a b (sub c r) (- p) z res)))
1134 ;; F(a,b;c;z) when a and b are integers (or declared to be integers)
1135 ;; and c is an integral number.
1136 (defun simpr2f1 (arg-l1 arg-l2 arg)
1137 (destructuring-bind (a b)
1138 arg-l1
1139 (destructuring-bind (c)
1140 arg-l2
1141 (let ((inl1p (hyp-integerp a))
1142 (inl1bp (hyp-integerp b))
1143 (inl2p (hyp-integerp c)))
1144 (cond (inl2p
1145 ;; c is an integer
1146 (cond ((and inl1p inl1bp)
1147 ;; a, b, c are (numerical) integers
1148 (derivint a b c arg))
1149 (inl1p
1150 ;; a and c are integers
1151 (geredno2 b a c arg))
1152 (inl1bp
1153 ;; b and c are integers.
1154 (geredno2 a b c arg))
1155 (t 'fail1)))
1156 ;; Can't really do anything else if c is not an integer.
1157 (inl1p
1158 (cond (inl1bp
1161 'c)))
1162 ((eq (caaar arg-l1) 'rat)
1163 ;; How do we ever get here?
1164 (cond (inl1bp
1167 'd)))
1169 'failg))))))
1171 ;; This isn't called from anywhere?
1172 #+nil
1173 (defun geredno1
1174 (arg-l1 arg-l2)
1175 (cond ((and (> (car arg-l2)(car arg-l1))
1176 (> (car arg-l2)(cadr arg-l1)))
1177 (geredf (car arg-l1) (cadr arg-l1) (car arg-l2) arg))
1178 (t (gered1 arg-l1 arg-l2 #'hgfsimp arg))))
1180 (defun geredno2 (a b c arg)
1181 (cond ((> c b) (geredf b a c arg))
1182 (t (gered2 a b c arg))))
1184 ;; Consider F(1,1;2;z). A&S 15.1.3 says this is equal to -log(1-z)/z.
1186 ;; Apply A&S 15.2.2:
1188 ;; diff(F(1,1;2;z),z,ell) = poch(1,ell)*poch(1,ell)/poch(2,ell)*F(1+ell,1+ell;2+ell;z)
1190 ;; A&S 15.2.7 says:
1192 ;; diff((1-z)^(m+ell)*F(1+ell;1+ell;2+ell;z),z,m)
1193 ;; = (-1)^m*poch(1+ell,m)*poch(1,m)/poch(2+ell,m)*(1-z)^ell*F(1+ell+m,1+ell;2+ell+m;z)
1195 ;; A&S 15.2.6 gives
1197 ;; diff((1-z)^ell*F(1+ell+m,1+ell;2+ell+m;z),z,n)
1198 ;; = poch(1,n)*poch(1+m,n)/poch(2+ell+m,n)*(1-z)^(ell-n)*F(1+ell+m,1+ell;2+ell+m+n;z)
1200 ;; The derivation above assumes that ell, m, and n are all
1201 ;; non-negative integers. Thus, F(a,b;c;z), where a, b, and c are
1202 ;; integers and a <= b <= c, can be written in terms of F(1,1;2;z).
1203 ;; The result also holds for b <= a <= c, of course.
1205 ;; Also note that the last two differentiations can be combined into
1206 ;; one differention since the result of the first is in exactly the
1207 ;; form required for the second. The code below does one
1208 ;; differentiation.
1210 ;; So if a = 1+ell, b = 1+ell+m, and c = 2+ell+m+n, we have ell = a-1,
1211 ;; m = b - a, and n = c - ell - m - 2 = c - b - 1.
1213 (defun derivint (a b c arg)
1214 (if (> a b)
1215 (derivint b a c arg)
1216 (let ((l (- a 1))
1217 (m (- b a))
1218 (n (- c b 1))
1219 (psey (gensym))
1220 result)
1222 (setq result
1223 (mul (power -1 m)
1224 (factorial (+ n m l 1))
1225 (inv (factorial n))
1226 (inv (factorial l))
1227 (inv (factorial (+ n m)))
1228 (inv (factorial (+ m l)))
1229 (power (sub 1 psey) (sub n l))
1230 ($diff (mul (power (sub 1 psey) (+ m l))
1231 ($diff (mul (power psey -1)
1233 (take '(%log) (sub 1 psey)))
1234 psey
1236 psey
1237 (+ n m))))
1238 (if (onep1 arg)
1239 ($limit result psey arg)
1240 (maxima-substitute arg psey result)))))
1242 ;; Handle F(a, b; c; z) for certain values of a, b, and c. See the
1243 ;; comments below for these special values.
1244 (defun hyp-cos (a b c z)
1245 (let ((a1 (div (sub (add a b) (div 1 2)) 2))
1246 (z1 (sub 1 z)))
1247 ;; a1 = (a+b-1/2)/2
1248 ;; z1 = 1-z
1249 (cond ((eql 0 ($ratsimp (sub (sub (add a b)
1250 (div 1 2))
1251 c)))
1252 ;; a+b-1/2 = c
1254 ;; A&S 15.1.14
1256 ;; F(a,a+1/2;2*a;z)
1257 ;; = 2^(2*a-1)*(1-z)^(-1/2)*(1+sqrt(1-z))^(1-2*a)
1259 ;; But if 1-2*a is a negative integer, let's rationalize the answer to give
1261 ;; F(a,a+1/2;2*a;z)
1262 ;; = 2^(2*a-1)*(1-z)^(-1/2)*(1-sqrt(1-z))^(2*a-1)/z^(2*a-1)
1263 (when $trace2f1
1264 (format t " Case a+b-1/2=c~%"))
1265 (let ((2a-1 (sub (mul a1 2) 1)))
1266 (cond ((and (integerp 2a-1) (plusp 2a-1))
1267 ;; 2^(2*a-1)*(1-z)^(-1/2)*(1-sqrt(1-z))^(2*a-1)/z^(2*a-1)
1268 (mul (power 2 2a-1)
1269 (inv (power z1 1//2))
1270 (power (sub 1 (power z1 1//2)) 2a-1)
1271 (inv (power z 2a-1))))
1273 ;; 2^(2*a-1)*(1-z)^(-1/2)*(1+sqrt(1-z))^(1-2*a)
1274 (mul (power 2 (sub (mul a1 2) 1))
1275 (inv (power z1 (div 1 2)))
1276 (power (add 1
1277 (power z1
1278 (div 1
1279 2)))
1280 (sub 1 (mul 2 a1))))))))
1281 ((eql 0 ($ratsimp (sub (add 1 (mul 2 a1)) c)))
1282 ;; c = 1+2*a1 = a+b+1/2
1284 ;; A&S 15.1.13:
1286 ;; F(a,1/2+a;1+2*a;z) = 2^(2*a)*(1+sqrt(1-z))^(-2*a)
1288 ;; But if 2*a is a positive integer, let's rationalize the answer to give
1290 ;; F(a,1/2+a;1+2*a;z) = 2^(2*a)*(1-sqrt(1-z))^(2*a)/z^(2*a)
1291 (when $trace2f1
1292 (format t " Case c=1+2*a=a+b+1/2~%"))
1293 (let ((2a (sub c 1)))
1294 (cond ((and (integerp 2a) (plusp 2a))
1295 ;; 2^(2*a)*(1-sqrt(1-z))^(2*a)/z^(2*a)
1296 (mul (power 2 2a)
1297 (power (sub 1 (power z1 1//2))
1299 (power z (mul -1 2a))))
1301 ;; 2^(2*a)*(1+sqrt(1-z))^(-2*a)
1302 (mul (power 2 2a)
1303 (power (add 1 (power z1 1//2))
1304 (mul -1 2a))))))))))
1306 ;; Is A a non-negative integer?
1307 (defun nni (a)
1308 (cond ((hyp-integerp a)
1309 (not (minusp a)))))
1312 ;;; The following code doesn't appear to be used at all. Comment it all out for now.
1314 (defun degen2f1
1315 (a b c)
1316 (cond ((eq (quest (sub c b)) '$negative)
1317 (cond ((eq (quest (sub c a)) '$negative)
1318 (gered1 (list a b)(list c) #'hgfsimp))
1319 (t (gered2 a b c))))
1320 ((eq (quest (sub c a)) '$negative)(gered2 b a c))
1321 (t (rest-degen a b c))))
1324 (defun rest-degen
1325 (a b c)
1326 (prog(m n l)
1327 (cond ((nni (setq m (sub a 1)))
1328 (return (rest-degen-1 a b c m))))
1329 (cond ((ni b)(return (rest-degen-2 a b c))))
1330 (cond ((and (nni (setq n (sub c 1)))
1331 (nni (setq m (sub (sub a n) 1)))
1332 (nni (setq l (sub b a)))
1333 (eq (sub (sub c a) b)
1334 (mul -1 (add m m n l 1))))
1335 (return (gered1 (list a b)
1336 (list c)
1337 #'hgfsimp))))
1338 (return (hyp-deg b a c))))
1341 (defun rest-degen-1
1342 (a b c m)
1343 (prog(n l)
1344 (cond ((and (ni b)
1345 (ni (sub (sub c a) b))
1346 (nni (sub (sub c a) 1)))
1347 (return (deg299 a b c))))
1348 (cond ((and (nni (setq n (sub (sub c m) 2)))
1349 (nni (setq l (sub b c)))
1350 (equal (sub (sub c a) b)
1351 (mul -1 (add l m 1))))
1352 (return (gered1 (list a b)
1353 (list c)
1354 #'hgfsimp))))
1355 (cond ((nni (setq l (sub (sub b m) 1)))
1356 (return (rest-degen-1a a b c m l))))
1357 (return (hyp-deg b a c))))
1360 (defun rest-degen-1a
1361 (a b c m l)
1362 (prog(n)
1363 (cond ((and (nni (setq n
1364 (sub (sub (sub c m) l) 2)))
1365 (equal (sub n m)(sub (sub c a) b)))
1366 (return (deg2913 a b c))))
1367 (cond ((and (equal c (mul -1 n))
1368 (equal (sub (sub c a) b)
1369 (mul -1 (add m m l n 2))))
1370 (return (deg2918 a b c))))
1371 (return (hyp-deg b a c))))
1374 (defun rest-degen-2
1375 (a b c)
1376 (prog(m l)
1377 (cond ((and (ni c)(ni (sub (sub c a) b)))
1378 (return (rest-degen-2a a b c))))
1379 (cond ((and (nni (setq m (sub c 1)))
1380 (nni (setq l (sub a c)))
1381 (ni (sub (sub c a) b)))
1382 (return (deg292 a b c))))
1383 (return (hyp-deg b a c))))
1386 (defun rest-degen-2a
1387 (a b c)
1388 (prog()
1389 (cond ((nni (sub a c))
1390 (return (gered1 (list a b)
1391 (list c)
1392 #'hgfsimp))))
1393 (cond ((nni (sub (sub c a) 1))
1394 (return (deg2917 a b c))))
1395 (return (hyp-deg b a c))))
1397 (defun quest
1399 (cond ((numberp a)(checksigntm a))
1400 ((equal (caar a) 'rat)(checksigntm a))
1401 (t nil)))
1403 (defun ni(a)(not (hyp-integerp a)))
1406 (defun hyp-deg
1407 (a b c)
1408 (prog()
1409 (cond (fldeg (setq fldeg nil)
1410 (return (hgfsimp (list a b)
1411 (list c)
1412 var))))
1413 (setq fldeg t)
1414 (return (fpqform (list a b)(list c) var))))
1417 (defun deg2913
1418 (a b c)
1419 (mul (power (mul -1 var)(mul -1 b))
1420 (hgfsimp (list (add b 1 (mul -1 c)) b)
1421 (list (add b 1 (mul -1 a)))
1422 (inv var))))
1425 (defun deg2918
1426 (a b c)
1427 (mul (power var (sub 1 c))
1428 (power (sub 1 var)(add c (mul -1 a)(mul -1 b)))
1429 (hgfsimp (list (sub 1 a)(sub 1 b))
1430 (list (sub 2 c))
1431 var)))
1434 (defun deg2917
1435 (a b c)
1436 (mul (power var (sub 1 c))
1437 (hgfsimp (list (add a 1 (mul -1 c))
1438 (add b 1 (mul -1 c)))
1439 (list (sub 2 c))
1440 var)))
1443 (defun deg299
1444 (a b c)
1445 (mul (power (mul -1 var)(mul -1 a))
1446 (hgfsimp (list a (add a 1 (mul -1 c)))
1447 (list (add a 1 (mul -1 b)))
1448 (inv var))))
1453 ;; Is F(a, b; c; z) is Legendre function?
1455 ;; Lemma 29 (see ref) says F(a, b; c; z) can be reduced to a Legendre
1456 ;; function if two of the numbers 1-c, +/-(a-b), and +/- (c-a-b) are
1457 ;; equal to each other or one of them equals +/- 1/2.
1459 ;; This routine checks for each of the possibilities.
1460 (defun legfun (a b c arg)
1461 (let ((1-c (sub 1 c))
1462 (a-b (sub a b))
1463 (c-a-b (sub (sub c a) b))
1464 (inv2 (inv 2)))
1465 (cond ((alike1 a-b inv2)
1466 ;; a-b = 1/2
1467 (when $trace2f1
1468 (format t "Legendre a-b = 1/2~%"))
1469 (gered1 (list a b) (list c) #'legf24 arg))
1471 ((alike1 a-b (mul -1 inv2))
1472 ;; a-b = -1/2
1474 ;; For example F(a,a+1/2;c;x)
1475 (when $trace2f1
1476 (format t "Legendre a-b = -1/2~%"))
1477 (legf24 (list a b) (list c) arg))
1479 ((alike1 c-a-b '((rat simp) 1 2))
1480 ;; c-a-b = 1/2
1481 ;; For example F(a,b;a+b+1/2;z)
1482 (when $trace2f1
1483 (format t "Legendre c-a-b = 1/2~%"))
1484 (legf20 (list a b) (list c) arg))
1486 ((and (alike1 c-a-b '((rat simp) 3 2))
1487 (not (alike1 c 1))
1488 (not (alike1 a -1//2))
1489 (not (alike1 b -1//2)))
1490 ;; c-a-b = 3/2 e.g. F(a,b;a+b+3/2;z) Reduce to
1491 ;; F(a,b;a+b+1/2) and use A&S 15.2.6. But if c = 1, we
1492 ;; don't want to reduce c to 0! Problem: The derivative of
1493 ;; assoc_legendre_p introduces a unit_step function and the
1494 ;; result looks very complicated. And this doesn't work if
1495 ;; a+1/2 or b+1/2 is zero, so skip that too.
1496 (when $trace2f1
1497 (format t "Legendre c-a-b = 3/2~%")
1498 (mformat t " : a = ~A~%" a)
1499 (mformat t " : b = ~A~%" b)
1500 (mformat t " : c = ~A~%" c))
1501 (let ((psey (gensym)))
1502 (maxima-substitute
1503 *par* psey
1504 (mul (power (sub 1 psey) '((rat simp) 3 2))
1505 (add a b '((rat simp) 1 2))
1506 (inv (add b '((rat simp) 1 2)))
1507 (inv (add a '((rat simp) 1 2)))
1508 ($diff (mul (power (sub 1 psey) '((rat simp) -1 2))
1509 (hgfsimp (list a b)
1510 (list (add a b '((rat simp) 1 2)))
1511 psey))
1512 psey
1513 1)))))
1515 ((alike1 c-a-b '((rat simp) -1 2))
1516 ;; c-a-b = -1/2, e.g. F(a,b; a+b-1/2; z)
1517 ;; This case is reduce to F(a,b; a+b+1/2; z) with
1518 ;; F(a,b;c;z) = (1-z)^(c-a-b)*F(c-a,c-b;c;z)
1519 (when $trace2f1
1520 (format t "Legendre c-a-b = -1/2~%"))
1521 (gered1 (list a b) (list c) #'legf20 arg))
1523 ((alike1 1-c a-b)
1524 ;; 1-c = a-b, F(a,b; b-a+1; z)
1525 (when $trace2f1
1526 (format t "Legendre 1-c = a-b~%"))
1527 (gered1 (list a b) (list c) #'legf16 arg))
1529 ((alike1 1-c (mul -1 a-b))
1530 ;; 1-c = b-a, e.g. F(a,b; a-b+1; z)
1531 (when $trace2f1
1532 (format t "Legendre 1-c = b-a~%"))
1533 (legf16 (list a b) (list c) arg))
1535 ((alike1 1-c c-a-b)
1536 ;; 1-c = c-a-b, e.g. F(a,b; (a+b+1)/2; z)
1537 (when $trace2f1
1538 (format t "Legendre 1-c = c-a-b~%"))
1539 (gered1 (list a b) (list c) #'legf14 arg))
1541 ((alike1 1-c (mul -1 c-a-b))
1542 ;; 1-c = a+b-c
1544 ;; For example F(a,1-a;c;x)
1545 (when $trace2f1
1546 (format t "Legendre 1-c = a+b-c~%"))
1547 (legf14 (list a b) (list c) arg))
1549 ((alike1 a-b (mul -1 c-a-b))
1550 ;; a-b = a+b-c, e.g. F(a,b;2*b;z)
1551 (when $trace2f1
1552 (format t "Legendre a-b = a+b-c~%"))
1553 (legf36 (list a b) (list c) arg))
1555 ((or (alike1 1-c inv2)
1556 (alike1 1-c (mul -1 inv2)))
1557 ;; 1-c = 1/2 or 1-c = -1/2
1558 ;; For example F(a,b;1/2;z) or F(a,b;3/2;z)
1559 (when $trace2f1
1560 (format t "Legendre |1-c| = 1/2~%"))
1561 ;; At this point it does not make sense to call legpol. legpol can
1562 ;; handle only cases with a negative integer in the first argument
1563 ;; list. This has been tested already. Therefore we can not get
1564 ;; a result from legpol. For this case a special algorithm is needed.
1565 ;; At this time we return nil.
1566 ;(legpol a b c)
1567 nil)
1569 ((alike1 a-b c-a-b)
1570 ;; a-b = c-a-b
1571 (when $trace2f1
1572 (format t "Legendre a-b = c-a-b~%"))
1573 'legendre-funct-to-be-discovered)
1575 nil))))
1577 ;;; The following legf<n> functions correspond to formulas in Higher
1578 ;;; Transcendental Functions. See the chapter on Legendre functions,
1579 ;;; in particular the table on page 124ff,
1581 ;; Handle the case c-a-b = 1/2
1583 ;; Formula 20:
1585 ;; P(n,m,z) = 2^m*(z^2-1)^(-m/2)/gamma(1-m)*F(1/2+n/2-m/2, -n/2-m/2; 1-m; 1-z^2)
1587 ;; See also A&S 15.4.12 and 15.4.13.
1589 ;; Let a = 1/2+n/2-m/2, b = -n/2-m/2, c = 1-m. Then, m = 1-c. And we
1590 ;; have two equivalent expressions for n:
1592 ;; n = c - 2*b - 1 or n = 2*a - c
1594 ;; The code below chooses the first solution. A&S chooses second.
1596 ;; F(a,b;c;w) = 2^(c-1)*gamma(c)*(-w)^((1-c)/2)*P(c-2*b-1,1-c,sqrt(1-w))
1599 (defun legf20 (arg-l1 arg-l2 arg)
1600 ;; F(a,b;a+b+1/2;x)
1601 (let* (($radexpand nil)
1602 (b (cadr arg-l1))
1603 (c (car arg-l2))
1604 (a (sub (sub c b) '((rat simp) 1 2)))
1605 (m (sub 1 c))
1606 (n (mul -1 (add b b m))))
1607 ;; m = 1 - c
1608 ;; n = -(2*b+1-c) = c - 1 - 2*b
1609 ;; A&S 15.4.13
1611 ;; F(a,b;a+b+1/2;x) = 2^(a+b-1/2)*gamma(a+b+1/2)*x^((1/2-a-b)/2)
1612 ;; *assoc_legendre_p(a-b-1/2,1/2-a-b,sqrt(1-x))
1613 ;; This formula is correct for all arguments x.
1614 (mul (power 2 (add a b '((rat simp) -1 2)))
1615 (take '(%gamma) (add a b '((rat simp) 1 2)))
1616 (power arg
1617 (div (sub '((rat simp) 1 2) (add a b))
1619 (legen n
1621 (power (sub 1 arg) '((rat simp) 1 2))
1622 '$p))))
1624 ;; Handle the case a-b = -1/2.
1626 ;; Formula 24:
1628 ;; P(n,m,z) = 2^m*(z^2-1)^(-m/2)*z^(n+m)/gamma(1-m)*F(-n/2-m/2,1/2-n/2-m/2;1-m;1-1/z^2)
1630 ;; See also A&S 15.4.10 and 15.4.11.
1632 ;; Let a = -n/2-m/2, b = 1/2-n/2-m/2, c = 1-m. Then m = 1-c. Again,
1633 ;; we have 2 possible (equivalent) values for n:
1635 ;; n = -(2*a + 1 - c) or n = c-2*b
1637 ;; The code below chooses the first solution.
1639 ;; F(a,b;c;w) = 2^(c-1)*w^(1/2-c/2)*(1-w)^(c/2-a-1/2)*P(c-2*a-1,1-c,1/sqrt(1-w))
1641 ;; F(a,b;c;w) = 2^(c-1)*w^(1/2-c/2)*(1-w)^(c/2-b)*P(c-2*b,1-c,sqrt(1-w))
1643 ;; Is there a mistake in 15.4.10 and 15.4.11?
1645 (defun legf24 (arg-l1 arg-l2 arg)
1646 (let* (($radexpand nil)
1647 (a (car arg-l1))
1648 (c (car arg-l2))
1649 (m (sub 1 c))
1650 ; (n (mul -1 (add a a m))) ; This is not 2*a-c
1651 (n (sub (add a a) c)) ; but this.
1652 (z (inv (power (sub 1 arg) (inv 2)))))
1653 ;; A&S 15.4.10, 15.4.11
1654 (cond ((eq (asksign arg) '$negative)
1655 ;; A&S 15.4.11
1657 ;; F(a,a+1/2;c;x) = 2^(c-1)*gamma(c)(-x)^(1/2-c/2)*(1-x)^(c/2-a-1/2)
1658 ;; *assoc_legendre_p(2*a-c,1-c,1/sqrt(1-x))
1659 (mul (inv (power 2 m))
1660 (gm (sub 1 m))
1661 (power (mul -1 arg) (div m 2))
1662 (power (sub 1 arg) (sub (div m -2) a))
1663 (legen n
1666 '$p)))
1668 (mul (inv (power 2 m))
1669 (gm (sub 1 m))
1670 (power arg (div m 2))
1671 (power (sub 1 arg) (sub (div m -2) a))
1672 (legen n
1675 '$p))))))
1677 ;; Handle 1-c = a-b
1679 ;; Formula 16:
1681 ;; P(n,m,z) = 2^(-n)*(z+1)^(m/2+n)*(z-1)^(-m/2)/gamma(1-m)*F(-n,-n-m;1-m;(z-1)/(z+1))
1683 ;; See also A&S 15.4.14 and 15.4.15.
1685 ;; Let a = -n, b = -n-m, c = 1-m. Then m = 1-c. We have 2 solutions
1686 ;; for n:
1688 ;; n = -a or n = c-b-1.
1690 ;; The code below chooses the first solution.
1692 ;; F(a,b;c;w) = gamma(c)*w^(1/2-c/2)*(1-w)^(-a)*P(-a,1-c,(1+w)/(1-w));
1694 ;; FIXME: We don't correctly handle the branch cut here!
1695 (defun legf16 (arg-l1 arg-l2 arg)
1696 (let* (($radexpand nil)
1697 (a (car arg-l1))
1698 (c (car arg-l2))
1699 ;; m = 1-c = b-a
1700 (m (sub 1 c))
1701 ;; n = -b
1702 ;; m = b - a, so b = a + m
1703 (b (add a m))
1704 (n (mul -1 b))
1705 (z (div (add 1 arg)
1706 (sub 1 arg))))
1707 (when $trace2f1
1708 (format t "a, c = ~A ~A~%" a c)
1709 (format t "b = ~A~%" b))
1710 ;; A&S 15.4.14, 15.4.15
1711 (cond ((eq (asksign arg) '$negative)
1712 ;; A&S 15.4.15
1714 ;; F(a,b;a-b+1,x) = gamma(a-b+1)*(1-x)^(-b)*(-x)^(b/2-a/2)
1715 ;; * assoc_legendre_p(-b,b-a,(1+x)/(1-x))
1717 ;; for x < 0
1718 (mul (take '(%gamma) c)
1719 (power (sub 1 arg) (mul -1 b))
1720 (power (mul -1 arg) (div m 2))
1721 (legen n m z '$p)))
1723 (mul (take '(%gamma) c)
1724 (power (sub 1 arg) (mul -1 b))
1725 (power arg (div m 2))
1726 (legen n m z '$p))))))
1728 ;; Handle the case 1-c = a+b-c.
1730 ;; See, for example, A&S 8.1.2 (which
1731 ;; might have a bug?) or
1732 ;; http://functions.wolfram.com/HypergeometricFunctions/LegendreP2General/26/01/02/
1734 ;; Formula 14:
1736 ;; P(n,m,z) = (z+1)^(m/2)*(z-1)^(-m/2)/gamma(1-m)*F(-n,1+n;1-m;(1-z)/2)
1738 ;; See also A&S 8.1.2, 15.4.16, 15.4.17
1740 ;; Let a=-n, b = 1+n, c = 1-m. Then m = 1-c and n has 2 solutions:
1742 ;; n = -a or n = b - 1.
1744 ;; The code belows chooses the first solution.
1746 ;; F(a,b;c;w) = gamma(c)*(-w)^(1/2-c/2)*(1-w)^(c/2-1/2)*P(-a,1-c,1-2*w)
1747 (defun legf14 (arg-l1 arg-l2 arg)
1748 ;; Set $radexpand to NIL, because we don't want (-z)^x getting
1749 ;; expanded to (-1)^x*z^x because that's wrong for this.
1750 (let* (($radexpand nil)
1751 (a (first arg-l1))
1752 (b (second arg-l1))
1753 (c (first arg-l2))
1754 (m (sub 1 c))
1755 (n (mul -1 a))
1756 (z (sub 1 (mul 2 arg))))
1757 (when $trace2f1
1758 (format t "~&legf14~%"))
1759 ;; A&S 15.4.16, 15.4.17
1760 (cond ((not (alike1 (add a b) 1))
1761 ;; I think 15.4.16 and 15.4.17 require the form
1762 ;; F(a,1-a;c;x). That is, a+b = 1. If we don't have it
1763 ;; exit now.
1764 nil)
1765 ((and (eq (asksign arg) '$positive)
1766 (eq (asksign (sub 1 arg)) '$positive))
1767 (when $trace2f1
1768 (format t " A&S 15.4.17~%"))
1769 ;; A&S 15.4.17
1771 ;; F(a,1-a;c;x) = gamma(c)*x^(1/2-c/2)*(1-x)^(c/2-1/2)*
1772 ;; assoc_legendre_p(-a,1-c,1-2*x)
1774 ;; for 0 < x < 1
1775 (mul (gm c)
1776 (power arg (div (sub 1 c) 2))
1777 (power (sub 1 arg) (div (sub c 1) 2))
1778 (legen n m z '$p)))
1780 ;; A&S 15.4.16
1782 ;; F(a,1-a;c;z) = gamma(c)*(-z)^(1/2-c/2)*(1-z)^(c/2-1/2)*
1783 ;; assoc_legendre_p(-a,1-c,1-2*z)
1784 (when $trace2f1
1785 (format t " A&S 15.4.17~%"))
1786 (mul (gm c)
1787 (power (mul -1 arg) (div (sub 1 c) 2))
1788 (power (sub 1 arg) (div (sub c 1) 2))
1789 (legen n m z '$p))))))
1791 ;; Handle a-b = a+b-c
1793 ;; Formula 36:
1795 ;; exp(-%i*m*%pi)*Q(n,m,z) =
1796 ;; 2^n*gamma(1+n)*gamma(1+n+m)*(z+1)^(m/2-n-1)*(z-1)^(-m/2)/gamma(2+2*n)
1797 ;; * hgfred([1+n-m,1+n],[2+2*n],2/(1+z))
1799 ;; Let a = 1+n-m, b = 1+n, c = 2+2*n. then n = b-1 and m = b - a.
1800 ;; (There are other solutions.)
1802 ;; F(a,b;c;z) = 2*gamma(2*b)/gamma(b)/gamma(2*b-a)*w^(-b)*(1-w)^((b-a)/2)
1803 ;; *Q(b-1,b-a,2/w-1)*exp(-%i*%pi*(b-a))
1805 (defun legf36 (arg-l1 arg-l2 arg)
1806 (declare (ignore arg-l2))
1807 (let* ((a (car arg-l1))
1808 (b (cadr arg-l1))
1809 (n (sub b 1))
1810 (m (sub b a))
1811 ;;z (div (sub 2 arg) arg)
1812 (z (sub (div 2 arg) 1)))
1813 (mul (inv (power 2 n))
1814 (inv (gm (add 1 n)))
1815 (inv (gm (add 1 n m)))
1816 (inv (power (add z 1)
1817 (add (div m 2)
1818 (mul -1 n)
1819 -1)))
1820 (inv (power (sub z 1) (div m -2)))
1821 (gm (add 2 n n))
1822 (power '$%e (mul -1 '$%i m '$%pi))
1823 (legen n m z '$q))))
1825 (defun legen (n m x pq)
1826 ;; A&S 8.2.1: P(-n-1,m,z) = P(n,m,z)
1827 (let ((n (if (or (member ($sign n) '($neg $nz)) ; negative sign or
1828 (mminusp n)) ; negative form like -n-1
1829 (mul -1 (add 1 n))
1830 n)))
1831 (cond ((equal m 0)
1832 (list (if (eq pq '$q)
1833 '($legendre_q simp)
1834 '($legendre_p simp))
1835 n x))
1837 (list (if (eq pq '$q)
1838 '($assoc_legendre_q simp)
1839 '($assoc_legendre_p simp))
1840 n m x)))))
1842 (defun legpol-core (a b c arg)
1843 ;; I think for this to be correct, we need a to be a negative integer.
1844 (unless (and (eq '$yes (ask-integerp a))
1845 (eq (asksign a) '$negative))
1846 (return-from legpol-core nil))
1847 (let* ((l (vfvp (div (add b a) 2) arg))
1848 (v (cdr (assoc 'v l :test #'equal))))
1849 ;; v is (a+b)/2
1850 (cond
1851 ((and (alike1 v '((rat simp) 1 2))
1852 (alike1 c 1))
1853 ;; A&S 22.5.49:
1854 ;; P(n,x) = F(-n,n+1;1;(1-x)/2)
1855 (legenpol (mul -1 a)
1856 (sub 1 (mul 2 arg))))
1858 ((and (alike1 c '((rat simp) 1 2))
1859 (alike1 (add b a) '((rat simp) 1 2)))
1860 ;; c = 1/2, a+b = 1/2
1862 ;; A&S 22.5.52
1863 ;; P(2*n,x) = (-1)^n*(2*n)!/2^(2*n)/(n!)^2*F(-n,n+1/2;1/2;x^2)
1865 ;; F(-n,n+1/2;1/2;x^2) = P(2*n,x)*(-1)^n*(n!)^2/(2*n)!*2^(2*n)
1867 (let ((n (mul -1 a)))
1868 (mul (power -1 n)
1869 (power (gm (add n 1)) 2)
1870 (inv (gm (add 1 (mul 2 n))))
1871 (power 2 (mul 2 n))
1872 (legenpol (mul 2 n)
1873 (power arg (div 1 2))))))
1875 ((and (alike1 c '((rat simp) 3 2))
1876 (alike1 (add b a) '((rat simp) 3 2)))
1877 ;; c = 3/2, a+b = 3/2
1879 ;; A&S 22.5.53
1880 ;; P(2*n+1,x) = (-1)^n*(2*n+1)!/2^(2*n)/(n!)^2*F(-n,n+3/2;3/2;x^2)*x
1882 ;; F(-n,n+3/2;3/2;x^2) = P(2*n+1,x)*(-1)^n*(n!)^2/(2*n+1)!*2^(2*n)/x
1884 (let ((n (mul -1 a)))
1885 (mul (power -1 n)
1886 (power (gm (add 1 n)) 2)
1887 (inv (gm (add 2 (mul 2 n))))
1888 (power 2 (mul 2 n))
1889 (legenpol (add 1 (mul 2 n))
1890 (power arg (div 1 2)))
1891 (inv (power arg (div 1 2))))))
1893 ((and (zerp (sub b a))
1894 (zerp (sub c (add a b))))
1895 ;; a = b, c = a + b
1897 ;; A&S 22.5.50
1898 ;; P(n,x) = binomial(2*n,n)*((x-1)/2)^n*F(-n,-n;-2*n;2/(1-x))
1900 ;; F(-n,-n;-2*n;x) = P(n,1-2/x)/binomial(2*n,n)(-1/x)^(-n)
1901 (mul (power (gm (add 1 (mul -1 a))) 2)
1902 (inv (gm (add 1 (mul -2 a))))
1903 (power (mul -1 arg) (mul -1 a))
1904 (legenpol (mul -1 a)
1905 (add 1 (div -2 arg)))))
1906 ((and (alike1 (sub a b) '((rat simp) 1 2))
1907 (alike1 (sub c (mul 2 b)) '((rat simp) 1 2)))
1908 ;; a - b = 1/2, c - 2*b = 1/2
1910 ;; A&S 22.5.51
1911 ;; P(n,x) = binomial(2*n,n)*(x/2)^n*F(-n/2,(1-n)/2;1/2-n;1/x^2)
1913 ;; F(-n/2,(1-n)/2;1/2-n,1/x^2) = P(n,x)/binomial(2*n,n)*(x/2)^(-n)
1914 (mul (power (gm (add 1 (mul -2 b))) 2)
1915 (inv (gm (add 1 (mul -4 b))))
1916 (power (mul 2 (power arg (div 1 2))) (mul -2 b))
1917 (legenpol (mul -2 b)
1918 (power arg (div -1 2)))))
1919 ((and (alike1 (sub b a) '((rat simp) 1 2))
1920 (alike1 (sub c (mul 2 a)) '((rat simp) 1 2)))
1921 ;; b - a = 1/2, c + 2*a = 1/2
1923 ;; A&S 22.5.51
1924 ;; P(n,x) = binomial(2*n,n)*(x/2)^n*F(-n/2,(1-n)/2;1/2-n;1/x^2)
1926 ;; F(-n/2,(1-n)/2;1/2-n,1/x^2) = P(n,x)/binomial(2*n,n)*(x/2)^(-n)
1927 (mul (power (gm (add 1 (mul -2 a))) 2)
1928 (inv (gm (add 1 (mul -4 a))))
1929 (power (mul 2 (power arg (div 1 2))) (mul -2 a))
1930 (legenpol (mul -2 a)
1931 (power arg (div -1 2)))))
1933 nil))))
1935 ;; This seems not be be called at all? There is a commented-out call
1936 ;; in LEGFUN that says it doesn't make sense to call LEGPOL there.
1937 #+nil
1938 (defun legpol (a b c)
1939 ;; See if F(a,b;c;z) is a Legendre polynomial. If not, try
1940 ;; F(b,a;c;z).
1941 (or (legpol-core a b c var)
1942 (legpol-core b a c var)))
1944 ;; See A&S 15.3.3:
1946 ;; F(a,b;c;z) = (1-z)^(c-a-b)*F(c-a,c-b;c;z)
1947 (defun gered1 (arg-l1 arg-l2 simpflg arg)
1948 (destructuring-bind (a b)
1949 arg-l1
1950 (destructuring-bind (c)
1951 arg-l2
1952 (mul (power (sub 1 arg)
1953 (add c
1954 (mul -1 a)
1955 (mul -1 b)))
1956 (funcall simpflg
1957 (list (sub c a)
1958 (sub c b))
1959 arg-l2
1960 arg)))))
1962 ;; See A&S 15.3.4
1964 ;; F(a,b;c;z) = (1-z)^(-a)*F(a,c-b;c;z/(z-1))
1965 (defun gered2 (a b c arg)
1966 (mul (power (sub 1 arg) (mul -1 a))
1967 (hgfsimp (list a (sub c b))
1968 (list c)
1969 (div arg (sub arg 1)))))
1971 ;; See A&S 15.3.9:
1973 ;; F(a,b;c;z) = A*z^(-a)*F(a,a-c+1;a+b-c+1;1-1/z)
1974 ;; + B*(1-z)^(c-a-b)*z^(a-c)*F(c-a,1-a;c-a-b+1,1-1/z)
1976 ;; where A = gamma(c)*gamma(c-a-b)/gamma(c-a)/gamma(c-b)
1977 ;; B = gamma(c)*gamma(a+b-c)/gamma(a)/gamma(b)
1979 (defun geredf (a b c arg)
1980 (let (($gamma_expand t))
1981 (add (div (mul (take '(%gamma) c)
1982 (take '(%gamma) (add c (mul -1 a) (mul -1 b)))
1983 (power arg (mul -1 a))
1984 ($hgfred `((mlist) ,a ,(add a 1 (mul -1 c)))
1985 `((mlist) ,(add a b (mul -1 c) 1))
1986 (sub 1 (div 1 arg))))
1987 (mul (take '(%gamma) (sub c a))
1988 (take '(%gamma) (sub c b))))
1989 (div (mul (take '(%gamma) c)
1990 (take '(%gamma) (add a b (mul -1 c)))
1991 (power (sub 1 arg)
1992 (add c (mul -1 a) (mul -1 b)))
1993 (power arg (sub a c))
1994 ($hgfred `((mlist) ,(sub c a) ,(sub 1 a))
1995 `((mlist) ,(add c (mul -1 a) (mul -1 b) 1))
1996 (sub 1 (div 1 arg))))
1997 (mul (take '(%gamma) a) (take '(%gamma) b))))))
1999 (defun trig-log (arg-l1 arg-l2 arg)
2000 (cond ((equal (simplifya (car arg-l2) nil) '((rat simp) 3 2))
2001 ;; c = 3/2
2002 (when $trace2f1
2003 (format t " trig-log: Test c=3/2~%"))
2004 (trig-log-3 arg-l1 arg-l2 arg))
2005 ((equal (simplifya (car arg-l2) nil) '((rat simp) 1 2))
2006 ;; c = 1/2
2007 (when $trace2f1
2008 (format t " trig-log: Test c=1/2~%"))
2009 (trig-log-1 arg-l1 arg-l2 arg))
2010 (t nil)))
2012 (defun trig-log-3 (arg-l1 arg-l2 arg)
2013 (cond ((and (or (equal (car arg-l1) 1) (equal (cadr arg-l1) 1))
2014 (or (equal (car arg-l1) (div 1 2))
2015 (equal (cadr arg-l1) (div 1 2))))
2016 ;; (a = 1 or b = 1) and (a = 1/2 or b = 1/2)
2017 (when $trace2f1
2018 (format t " Case a or b is 1 and the other is 1/2~%"))
2019 (trig-log-3-exec arg-l1 arg-l2 arg))
2020 ((and (equal (car arg-l1) (cadr arg-l1))
2021 (or (equal 1 (car arg-l1))
2022 (equal (div 1 2) (car arg-l1))))
2023 ;; a = b and (a = 1 or a = 1/2)
2024 (when $trace2f1
2025 (format t " Case a=b and a is 1 or 1/2~%"))
2026 (trig-log-3a-exec arg-l1 arg-l2 arg))
2027 ((or (equal (add (car arg-l1) (cadr arg-l1)) 1)
2028 (equal (add (car arg-l1) (cadr arg-l1)) 2))
2029 ;; a + b = 1 or a + b = 2
2030 (when $trace2f1
2031 (format t " Case a+b is 1 or 2~%"))
2032 (trig-sin arg-l1 arg-l2 arg))
2033 ((or (equal (sub (car arg-l1) (cadr arg-l1)) (div 1 2))
2034 (equal (sub (cadr arg-l1) (car arg-l1)) (div 1 2)))
2035 ;; a - b = 1/2 or b - a = 1/2
2036 (when $trace2f1
2037 (format t " Case a-b=1/2 or b-a=1/2~%"))
2038 (trig-3 arg-l1 arg-l2 arg))
2039 (t nil)))
2041 (defun trig-3 (arg-l1 arg-l2 arg)
2042 (declare (ignore arg-l2))
2043 ;; A&S 15.1.10
2045 ;; F(a,a+1/2,3/2,z^2) =
2046 ;; ((1+z)^(1-2*a) - (1-z)^(1-2*a))/2/z/(1-2*a)
2047 (let* (($radexpand '$all)
2048 (a (sub 1
2049 (sub (add (car arg-l1)
2050 (cadr arg-l1))
2051 (div 1 2))))
2052 (z (power arg (div 1 2))))
2053 (mul (inv z)
2054 (inv 2)
2055 (inv a)
2056 (sub (power (add 1 z) a)
2057 (power (sub 1 z) a)))))
2059 (defun trig-sin (arg-l1 arg-l2 arg)
2060 (declare (ignore arg-l2))
2061 ;; A&S 15.1.15, 15.1.16
2062 (destructuring-bind (a b)
2063 arg-l1
2064 ;; I think it's ok to convert sqrt(z^2) to z here, so $radexpand
2065 ;; is $all.
2066 (let (($radexpand '$all)
2067 a1 z1)
2068 (cond ((equal (add a b) 1)
2069 ;; A&S 15.1.15
2071 ;; F(a,1-a;3/2;sin(z)^2) =
2073 ;; sin((2*a-1)*z)/(2*a-1)/sin(z)
2074 (mul (inv (mul (mul -1 (sub a b))
2075 (msin (masin (msqrt arg)))))
2076 (msin (mul (mul -1
2077 (sub a b))
2078 (masin (msqrt arg))))))
2079 ((equal (add a b) 2)
2080 ;; A&S 15.1.16
2082 ;; F(a, 2-a; 3/2; sin(z)^2) =
2084 ;; sin((2*a-2)*z)/(a-1)/sin(2*z)
2085 (mul (msin (mul (setq z1
2086 (masin (msqrt
2087 arg)))
2088 (setq a1
2089 (mul -1
2090 (sub a
2091 b)))))
2092 (inv (mul a1
2093 (msin z1)
2094 (mcos z1)))))
2096 nil)))))
2099 ;;Generates atan if arg positive else log
2100 (defun trig-log-3-exec (arg-l1 arg-l2 arg)
2101 (declare (ignore arg-l1 arg-l2))
2102 ;; See A&S 15.1.4 and 15.1.5
2104 ;; F(a,b;3/2;z) where a = 1/2 and b = 1 (or vice versa).
2106 ;; I think it's ok to convert sqrt(z^2) to z here, so $radexpand is
2107 ;; $all.
2108 (let (($radexpand '$all))
2109 (cond ((equal (checksigntm arg) '$positive)
2110 ;; A&S 15.1.4
2112 ;; F(1/2,1;3/2,z^2) =
2114 ;; log((1+z)/(1-z))/z/2
2116 ;; This is the same as atanh(z)/z. Should we return that
2117 ;; instead? This would make this match what hyp-atanh
2118 ;; returns.
2119 (let ((z (power arg (div 1 2))))
2120 (mul (power z -1)
2121 (inv 2)
2122 (mlog (div (add 1 z)
2123 (sub 1 z))))))
2124 ((equal (checksigntm arg) '$negative)
2125 ;; A&S 15.1.5
2127 ;; F(1/2,1;3/2,z^2) =
2128 ;; atan(z)/z
2129 (let ((z (power (mul -1 arg)
2130 (div 1 2))))
2131 (mul (power z -1)
2132 (matan z)))))))
2134 (defun trig-log-3a-exec (arg-l1 arg-l2 arg)
2135 ;; See A&S 15.1.6 and 15.1.7
2137 ;; F(a,b;3/2,z) where a = b and a = 1/2 or a = 1.
2139 ;; I think it's ok to convert sqrt(z^2) to z here, so $radexpand is
2140 ;; $all.
2141 (let ((a (first arg-l1))
2142 ($radexpand '$all))
2143 (cond ((equal (checksigntm arg) '$positive)
2144 ;; A&S 15.1.6
2146 ;; F(1/2,1/2; 3/2; z^2) = sqrt(1-z^2)*F(1,1;3/2;z^2) =
2147 ;; asin(z)/z
2148 (let ((z (power arg (div 1 2))))
2149 (if (equal a 1)
2150 (div (trig-log-3a-exec (list (div 1 2) (div 1 2)) arg-l2 arg)
2151 (power (sub 1 (power z 2)) (div 1 2)))
2152 (div (masin z) z))))
2153 ((equal (checksigntm arg) '$negative)
2154 ;; A&S 15.1.7
2156 ;; F(1/2,1/2; 3/2; -z^2) = sqrt(1+z^2)*F(1,1,3/2; -z^2) =
2157 ;;log(z + sqrt(1+z^2))/z
2158 (let* ((z (power (mul -1 arg)
2159 (div 1 2)))
2160 (1+z^2 (add 1 (power z 2))))
2161 (if (equal a 1)
2162 (div (trig-log-3a-exec (list (div 1 2) (div 1 2))
2163 arg-l2
2164 arg)
2165 (power 1+z^2
2166 (div 1 2)))
2167 (div (mlog (add z (power 1+z^2
2168 (div 1 2))))
2169 z)))))))
2172 (defun trig-log-1 (arg-l1 arg-l2 arg) ;; 2F1's with C = 1/2
2173 (declare (ignore arg-l2))
2174 ;; 15.1.17, 11, 18, 12, 9, and 19
2176 ;; I think it's ok to convert sqrt(z^2) to z here, so $radexpand is
2177 ;; $all.
2178 (let (($radexpand '$all)
2179 x z $exponentialize a b)
2180 (setq a (car arg-l1) b (cadr arg-l1))
2181 (cond ((=0 (m+t a b))
2182 ;; F(-a,a;1/2,z)
2184 (cond ((equal (checksigntm arg) '$positive)
2185 ;; A&S 15.1.17:
2186 ;; F(-a,a;1/2;sin(z)^2) = cos(2*a*z)
2187 (trig-log-1-pos a arg))
2188 ((equal (checksigntm arg) '$negative)
2189 ;; A&X 15.1.11:
2190 ;; F(-a,a;1/2;-z^2) = 1/2*((sqrt(1+z^2)+z)^(2*a)
2191 ;; +(sqrt(1+z^2)-z)^(2*a))
2193 (trig-log-1-neg a b arg))
2194 (t ())))
2195 ((equal (m+t a b) 1.)
2196 ;; F(a,1-a;1/2,z)
2197 (cond ((equal (checksigntm arg) '$positive)
2198 ;; A&S 15.1.18:
2199 ;; F(a,1-a;1/2;sin(z)^2) = cos((2*a-1)*z)/cos(z)
2200 (m//t (mcos (m*t (m-t a b) (setq z (masin (msqrt arg)))))
2201 (mcos z)))
2202 ((equal (checksigntm arg) '$negative)
2203 ;; A&S 15.1.12
2204 ;; F(a,1-a;1/2;-z^2) = 1/2*(1+z^2)^(-1/2)*
2205 ;; {[(sqrt(1+z^2)+z]^(2*a-1)
2206 ;; +[sqrt(1+z^2)-z]^(2*a-1)}
2207 (m*t 1//2 (m//t (setq x (msqrt (m-t 1. arg))))
2208 (m+t (m^t (m+t x (setq z (msqrt (m-t arg))))
2209 (setq b (m-t a b)))
2210 (m^t (m-t x z) b))))
2211 (t ())))
2212 ((=1//2 (hyp-mabs (m-t b a)))
2213 ;; F(a, a+1/2; 1/2; z)
2214 (cond ((equal (checksigntm arg) '$positive)
2215 ;; A&S 15.1.9
2216 ;; F(a,1/2+a;1/2;z^2) = ((1+z)^(-2*a)+(1-z)^(-2*a))/2
2217 (m*t 1//2
2218 (m+t (m^t (m1+t (setq z (msqrt arg)))
2219 (setq b (m-t 1//2 (m+t a b))))
2220 (m^t (m-t 1. z) b))))
2221 ((equal (checksigntm arg) '$negative)
2222 ;; A&S 15.1.19
2223 ;; F(a,1/2+a;1/2;-tan(z)^2) = cos(z)^(2*a)*cos(2*a*z)
2224 (m*t (m^t (mcos (setq z (matan (msqrt (m-t arg)))))
2225 (setq b (m+t a b -1//2)))
2226 (mcos (m*t b z))))
2227 (t ())))
2228 (t ()))))
2230 (defun trig-log-1-pos (a z)
2231 ;; I think it's ok to convert sqrt(z^2) to z here, so $radexpand is
2232 ;; $all.
2233 (let (($radexpand '$all))
2234 (mcos (m*t 2. a (masin (msqrt z))))))
2236 (defun trig-log-1-neg (a b v)
2237 ;; Look to see a is of the form m*s+c where m and c
2238 ;; are numbers. If m is positive, swap a and b.
2239 ;; Basically we want F(-a,a;1/2;-z^2) =
2240 ;; F(a,-a;1/2;-z^2), as they should be.
2241 (let* (($radexpand '$all)
2242 (match (m*s+c a))
2243 (m (cdras 'm match))
2244 (s (cdras 's match))
2245 (b (if s
2246 (if (and m (eq (checksigntm m) '$positive))
2249 (if (eq (checksigntm a) '$negative)
2251 a)))
2252 (x (msqrt (m-t 1. v)))
2253 (z (msqrt (m-t v))))
2254 (m*t 1//2
2255 (m+t (m^t (m+t x z)
2256 (setq b (m*t 2. b)))
2257 (m^t (m-t x z) b)))))
2259 ;; Pattern match for m*s+c where a is a number, x is symbolic, and c
2260 ;; is a number.
2261 (defun m*s+c (exp)
2262 (m2 exp
2263 '((mplus) ((coeffpt) (m $numberp) (s nonnump))
2264 ((coeffpp) (c $numberp)))))
2266 ;; List L contains two elements first the numerator parameter that
2267 ;;exceeds the denumerator one and is called "C", second
2268 ;;the difference of the two parameters which is called "M".
2271 (defun diffintprop-gen-exec (l l1 l2)
2272 (prog (c m poly constfact )
2273 (setq c (car l)
2274 m (cadr l)
2275 l1 (delete c l1 :count 1 :test #'equal)
2276 c (sub c m)
2277 l2 (delete c l2 :count 1 :test equal)
2278 poly ($expand (constrpoly c m 'avgoustis))
2279 constfact (createconstfact c m))
2280 (return (yanmult constfact
2281 (diffintprop-exec poly l1 l2)))))
2283 (defun constrpoly (c m k)
2284 (cond ((zerop m) 1.)
2285 (t (mul (add c k (1- m)) (constrpoly c (1- m) k)))))
2287 (defun createconstfact (c m)
2288 (cond ((zerop m) 1.)
2289 (t (mul (inv (add c (1- m)))
2290 (createconstfact c (1- m))))))
2292 (defun diffintprop-exec (poly l1 l2)
2293 (distrdiffintprop (createcoefpowlist-exec poly) l1 l2))
2295 (defun distrdiffintprop (l l1 l2)
2296 (cond ((null l) 0.)
2297 (t (add (yanmult ($factor (cadar l))
2298 (diffintprop (caar l) l1 l2))
2299 (distrdiffintprop (cdr l) l1 l2)))))
2301 (defun diffintprop (pow l1 l2)
2302 (cond ((zerop pow) (hgfsimp l1 l2 var))
2303 ((equal pow 1.)
2304 (yanmult (mul (div (multpl l1) (multpl l2)) var)
2305 (hgfsimp (incr1 l1) (incr1 l2) var)))
2306 (t (searchaddserieslist pow l1 l2))))
2308 (defun searchaddserieslist (pow l1 l2)
2309 (prog (series res)
2310 (cond ((setq series (searchserieslistp serieslist pow))
2311 (return (eval series))))
2312 (setq
2313 serieslist
2314 (append
2315 serieslist
2316 (list
2317 (list
2319 (setq res
2320 '(yanmult (mul (div (multpl l1) (multpl l2))
2321 var)
2322 (diffintproprecurse (1- pow)
2323 (incr1 l1)
2324 (incr1 l2))))))))
2325 (return (eval res))))
2327 (defun diffintproprecurse (pow l1 l2)
2328 (prog (poly)
2329 (setq poly
2330 ($expand (power (add 'avgoustis 1.) pow)))
2331 (return (diffintprop-exec poly l1 l2))))
2333 (defun multpl (l)
2334 (cond ((null l) 1.) (t (mul (car l) (multpl (cdr l))))))
2336 (defun createcoefpowlist-exec (poly)
2337 (prog (hp conster)
2338 (setq conster (consterminit poly 'avgoustis)
2339 hp ($hipow poly 'avgoustis))
2340 (return (append (list (list 0. conster))
2341 (createcoefpowlist poly hp)))))
2343 (defun createcoefpowlist (poly hp)
2344 (cond ((equal hp 1.)
2345 (list (list 1. ($coeff poly 'avgoustis))))
2346 (t (append (createcoefpowlist poly (1- hp))
2347 (list (list hp
2348 ($coeff poly
2349 (power 'avgoustis
2350 hp))))))))
2352 (defun consterminit (fun var)
2353 (cond ((eq (caar fun) 'mplus) (consterm (cdr fun) var))
2354 (t (cond ((freevar fun) fun) (t 0.)))))
2356 (defun searchserieslistp (serieslist pow)
2357 (cond ((null serieslist) nil)
2358 ((equal (caar serieslist) pow) (cadar serieslist))
2359 (t (searchserieslistp (cdr serieslist) pow))))
2361 (defun yanmult (a b)
2362 (cond ((eq (caar b) 'mplus) (yanmul a (cdr b)))
2363 (t (mul a b))))
2365 (defun yanmul (a b)
2366 (cond ((null b) 0.)
2367 (t (add (mul a (car b)) (yanmul a (cdr b))))))
2371 (defun freevarpar (exp)
2372 (cond ((freevar exp) (freepar exp))
2373 (t nil)))
2375 (defun freevarpar2 (exp arg)
2376 (cond ((freevar2 exp arg)
2377 (freepar exp))
2378 (t nil)))
2380 ;; Why is this needed?
2381 (setq *par* '$p)
2383 (defun freepar (exp)
2384 (cond ((atom exp)
2385 (not (eq exp *par*)))
2386 (t (and (freepar (car exp))
2387 (freepar (cdr exp))))))
2389 ;; Confluent hypergeometric function.
2391 ;; F(a;c;z)
2392 (defun confl (arg-l1 arg-l2 arg)
2393 (let* ((a (car arg-l1))
2394 (c (car arg-l2))
2395 (a-c (sub a c))
2397 (cond ((zerop1 arg)
2398 ;; F(a;c;0) = 1
2399 (add 1 arg))
2401 ((and (equal 1 c)
2402 (not (integerp a)) ; Do not handle an integer or
2403 (not (integerp (add a a)))) ; an half integral value
2404 ;; F(a;1;z) = laguerre(-a,z)
2405 (lagpol (neg a) 0 arg))
2407 ((and (maxima-integerp a)
2408 (member ($sign a) '($neg nz)))
2409 ;; F(-n; c; z) and n a positive integer
2410 (1f1polys (list c) a arg))
2412 ((alike1 c (add a a))
2413 ;; F(a;2a;z)
2414 ;; A&S 13.6.6
2416 ;; F(n+1;2*n+1;2*z) =
2417 ;; gamma(3/2+n)*exp(z)*(z/2)^(-n-1/2)*bessel_i(n+1/2,z).
2419 ;; So
2421 ;; F(n,2*n,z) =
2422 ;; gamma(n+1/2)*exp(z/2)*(z/4)^(-n-3/2)*bessel_i(n-1/2,z/2);
2424 ;; If n is a negative integer, we use laguerre polynomial.
2425 (if (and (maxima-integerp a)
2426 (eq (asksign a) '$negative))
2427 ;; We have already checked a negative integer. This algorithm
2428 ;; is now present in 1f1polys and at this place never called.
2429 (let ((n (neg a)))
2430 (mul (power -1 n)
2431 (inv (take '(%binomial) (add n n) n))
2432 (lagpol n (sub c 1) arg)))
2433 (let ((z (div arg 2)))
2434 (mul (power '$%e z)
2435 (bestrig (add a '((rat simp) 1 2))
2436 (div (mul z z) 4))))))
2438 ((and (integerp (setq n (sub (add a a) c)))
2439 (plusp n)
2440 (not (integerp a))
2441 (not (integerp (add a a))))
2442 ;; F(a,2*a-n,z) and a not an integer or a half integral value
2443 (when *debug-hyp*
2444 (format t "~&Case 1F1(a,2*a-n,x):")
2445 (format t "~& ; a = ~A~%" a)
2446 (format t "~& ; c = ~A~%" c)
2447 (format t "~& : n = ~A~%" n))
2448 (sratsimp
2449 (mul (take '(%gamma) (sub a (add n '((rat simp) 1 2))))
2450 (power (div arg 4)
2451 (sub (add '((rat simp) 1 2) n) a))
2452 (power '$%e (div arg 2))
2453 (let ((index (gensym)))
2454 (dosum
2455 (mul (power -1 index)
2456 (take '($pochhammer) (- n) index)
2457 (take '($pochhammer) (add a a (* -2 n) -1) index)
2458 (add a index (- n) '((rat simp) -1 2))
2459 (inv (take '($pochhammer) (sub (add a a) n) index))
2460 (inv (take '(mfactorial) index))
2461 (take '(%bessel_i)
2462 (add a index (- n) '((rat simp) -1 2))
2463 (div arg 2)))
2464 index 0 n t)))))
2466 ((and (integerp (setq n (sub c (add a a))))
2467 (plusp n)
2468 (not (integerp a))
2469 (not (integerp (add a a))))
2470 ;; F(a,2*a+n,z) and a not an integer or a half integral value
2471 (when *debug-hyp*
2472 (format t "~&Case 1F1(a,2*a+n,x):")
2473 (format t "~& ; a = ~A~%" a)
2474 (format t "~& ; c = ~A~%" c)
2475 (format t "~& : n = ~A~%" n))
2476 (sratsimp
2477 (mul (take '(%gamma) (sub a '((rat simp) 1 2)))
2478 (power (div arg 4)
2479 (sub '((rat simp) 1 2) a))
2480 (power '$%e (div arg 2))
2481 (let ((index (gensym)))
2482 (dosum
2483 (mul (take '($pochhammer) (- n) index)
2484 (take '($pochhammer) (add a a -1) index)
2485 (add a index '((rat simp) -1 2))
2486 (inv (take '($pochhammer) (add a a n) index))
2487 (inv (take '(mfactorial) index))
2488 (take '(%bessel_i)
2489 (add a index '((rat simp) -1 2))
2490 (div arg 2)))
2491 index 0 n t)))))
2493 ((and (integerp (setq n (sub a '((rat simp) 1 2))))
2494 (>= n 0)
2495 (integerp c)
2496 (plusp c))
2497 (let ((m c)
2498 ($simpsum t)
2499 ($bessel_reduce t))
2500 (when *debug-hyp*
2501 (format t "~&Case 1F1(n+1/2,m,x):")
2502 (format t "~& ; a = ~A~%" a)
2503 (format t "~& ; c = ~A~%" c)
2504 (format t "~& : n = ~A~%" n)
2505 (format t "~& : m = ~A~%" m))
2506 (sratsimp
2507 (mul (power 2 (- 1 m))
2508 (power '$%e (div arg 2))
2509 (factorial (- m 1))
2510 (factorial n)
2511 (inv (take '($pochhammer) '((rat simp) 1 2) (- m 1)))
2512 (inv (take '($pochhammer) '((rat simp) 1 2) n))
2513 (let ((index1 (gensumindex)))
2514 (dosum
2515 (mul (power 2 (neg index1))
2516 (power (neg arg) index1)
2517 (inv (take '(mfactorial) index1))
2518 (mfuncall '$gen_laguerre
2519 (sub n index1)
2520 (sub index1 '((rat simp) 1 2))
2521 (neg arg))
2522 (let ((index2 (gensumindex)))
2523 (dosum
2524 (mul (power -1 index2)
2525 (power 2 (neg index2))
2526 (take '(%binomial)
2527 (add index1 m -1)
2528 index2)
2529 (let ((index3 (gensumindex)))
2530 (dosum
2531 (mul (take '(%binomial) index2 index3)
2532 (take '(%bessel_i)
2533 (sub index2 (mul 2 index3))
2534 (div arg 2)))
2535 index3 0 index2 t)))
2536 index2 0 (add index1 m -1) t)))
2537 index1 0 n t))))))
2539 ((and (integerp (setq n (sub a '((rat simp) 1 2))))
2540 (< n 0)
2541 (integerp c)
2542 (plusp c))
2543 (let ((n (- n))
2544 (m c)
2545 ($simpsum t)
2546 ($bessel_reduce t))
2547 (when *debug-hyp*
2548 (format t "~&Case 1F1(1/2-n,m,x):")
2549 (format t "~& ; a = ~A~%" a)
2550 (format t "~& ; c = ~A~%" c)
2551 (format t "~& : n = ~A~%" n)
2552 (format t "~& : m = ~A~%" m))
2553 (sratsimp
2554 (mul (power -1 n)
2555 (power 2 (- 1 m))
2556 (power '$%e (div arg 2))
2557 (factorial (- m 1))
2558 (inv (take '($pochhammer) '((rat simp) 1 2) (- m 1)))
2559 (inv (take '($pochhammer) (sub m '((rat simp) 1 2)) n))
2560 (let ((index1 (gensumindex)))
2561 (dosum
2562 (mul (power 2 (neg index1))
2563 (power arg index1)
2564 (take '(%binomial) n index1)
2565 (take '($pochhammer)
2566 (sub '((rat simp) 3 2) (+ m n))
2567 (sub n index1))
2568 (let ((index2 (gensumindex)))
2569 (dosum
2570 (mul (power '((rat simp) -1 2) index2)
2571 (take '(%binomial)
2572 (add index1 m -1)
2573 index2)
2574 (let ((index3 (gensumindex)))
2575 (dosum
2576 (mul (take '(%binomial) index2 index3)
2577 (take '(%bessel_i)
2578 (sub index2 (mul 2 index3))
2579 (div arg 2)))
2580 index3 0 index2 t)))
2581 index2 0 (add index1 m -1) t)))
2582 index1 0 n t))))))
2584 ((not (hyp-integerp a-c))
2585 (cond ((hyp-integerp a)
2586 (kummer arg-l1 arg-l2 arg))
2587 ($prefer_whittaker
2588 ;; A&S 15.1.32:
2590 ;; %m[k,u](z) = exp(-z/2)*z^(u+1/2)*M(1/2+u-k,1+2*u,z)
2592 ;; So
2594 ;; M(a,c,z) = exp(z/2)*z^(-c/2)*%m[c/2-a,c/2-1/2](z)
2596 ;; But if we apply Kummer's transformation, we can also
2597 ;; derive the expression
2599 ;; %m[k,u](z) = exp(z/2)*z^(u+1/2)*M(1/2+u+k,1+2*u,-z)
2601 ;; or
2603 ;; M(a,c,z) = exp(-z/2)*(-z)^(-c/2)*%m[a-c/2,c/2-1/2](-z)
2605 ;; For Laplace transforms it might be more beneficial to
2606 ;; return this second form instead.
2607 (let* ((m (div (sub c 1) 2))
2608 (k (add '((rat simp) 1 2) m (mul -1 a))))
2609 (mul (power arg (mul -1 (add '((rat simp) 1 2) m)))
2610 (power '$%e (div arg 2))
2611 (whitfun k m arg))))
2613 (fpqform arg-l1 arg-l2 arg))))
2614 ((minusp a-c)
2615 (sratsimp (erfgammared a c arg)))
2617 (kummer arg-l1 arg-l2 arg)))))
2619 ;; A&S 13.6.19:
2620 ;; M(1/2,3/2,-z^2) = sqrt(%pi)*erf(z)/2/sqrt(z)
2622 ;; So M(1/2,3/2,z) = sqrt(%pi)*erf(sqrt(-z))/2/sqrt(-z)
2623 ;; = sqrt(%pi)*erf(%i*sqrt(z))/2/(%i*sqrt(z))
2624 (defun hyprederf (x)
2625 (let ((x (mul '$%i (power x '((rat simp) 1 2 )))))
2626 (mul (power '$%pi '((rat simp) 1 2))
2627 '((rat simp) 1 2)
2628 (inv x)
2629 (take '(%erf) x))))
2631 ;; M(a,c,z), where a-c is a negative integer.
2632 (defun erfgammared (a c z)
2633 (cond ((and (nump a) (nump c))
2634 (erfgamnumred a c z))
2635 (t (gammareds a c z))))
2637 ;; I (rtoy) think this is what the function above is doing, but I'm
2638 ;; not sure. Plus, I think it's wrong.
2640 ;; For hgfred([n],[2+n],-z), the above returns
2642 ;; 2*n*(n+1)*z^(-n-1)*(gamma_incomplete_lower(n,z)*z-gamma_incomplete_lower(n+1,z))
2644 ;; But from A&S 13.4.3
2646 ;; -M(n,2+n,z) - n*M(n+1,n+2,z) + (n+1)*M(n,n+1,z) = 0
2648 ;; so M(n,2+n,z) = (n+1)*M(n,n+1,z)-n*M(n+1,n+2,z)
2650 ;; And M(n,n+1,-z) = n*z^(-n)*gamma_incomplete_lower(n,z)
2652 ;; This gives
2654 ;; M(n,2+n,z) = (n+1)*n*z^(-n)*gamma_incomplete_lower(n,z) - n*(n+1)*z^(-n-1)*gamma_incomplete_lower(n+1,z)
2655 ;; = n*(n+1)*z^(-n-1)*(gamma_incomplete_lower(n,z)*n-gamma_incomplete_lower(n+1,z))
2657 ;; So the version above is off by a factor of 2. But I think it's more than that.
2658 ;; Using A&S 13.4.3 again,
2660 ;; M(n,n+3,-z) = [n*M(n+1,n+3,-z) - (n+2)*M(n,n+2,-z)]/(-2);
2662 ;; The version above doesn't produce anything like this equation would
2663 ;; produce, given the value of M(n,n+2,-z) derived above.
2664 (defun gammareds (a c z)
2665 ;; M(a,c,z) where a-c is a negative integer.
2666 (let ((diff (sub c a)))
2667 (cond ((eql diff 1)
2668 ;; We have M(a,a+1,z).
2669 (hypredincgm a z))
2670 ((eql a 1)
2671 ;; We have M(1,a,z)
2672 ;; Apply Kummer's transformation to get the form M(a-1,a,z)
2674 ;; (I don't think we ever get here, but just in case, we leave it.)
2675 (kummer (list a) (list c) z))
2677 ;; We have M(a, a+n, z)
2679 ;; A&S 13.4.3 says
2680 ;; (1+a-b)*M(a,b,z) - a*M(a+1,b,z)+(b-1)*M(a,b-1,z) = 0
2682 ;; So
2684 ;; M(a,b,z) = [a*M(a+1,b,z) - (b-1)*M(a,b-1,z)]/(1+a-b);
2686 ;; Thus, the difference between b and a is reduced, until
2687 ;; b-a=1, which we handle above.
2688 (mul (sub (mul a
2689 (gammareds (add 1 a) c z))
2690 (mul (sub c 1)
2691 (gammareds a (sub c 1) z)))
2692 (inv (sub (add 1 a) c)))))))
2694 ;; A&S 6.5.12:
2695 ;; gamma_incomplete_lower(a,x) = x^a/a*M(a,1+a,-x)
2696 ;; = x^a/a*exp(-x)*M(1,1+a,x)
2698 ;; where gamma_incomplete_lower(a,x) is the lower incomplete gamma function.
2700 ;; M(a,1+a,x) = a*(-x)^(-a)*gamma_incomplete_lower(a,-x)
2701 (defun hypredincgm (a z)
2702 (let ((-z (mul -1 z)))
2703 (if (not $prefer_gamma_incomplete)
2704 (mul a
2705 (power -z (mul -1 a))
2706 (take '(%gamma_incomplete_lower) a -z))
2707 (mul a
2708 (power -z (mul -1 a))
2709 (sub (take '(%gamma) a)
2710 (take '(%gamma_incomplete) a -z))))))
2712 ;; M(a,c,z), when a and c are numbers, and a-c is a negative integer
2713 (defun erfgamnumred (a c z)
2714 (cond ((hyp-integerp (sub c '((rat simp) 1 2)))
2715 (erfred a c z))
2716 (t (gammareds a c z))))
2718 ;; M(a,c,z) when a and c are numbers and c-1/2 is an integer and a-c
2719 ;; is a negative integer. Thus, we have M(p+1/2, q+1/2,z)
2720 (defun erfred (a c z)
2721 (prog (n m)
2722 (setq n (sub a '((rat simp) 1 2))
2723 m (sub c '((rat simp) 3 2)))
2724 ;; a = n + 1/2
2725 ;; c = m + 3/2
2726 ;; a - c < 0 so n - m - 1 < 0
2727 (cond ((not (or (> n m) (minusp n)))
2728 ;; 0 <= n <= m
2729 (return (thno33 n m z)))
2730 ((and (minusp n) (minusp m))
2731 ;; n < 0 and m < 0
2732 (return (thno35 (mul -1 n) (mul -1 m) z)))
2733 ((and (minusp n) (plusp m))
2734 ;; n < 0 and m > 0
2735 (return (thno34 (mul -1 n) m z)))
2737 ;; n = 0 or m = 0
2738 (return (gammareds (add n '((rat simp) 1 2))
2739 (add m '((rat simp) 3 2))
2740 z))))))
2742 ;; Compute M(n+1/2, m+3/2, z) with 0 <= n <= m.
2744 ;; I (rtoy) think this is what this routine is doing. (I'm guessing
2745 ;; that thno33 means theorem number 33 from Yannis Avgoustis' thesis.)
2747 ;; I don't have his thesis, but I see there are similar ways to derive
2748 ;; the result we want.
2750 ;; Method 1:
2751 ;; Use Kummer's transformation (A&S ) to get
2753 ;; M(n+1/2,m+3/2,z) = exp(z)*M(m-n+1,m+3/2,-z)
2755 ;; From A&S, we have
2757 ;; diff(M(1,n+3/2,z),z,m-n) = poch(1,m-n)/poch(n+3/2,m-n)*M(m-n+1,m+3/2,z)
2759 ;; Apply Kummer's transformation again:
2761 ;; M(1,n+3/2,z) = exp(z)*M(n+1/2,n+3/2,-z)
2763 ;; Apply the differentiation formula again:
2765 ;; diff(M(1/2,3/2,z),z,n) = poch(1/2,n)/poch(3/2,n)*M(n+1/2,n+3/2,z)
2767 ;; And we know that M(1/2,3/2,z) can be expressed in terms of erf.
2769 ;; Method 2:
2771 ;; Since n <= m, apply the differentiation formula:
2773 ;; diff(M(1/2,m-n+3/2,z),z,n) = poch(1/2,n)/poch(m-n+3/2,n)*M(n+1/2,m+3/2,z)
2775 ;; Apply Kummer's transformation:
2777 ;; M(1/2,m-n+3/2,z) = exp(z)*M(m-n+1,m-n+3/2,z)
2779 ;; Apply the differentiation formula again:
2781 ;; diff(M(1,3/2,z),z,m-n) = poch(1,m-n)/poch(3/2,m-n)*M(m-n+1,m-n+3/2,z)
2783 ;; I think this routine uses Method 2.
2784 (defun thno33 (n m x)
2785 ;; M(n+1/2,m+3/2,z) = diff(M(1/2,m-n+3/2,z),z,n)*poch(m-n+3/2,n)/poch(1/2,n)
2786 ;; M(1/2,m-n+3/2,z) = exp(z)*M(m-n+1,m-n+3/2,-z)
2787 ;; M(m-n+1,m-n+3/2,z) = diff(M(1,3/2,z),z,m-n)*poch(3/2,m-n)/poch(1,m-n)
2788 ;; diff(M(1,3/2,z),z,m-n) = (-1)^(m-n)*diff(M(1,3/2,-z),z,m-n)
2789 ;; M(1,3/2,-z) = exp(-z)*M(1/2,3/2,z)
2790 (let* ((yannis (gensym))
2791 (m-n (sub m n))
2792 ;; poch(m-n+3/2,n)/poch(1/2,n)
2793 (factor1 (div (take '($pochhammer) (add m-n '((rat simp) 3 2)) n)
2794 (take '($pochhammer) '((rat simp) 1 2) n)))
2795 ;; poch(3/2,m-n)/poch(1,m-n)
2796 (factor2 (div (take '($pochhammer) '((rat simp) 3 2) m-n)
2797 (take '($pochhammer) 1 m-n)))
2798 ;; M(1,3/2,-z) = exp(-z)*M(1/2,3/2,z)
2799 (hgferf (mul (power '$%e (mul -1 yannis))
2800 (hyprederf yannis)))
2801 ;; diff(M(1,3/2,z),z,m-n)
2802 (diff1 (meval (list '($diff) hgferf yannis m-n)))
2803 ;; exp(z)*M(m-n+1,m-n+3/2,-z)
2804 (kummer (mul (power '$%e yannis) diff1))
2805 ;; diff(M(1/2,m-n+3/2,z),z,n)
2806 (diff2 (meval (list '($diff) kummer yannis n))))
2807 ;; Multiply all the terms together.
2808 (sratsimp (mul (power -1 m-n)
2809 factor1
2810 factor2
2811 (maxima-substitute x yannis diff2)))))
2813 ;; M(n+1/2,m+3/2,z), with n < 0 and m > 0
2815 ;; Let's write it more explicitly as M(-n+1/2,m+3/2,z) with n > 0 and
2816 ;; m > 0.
2818 ;; First, use Kummer's transformation to get
2820 ;; M(-n+1/2,m+3/2,z) = exp(z)*M(m+n+1,m+3/2,-z)
2822 ;; We also have
2824 ;; diff(z^(n+m)*M(m+1,m+3/2,z),z,n) = poch(m+1,n)*z^m*M(m+n+1,m+3/2,z)
2826 ;; And finally
2828 ;; diff(M(1,3/2,z),z,m) = poch(1,m)/poch(3/2,m)*M(m+1,m+3/2,z)
2830 ;; Thus, we can compute M(-n+1/2,m+3/2,z) from M(1,3/2,z).
2832 ;; The second formula above can be derived easily by multiplying the
2833 ;; series for M(m+1,m+3/2,z) by z^(n+m) and differentiating n times.
2835 (defun thno34 (n m x)
2836 (let ((yannis (gensym)))
2837 (sratsimp
2838 (maxima-substitute
2840 yannis
2841 (mul (power -1 m)
2842 (div (mul (take '($pochhammer) '((rat simp) 3 2) m)
2843 (power '$%e yannis))
2844 (mul (take '($pochhammer) 1 m)
2845 (take '($pochhammer) (1+ m) n)
2846 (power yannis m)))
2847 (meval (list '($diff)
2848 (mul (power yannis (+ m n))
2849 (meval (list '($diff)
2850 (mul (power '$%e
2851 (mul -1 yannis))
2852 (hyprederf yannis))
2853 yannis
2854 m)))
2855 yannis
2856 n)))))))
2858 ;; M(n+1/2,m+3/2,z), with n < 0 and m < 0
2860 ;; Write it more explicitly as M(-n+1/2,-m+3/2,z) with n > 0 and m >
2861 ;; 0.
2863 ;; We know that
2865 ;; diff(sqrt(z)*M(-n+1/2,3/2,z),z,m) = poch(3/2-m,m)*M(-n+1/2,-m+3/2,z).
2867 ;; Apply Kummer's transformation:
2869 ;; M(-n+1/2,3/2,z) = exp(z) * M(n+1,3/2,-z)
2871 ;; Finally
2873 ;; diff(z^n*M(1,3/2,z),z,n) = n!*M(n+1,3/2,z)
2875 ;; So we can express M(-n+1/2,-m+3/2,z) in terms of M(1,3/2,z).
2877 ;; The first formula above follows from the more general formula
2879 ;; diff(z^(b-1)*M(a,b,z),z,n) = poch(b-n,n)*z^(b-n-1)*M(a,b-n,z)
2881 ;; The last formula follows from the general result
2883 ;; diff(z^(a+n-1)*M(a,b,z),z,n) = poch(a,n)*z^(a-1)*M(a+n,b,z)
2885 ;; Both of these are easily derived by using the series for M and
2886 ;; differentiating.
2888 (defun thno35 (n m x)
2889 (let ((yannis (gensym)))
2890 (sratsimp
2891 (maxima-substitute
2893 yannis
2894 (mul (div (power yannis (sub m '((rat simp) 1 2)))
2895 (mul (power -1 (* -1 m))
2896 (take '($pochhammer) 1 n)
2897 (take '($pochhammer) '((rat simp) -1 2) m)))
2898 (meval (list '($diff)
2899 (mul (power yannis '((rat simp) 1 2))
2900 (power '$%e yannis)
2901 (meval (list '($diff)
2902 (mul (power '$%e
2903 (mul -1 yannis))
2904 (power yannis n)
2905 (hyprederf yannis))
2906 yannis
2907 n)))
2908 yannis
2909 m)))))))
2911 ;; Pochhammer symbol. fctrl(a,n) = a*(a+1)*(a+2)*...*(a+n-1).
2913 ;; N must be a positive integer!
2915 ;; FIXME: This appears to be identical to factf below.
2916 (defun fctrl (a n)
2917 (cond ((zerop n)
2919 ((equal n 1)
2922 (mul (add a (1- n))
2923 (fctrl a (1- n))))))
2925 (setq *par* '$p)
2927 (defun vfvp (exp arg)
2928 (m2 exp `(v freevarpar2 ,arg)))
2931 (defun fpqform (arg-l1 arg-l2 arg)
2932 (list '(mqapply)
2933 (list '($%f simp array) (length arg-l1)(length arg-l2))
2934 (append (list '(mlist simp)) arg-l1)
2935 (append (list '(mlist simp)) arg-l2)
2936 arg))
2938 ;; Consider pFq([a_k]; [c_j]; z). If a_k = c_j + m for some k and j
2939 ;; and m >= 0, we can express pFq in terms of (p-1)F(q-1).
2941 ;; Here is a derivation for F(a,b;c;z), but it generalizes to the
2942 ;; generalized hypergeometric very easily.
2944 ;; From A&s 15.2.3:
2946 ;; diff(z^(a+n-1)*F(a,b;c;z), z, n) = poch(a,n)*z^(a-1)*F(a+n,b;c;z)
2948 ;; F(a+n,b;c;z) = diff(z^(a+n-1)*F(a,b;c;z), z, n)/poch(a,n)/z^(a-1)
2951 ;; So this expresses F(a+n,b;c;z) in terms of F(a,b;c;z). Let a = c +
2952 ;; n. This therefore gives F(c+n,b;c;z) in terms of F(c,b;c;z) =
2953 ;; 1F0(b;;z), which we know.
2955 ;; For simplicity, we will write F(z) for F(a,b;c;z).
2957 ;; Now,
2959 ;; n
2960 ;; diff(z^x*F(z),z,n) = sum binomial(n,k)*diff(z^x,z,n-k)*diff(F(z),z,k)
2961 ;; k=0
2963 ;; But diff(z^x,z,n-k) = x*(x-1)*...*(x-n+k+1)*z^(x-n+k)
2964 ;; = poch(x-n+k+1,n-k)*z^(x-n+k)
2966 ;; so
2968 ;; z^(-a+1)/poch(a,n)*diff(z^(a+n-1),z,n-k)
2969 ;; = poch(a+n-1-n+k+1,n-k)/poch(a,n)*z^(a+n-1-n+k)*z^(-a+1)
2970 ;; = poch(a+k,n-k)/poch(a,n)*z^k
2971 ;; = z^k/poch(a,k)
2973 ;; Combining these we have
2975 ;; n
2976 ;; F(a+n,b;c;z) = sum z^k/poch(a,k)*binomial(n,k)*diff(F(a,b;c;z),z,k)
2977 ;; k=0
2979 ;; Since a = c, we have
2981 ;; n
2982 ;; F(a+n,b;a;z) = sum z^k/poch(a,k)*binomial(n,k)*diff(F(a,b;a;z),z,k)
2983 ;; k=0
2985 ;; But F(a,b;a;z) = F(b;;z) and it's easy to see that A&S 15.2.2 can
2986 ;; be specialized to this case to give
2988 ;; diff(F(b;;z),z,k) = poch(b,k)*F(b+k;;z)
2990 ;; Finally, combining all of these, we have
2992 ;; n
2993 ;; F(a+n,b;c;z) = sum z^k/poch(a,k)*binomial(n,k)*poch(b,k)*F(b+k;;z)
2994 ;; k=0
2996 ;; Thus, F(a+n,b;c;z) is expressed in terms of 1F0(b+k;;z), as desired.
2997 (defun splitpfq (l arg-l1 arg-l2 arg)
2998 (destructuring-bind (a1 k)
3000 (let* ((result 0)
3001 (prodnum 1)
3002 (proden 1)
3003 (b1 (sub a1 k))
3004 (prod-b 1)
3005 (arg-l1 (delete a1 arg-l1 :count 1 :test #'equal))
3006 (arg-l2 (delete b1 arg-l2 :count 1 :test #'equal)))
3007 (loop for count from 0 upto k
3009 (when $trace2f1
3010 (format t "splitpfg term:~%")
3011 (maxima-display (mul (combin k count)
3012 (div prodnum proden)
3013 (inv prod-b)
3014 (power arg count)))
3015 (mtell "F(~:M, ~:M)~%"
3016 (cons '(mlist) arg-l1)
3017 (cons '(mlist) arg-l2)))
3018 (setq result (add result
3019 (mul (combin k count)
3020 (div prodnum proden)
3021 (inv prod-b)
3022 (power arg count)
3023 (hgfsimp arg-l1 arg-l2 arg))))
3024 (setq prod-b (mul prod-b (add b1 count)))
3025 (setq prodnum (mul prodnum (mull arg-l1))
3026 proden (mul proden (mull arg-l2)))
3027 (setq arg-l1 (incr1 arg-l1))
3028 (setq arg-l2 (incr1 arg-l2)))
3029 result)))
3031 ;; binomial(k,count)
3032 (defun combin (k count)
3033 (div (factorial k)
3034 (mul (factorial count)
3035 (factorial (sub k count)))))
3038 ;; We have something like F(s+m,-s+n;c;z)
3039 ;; Rewrite it like F(a'+d,-a';c;z) where a'=s-n=-b and d=m+n.
3041 (defun algii (a b)
3042 (let* ((sym (cdras 'f (s+c a)))
3043 (sign (cdras 'm (m2 sym '((mtimes) ((coefft) (m $numberp)) ((coefft) (s nonnump)))))))
3044 (when (and sign (minusp sign))
3045 (rotatef a b))
3046 (list nil (mul -1 b) (add a b))))
3049 ;;Algor. 2F1-RL from thesis:step 4:dispatch on a+m,-a+n,1/2+l cases
3050 (defun step4 (a b c arg)
3051 ;; F(a,b;c;z) where a+b is an integer and c+1/2 is an integer. If a
3052 ;; and b are not integers themselves, we can derive the result from
3053 ;; F(a1,-a1;1/2;z). However, if a and b are integers, we can't use
3054 ;; that because F(a1,-a1;1/2;z) is a polynomial. We need to derive
3055 ;; the result from F(1,1;3/2;z).
3056 (if (and (hyp-integerp a)
3057 (hyp-integerp b))
3058 (step4-int a b c arg)
3059 (step4-a a b c arg)))
3061 (defun step4-a (a b c arg)
3062 (let* ((alglist (algii a b))
3063 (aprime (cadr alglist))
3064 (m (caddr alglist))
3065 (n (sub c (inv 2)))
3066 ($ratsimpexpons $true)
3067 ($ratprint $false))
3068 ;; At this point, we have F(a'+m,-a';1/2+n;z) where m and n are
3069 ;; integers.
3070 (cond ((hyp-integerp (add aprime (inv 2)))
3071 ;; Ok. We have a problem if aprime + 1/2 is an integer.
3072 ;; We can't always use the algorithm below because we have
3073 ;; F(1/2,-1/2;1/2;z) which is 1F0(-1/2;;z) so the
3074 ;; derivation isn't quite right. Also, sometimes we'll end
3075 ;; up with a division by zero.
3077 ;; Thus, We need to do something else. So, use A&S 15.3.3
3078 ;; to change the problem:
3080 ;; F(a,b;c;z) = (1-z)^(c-a-b)*F(c-a, c-b; c; z)
3082 ;; which is
3084 ;; F('a+m,-a';1/2+n;z) = (1-z)^(1/2+n-m)*F(1/2+n-a'-m,1/2+n+a';1/2+n;z)
3086 ;; Recall that a' + 1/2 is an integer. Thus we have
3087 ;; F(<int>,<int>,1/2+n;z), which we know how to handle in
3088 ;; step4-int.
3089 (gered1 (list a b) (list c) #'hgfsimp arg))
3091 (let ((newf
3092 (cond ((equal (checksigntm arg) '$positive)
3093 (trig-log-1-pos aprime 'ell))
3094 ((equal (checksigntm arg) '$negative)
3095 (trig-log-1-neg (mul -1 aprime) aprime 'ell)))))
3096 ;; Ok, this uses F(a,-a;1/2;z). Since there are 2 possible
3097 ;; representations (A&S 15.1.11 and 15.1.17), we check the sign
3098 ;; of the arg (as done in trig-log-1) to select which form we
3099 ;; want to use. The original didn't and seemed to want to use
3100 ;; the negative form.
3102 ;; With this change, F(a,-a;3/2;z) matches what A&S 15.2.6 would
3103 ;; produce starting from F(a,-a;1/2;z), assuming z < 0.
3105 (subst arg 'ell
3106 (algiii newf
3107 m n aprime)))))))
3109 ;; F(a,b;c;z), where a and b are (positive) integers and c = 1/2+l.
3110 ;; This can be computed from F(1,1;3/2;z).
3112 ;; Assume a < b, without loss of generality.
3114 ;; F(m,n;3/2+L;z), m < n.
3116 ;; We start from F(1,1;3/2;z). Use A&S 15.2.3, differentiating m
3117 ;; times to get F(m,1;3/2;z). Swap a and b to get F(m,1;3/2;z) =
3118 ;; F(1,m;3/2;z) and use A&S 15.2.3 again to get F(n,m;3/2;z) by
3119 ;; differentiating n times. Finally, if L < 0, use A&S 15.2.4.
3120 ;; Otherwise use A&S 15.2.6.
3122 ;; I (rtoy) can't think of any way to do this with less than 3
3123 ;; differentiations.
3125 ;; Note that if d = (n-m)/2 is not an integer, we can write F(m,n;c;z)
3126 ;; as F(-d+u,d+u;c;z) where u = (n+m)/2. In this case, we could use
3127 ;; step4-a to compute the result.
3130 ;; Transform F(a,b;c;z) to F(a+n,b;c;z), given F(a,b;c;z)
3131 (defun as-15.2.3 (a bb cx n arg fun)
3132 (declare (ignore bb cx))
3133 (assert (>= n 0))
3134 ;; A&S 15.2.3:
3135 ;; F(a+n,b;c;z) = z^(1-a)/poch(a,n)*diff(z^(a+n-1)*fun,z,n)
3136 (mul (inv (factf a n))
3137 (power arg (sub 1 a))
3138 ($diff (mul (power arg (sub (add a n) 1))
3139 fun)
3140 arg n)))
3142 ;; Transform F(a,b;c;z) to F(a,b;c-n;z), given F(a,b;c;z)
3143 (defun as-15.2.4 (axax bb c n arg fun)
3144 (declare (ignore axax bb))
3145 (assert (>= n 0))
3146 ;; A&S 15.2.4
3147 ;; F(a,b;c-n;z) = 1/poch(c-n,n)/z^(c-n-1)*diff(z^(c-1)*fun,z,n)
3148 (mul (inv (factf (sub c n) n))
3149 (inv (power arg (sub (sub c n) 1)))
3150 ($diff (mul (power arg (sub c 1))
3151 fun)
3152 arg n)))
3154 ;; Transform F(a,b;c;z) to F(a-n,b;c;z), given F(a,b;c;z)
3155 (defun as-15.2.5 (a b c n arg fun)
3156 ;; A&S 15.2.5
3157 ;; F(a-n,b;c;z) = 1/poch(c-a,n)*z^(1+a-c)*(1-z)^(c+n-a-b)
3158 ;; *diff(z^(c-a+n-1)*(1-z)^(a+b-c)*F(a,b;c;z),z,n)
3159 (assert (>= n 0))
3160 (mul (inv (factf (sub c a) n))
3161 (power arg (sub (add a 1) c))
3162 (power (sub 1 arg)
3163 (sub (add c n) (add a b)))
3164 ($diff (mul (power arg (sub (add c n)
3165 (add a 1)))
3166 (power (sub 1 arg)
3167 (sub (add a b) c))
3168 fun)
3169 arg n)))
3171 ;; Transform F(a,b;c;z) to F(a,b;c+n;z), given F(a,b;c;z)
3172 (defun as-15.2.6 (a b c n arg fun)
3173 ;; A&S 15.2.6
3174 ;; F(a,b;c+n;z) = poch(c,n)/poch(c-a,n)/poch(c-b,n)*(1-z)^(c+n-a-b)
3175 ;; *diff((1-z)^(a+b-c)*fun,z,n)
3176 (assert (>= n 0))
3177 (mul (factf c n)
3178 (inv (factf (sub c a) n))
3179 (inv (factf (sub c b) n))
3180 (inv (power (sub 1 arg) (sub (add a b)
3181 (add c n))))
3182 ($diff (mul (power (sub 1 arg) (sub (add a b) c))
3183 fun)
3184 arg n)))
3186 ;; Transform F(a,b;c;z) to F(a+n, b; c+n; z)
3187 (defun as-15.2.7 (a b c n arg fun)
3188 ;; A&S 15.2.7
3189 ;; F(a+n,b;c+n;z) = (-1)^n*poch(c,n)/poch(a,n)/poch(c-b,n)*(1-z)^(1-a)
3190 ;; *diff((1-z)^(a+n-1)*fun, z, n)
3191 (assert (>= n 0))
3192 (mul (power -1 n)
3193 (factf c n)
3194 (inv (factf a n))
3195 (inv (factf (sub c b) n))
3196 (power (sub 1 arg) (sub 1 a))
3197 ($diff (mul (power (sub 1 arg) (sub (add a n) 1))
3198 fun)
3199 arg n)))
3201 ;; Transform F(a,b;c;z) to F(a-n, b; c-n; z)
3202 (defun as-15.2.8 (axax b c n arg fun)
3203 ;; A&S 15.2.8
3204 ;; F(a-n,b;c-n;z) = 1/poch(c-n,n)/(z^(c-n-1)*(1-z)^(b-c))
3205 ;; *diff(z^(c-1)*(1-z^(b-c+n)*fun, z, n))
3206 (declare (ignore axax))
3207 (assert (>= n 0))
3208 (mul (inv (factf (sub c n) n))
3209 (inv (mul (power arg (sub (sub c n) 1))
3210 (power (sub 1 arg) (sub b c))))
3211 ($diff (mul (power arg (sub c 1))
3212 (power (sub 1 arg) (add (sub b c) n))
3213 fun)
3214 arg n)))
3216 ;; Transform F(a,b;c;z) to F(a+n,b+n;c+n;z)
3217 (defun as-15.2.2 (a b c n arg fun)
3218 ;; A&S 15.2.2
3219 ;; F(a+n,b+n; c+n;z) = poch(c,n)/poch(a,n)/poch(b,n)
3220 ;; *diff(fun, z, n)
3221 (assert (>= n 0))
3222 (mul (factf c n)
3223 (inv (factf a n))
3224 (inv (factf b n))
3225 ($diff fun arg n)))
3227 ;; Transform F(a,b;c;z) to F(a-n,b-n;c-n;z)
3228 (defun as-15.2.9 (a b c n arg fun)
3229 ;; A&S 15.2.9
3230 ;; F(a-n,b-n; c-n;z) = 1/poch(c-n,n)/(z^(c-n-1)*(1-z)^(a+b-c-n))
3231 ;; *diff(z^(c-1)*(1-z)^(a+b-c)*fun, z, n)
3232 (assert (>= n 0))
3233 (mul (inv (factf (sub c n) n))
3234 (inv (mul (power arg (sub (sub c n) 1))
3235 (power (sub 1 arg) (sub (add a b)
3236 (add c n)))))
3237 ($diff (mul (power arg (sub c 1))
3238 (power (sub 1 arg) (sub (add a b) c))
3239 fun)
3240 arg n)))
3242 (defun step4-int (a b c arg)
3243 (if (> a b)
3244 (step4-int b a c arg)
3245 (let* ((s (gensym (symbol-name '#:step4-arg-)))
3246 (m (1- a))
3247 (n (1- b))
3248 (ell (sub c 3//2))
3249 (res (cond ((eq (checksigntm arg) '$negative)
3250 ;; F(1,1;3/2;z) =
3251 ;; -%i*log(%i*sqrt(zn)+sqrt(1-zn))/(sqrt(1-zn)*sqrt(zn))
3252 ;; for z < 0
3253 (let ((root1-z (power (sub 1 s) (inv 2)))
3254 (rootz (power s (inv 2))))
3255 (mul -1 '$%i
3256 (mlog (add (mul '$%i rootz)
3257 root1-z))
3258 (inv root1-z)
3259 (inv rootz))))
3261 ;; F(1,1;3/2;z) = asin(sqrt(zp))/(sqrt(1-zp)*sqrt(zp))
3262 ;; for z > 0
3263 (let ((root1-z (power (sub 1 s) (inv 2)))
3264 (rootz (power s (inv 2))))
3265 (mul (masin rootz)
3266 (inv root1-z)
3267 (inv rootz)))))))
3268 ;; Start with res = F(1,1;3/2;z). Compute F(m,1;3/2;z)
3269 (setf res (as-15.2.3 1 1 3//2 m s res))
3270 ;; We now have res = C*F(m,1;3/2;z). Compute F(m,n;3/2;z)
3271 (setf res (as-15.2.3 1 a 3//2 n s res))
3272 ;; We now have res = C*F(m,n;3/2;z). Now compute F(m,n;3/2+ell;z):
3273 (subst arg s
3274 (cond ((minusp ell)
3275 (as-15.2.4 a b 3//2 (- ell) s res))
3277 (as-15.2.6 a b 3//2 ell s res)))))))
3279 ;;Pattern match for s(ymbolic) + c(onstant) in parameter
3280 (defun s+c (exp)
3281 (m2 exp '((mplus) ((coeffpt)(f nonnump)) ((coeffpp)(c $numberp)))))
3283 (defun nonnump (z)
3284 (cond ((not ($numberp z)) t)
3285 (t nil)))
3287 ;;Algor. III from thesis:determines which Differ. Formula to use
3288 (defun algiii (fun m n aprime)
3289 (let ((mm (abs m))
3290 (nn (abs n)))
3291 (cond ((and (nni m) (nni n))
3292 (cond ((< m n)
3293 (f81 fun m n aprime))
3295 (f85 fun mm nn aprime))))
3296 ((and (hyp-negp n) (hyp-negp m))
3297 (cond ((> (abs m) (abs n))
3298 (f86 fun mm nn aprime))
3300 (f82 fun mm nn aprime))))
3301 ((and (hyp-negp m) (nni n))
3302 (f83 fun mm nn aprime))
3304 (f84 fun mm nn aprime)))))
3306 ;; Factorial function:x*(x+1)*(x+2)...(x+n-1)
3308 ;; FIXME: This appears to be identical to fctrl above
3309 (defun factf (x n)
3310 (cond ((zerop n) 1)
3311 (t (mul x (factf (add x 1) (sub n 1))))))
3313 ;;Formula #85 from Yannis thesis:finds by differentiating F[2,1](a,b,c,z)
3314 ;; given F[2,1](a+m,b,c+n,z) where b=-a and c=1/2, n,m integers
3316 ;; Like F81, except m > n.
3318 ;; F(a+m,-a;c+n;z), m > n, c = 1/2, m and n are non-negative integers
3320 ;; A&S 15.2.3
3321 ;; diff(z^(a+m-n-1)*F(a,-a;1/2;z),z,m-n) = poch(a,m-n)*z^(a-1)*F(a+m-n,-a;1/2;z)
3323 ;; A&S 15.2.7
3324 ;; diff((1-z)^(a+m-1)*F(a+m-n,-a;1/2;z),z,n)
3325 ;; = (-1)^n*poch(a+m-n,n)*poch(1/2+a,n)/poch(1/2,n)*(1-z)^(a+m-n-1)
3326 ;; * F(a+m,-a;1/2+n;z)
3328 (defun f85 (fun m n a)
3329 (mul (factf (inv 2) n)
3330 (inv (power -1 n))
3331 (inv (factf (sub (add a m)
3334 (inv (factf (sub (inv 2)
3335 (mul a -1))
3337 (inv (factf a (- m n)))
3338 (power (sub 1 'ell) (sub (sub (add 1 n) m) a))
3339 ($diff (mul (power (sub 1 'ell) (sub (add a m) 1))
3340 (power 'ell (sub 1 a))
3341 ($diff (mul (power 'ell (sub (add a m -1) n))
3342 fun)
3343 'ell (- m n)))
3344 'ell n)))
3346 ;;Used to find negative things that are not integers,eg RAT's
3347 (defun hyp-negp (x)
3348 (cond ((equal (asksign x) '$negative)
3350 (t nil)))
3352 ;; F(a,-a+m; c+n; z) where m,n are non-negative integers, m < n, c = 1/2.
3354 ;; A&S 15.2.6
3355 ;; diff((1-z)^(a+b-c)*F(a,b;c;z),z,n)
3356 ;; = poch(c-a,n)*poch(c-b,n)/poch(c,n)*(1-z)^(a+b-c-n)*F(a,b;c+n;z)
3358 ;; A&S 15.2.7:
3359 ;; diff((1-z)^(a+m-1))*F(a,b;c;z),z,m)
3360 ;; = (-1)^m*poch(a,m)*poch(c-b,m)/poch(c,m)*(1-z)^(a-1)*F(a+m,b;c+m;z)
3362 ;; Rewrite F(a,-a+m; c+n;z) as F(-a+m, a; c+n; z). Then apply 15.2.6
3363 ;; to F(-a,a;1/2;z), differentiating n-m times:
3365 ;; diff((1-z)^(-1/2)*F(-a,a;1/2;z),z,n-m)
3366 ;; = poch(1/2+a,n-m)*poch(1/2-a,n-m)/poch(1/2,n-m)*(1-z)^(-1/2-n+m)*F(-a,a;1/2+n-m;z)
3368 ;; Multiply this result by (1-z)^(n-a-1/2) and apply 15.2.7, differentiating m times:
3370 ;; diff((1-z)^(m-a-1)*F(-a,a;1/2+n-m;z),z,m)
3371 ;; = (-1)^m*poch(-a,m)*poch(1/2+n-m-a,m)/poch(1/2+n-m)*(1-z)^(-a-1)*F(-a+m,a;1/2+n;z)
3373 ;; Which gives F(-a+m,a;1/2+n;z), which is what we wanted.
3374 (defun f81 (fun m n a)
3375 (mul (factf (add (inv 2) (- n m)) m)
3376 (factf (inv 2) (- n m))
3377 (inv (power -1 m))
3378 (inv (factf a m))
3379 (inv (factf (add (inv 2) n (sub a m)) m))
3380 (inv (factf (sub (inv 2) a) (- n m)))
3381 (inv (factf (add (inv 2) a) (- n m)))
3382 (power (sub 1 'ell) (sub 1 a))
3383 ($diff (mul (power (sub 1 'ell) (add a n (inv -2)))
3384 ($diff (mul (power (sub 1 'ell) (inv -2))
3385 fun)
3386 'ell (- n m)))
3387 'ell m)))
3389 ;; Like f86, but |n|>=|m|
3391 ;; F(a-m,-a;1/2-n;z) where n >= m >0
3393 ;; A&S 15.2.4
3394 ;; diff(z^(c-1)*F(a,b;c;z),z,n)
3395 ;; = poch(c-n,n)*z^(c-n-1)*F(a;b;c-n;z)
3397 ;; A&S 15.2.8:
3398 ;; diff(z^(c-1)*(1-z)^(b-c+n)*F(a,b;c;z),z,n)
3399 ;; = poch(c-n,n)*z^(c-n-1)*(1-z)^(b-c)*F(a-n,b;c-n;z)
3401 ;; For our problem:
3403 ;; diff(z^(-1/2)*F(a,-a;1/2;z),z,n-m)
3404 ;; = poch(1/2-n+m,n-m)*z^(m-n-1/2)*F(a,-a;1/2-n+m;z)
3406 ;; diff(z^(m-n-1/2)*(1-z)^(n-a-1/2)*F(a,-a;1/2-n+m;z),z,m)
3407 ;; = poch(1/2-n,m)*z^(-1/2-n)*(1-z)^(n-m-a-1/2)*F(a-m,-a;1/2-n;z)
3409 ;; So
3411 ;; G(z) = z^(m-n-1/2)*F(a,-a;1/2-n+m;z)
3412 ;; = z^(n-m+1/2)/poch(1/2-n+m,n-m)*diff(z^(-1/2)*F(a,-a;1/2;z),z,n-m)
3414 ;; F(a-m,-a;1/2-n;z)
3415 ;; = z^(n+1/2)*(1-z)^(m+a-1/2-n)/poch(1/2-n,m)*diff((1-z)^(n-a-1/2)*G(z),z,m)
3416 (defun f82 (fun m n a)
3417 (mul (inv (factf (sub (inv 2) n) m))
3418 (inv (factf (sub (add (inv 2) m) n) (- n m)))
3419 (power 'ell (add n (inv 2)))
3420 (power (sub 1 'ell) (sub (add m (inv 2) a) n))
3421 ($diff (mul (power (sub 1 'ell)
3422 (sub (sub n a) (inv 2)))
3423 ($diff (mul (power 'ell (inv -2)) fun)
3424 'ell
3425 (- n m)))
3426 'ell
3427 m)))
3429 ;; F(a+m,-a;1/2+n;z) with m,n integers and m < 0, n >= 0
3431 ;; Write this more clearly as F(a-m,-a;1/2+n;z), m > 0, n >= 0
3432 ;; or equivalently F(a-m,-a;c+n;z)
3434 ;; A&S 15.2.6
3435 ;; diff((1-z)^(-1/2)*F(a,-a;1/2;z),z,n)
3436 ;; = poch((1/2+a,n)*poch(1/2-a,n)/poch(1/2,n)*(1-z)^(-1/2-n)
3437 ;; * F(a,-a;1/2+n;z)
3439 ;; A&S 15.2.5
3440 ;; diff(z^(n+m-a-1/2)*(1-z)^(-1/2-n)*F(a,-a;1/2+n;z),z,m)
3441 ;; = poch(1/2+n-a,m)*z^(1/2+n-a-1)*(1-z)^(-1/2-n-m)*F(a-m,-a;1/2+n;z)
3442 ;; = poch(1/2+n-a,m)*z^(n-a-1/2)*(1-z)^(-1/2-n-m)*F(a-m,-a;1/2+n;z)
3444 ;; (1-z)^(-1/2-n)*F(a,-a;1/2+n;z)
3445 ;; = poch(1/2,n)/poch(1/2-a,n)/poch(1/2+a,n)*diff((1-z)^(-1/2)*F(a,-a;1/2;z),z,n)
3447 ;; F(a-m,-a;1/2+n;z)
3448 ;; = (1-z)^(n+m+1/2)*z^(a-n+1/2)/poch(1/2+n-a,m)
3449 ;; *diff(z^(n+m-a-1/2)*(1-z)^(-1/2-n)*F(a,-a;1/2+n;z),z,m)
3450 (defun f83 (fun m n a)
3451 (mul (factf (inv 2) n)
3452 (inv (factf (sub (inv 2) a) n))
3453 (inv (factf (sub (add n (inv 2)) a) m))
3454 (inv (factf (add (inv 2) a) n))
3455 (power (sub 1 'ell) (add m n (inv 2)))
3456 (power 'ell (add (sub a n) (inv 2)))
3457 ($diff (mul (power 'ell (sub (sub (+ m n) a) (inv 2)))
3458 ($diff (mul (power (sub 1 'ell)
3459 (inv -2))
3460 fun)
3461 'ell
3463 'ell
3464 m)))
3466 ;; The last case F(a+m,-a;c+n;z), m,n integers, m >= 0, n < 0
3468 ;; F(a+m,-a;1/2-n;z)
3470 ;; A&S 15.2.4:
3471 ;; diff(z^(c-1)*F(a,b;c;z),z,n) = poch(c-n,n)*z^(c-n-1)*F(a,b;c-n;z)
3473 ;; A&S 15.2.3:
3474 ;; diff(z^(a+m-1)*F(a,b;c;z),z,m) = poch(a,n)*z^(a-1)*F(a+n,b;c;z)
3476 ;; For our problem:
3478 ;; diff(z^(-1/2)*F(a,-a;1/2;z),z,n) = poch(1/2-n,n)*z^(-n-1/2)*F(a,-a;1/2-n;z)
3480 ;; diff(z^(a+m-1)*F(a,-a;1/2-n;z),z,m) = poch(a,m)*z^(a-1)*F(a+m,-a;1/2-n;z)
3481 (defun f84 (fun m n a)
3482 (mul (inv (mul (factf a m)
3483 (factf (sub (inv 2) n) n)))
3484 (power 'ell (sub 1 a))
3485 ($diff (mul (power 'ell (sub (add a m n) (inv 2)))
3486 ($diff (mul (power 'ell (inv -2)) fun)
3487 'ell
3489 'ell
3490 m)))
3492 ;; Like f82, but |n|<|m|
3494 ;; F(a-m,-a;1/2-n;z), 0 < n < m
3496 ;; A&S 15.2.5
3497 ;; diff(z^(c-a+n-1)*(1-z)^(a+b-c)*F(a,b;c;z),z,n)
3498 ;; = poch(c-a,n)*z^(c-a-1)*(1-z)^(a+b-c-n)*F(a-n,b;c;z)
3500 ;; A&S 15.2.8:
3501 ;; diff(z^(c-1)*(1-z)^(b-c+n)*F(a,b;c;z),z,n)
3502 ;; = poch(c-n,n)*z^(c-n-1)*(1-z)^(b-c)*F(a-n,b;c-n;z)
3504 ;; For our problem:
3506 ;; diff(z^(-a+m-n-1/2)*(1-z)^(-1/2)*F(a,-a;1/2;z),z,m-n)
3507 ;; = poch(1/2-a,m-n)*z^(-a-1/2)*(1-z)^(-1/2-m+n)*F(a-m+n,-a;1/2;z)
3509 ;; diff(z^(-1/2)*(1-z)^(-a-1/2+n)*F(a-m+n,-a;1/2;z),z,n)
3510 ;; = poch(1/2-n,n)*z^(-n-1/2)*(1-z)^(-a-1/2)*F(a-m,-a;1/2-n;z)
3512 ;; G(z) = z^(-a-1/2)*(1-z)^(-1/2-m+n)*F(a-m+n,-a;1/2;z)
3513 ;; = 1/poch(1/2-a,m-n)*diff(z^(-a+m-n-1/2)*(1-z)^(-1/2)*F(a,-a;1/2;z),z,m-n)
3515 ;; F(a-m,-a;1/2-n;z)
3516 ;; = z^(n+1/2)*(1-z)^(a+1/2)/poch(1/2-n,n)
3517 ;; *diff(z^(-1/2)*(1-z)^(-a-1/2+n)*F(a-m+n,-a;1/2;z),z,n)
3518 ;; = z^(n+1/2)*(1-z)^(a+1/2)/poch(1/2-n,n)
3519 ;; *diff(z^a*(1-z)^(m-a)*G(z),z,n)
3521 (defun f86 (fun m n a)
3522 (mul (inv (mul (factf (sub (inv 2) n) n)
3523 (factf (sub (inv 2) a) (- m n))))
3524 (power 'ell (add n (inv 2)))
3525 (power (sub 1 'ell) (add (inv 2) a))
3526 ($diff (mul (power 'ell a)
3527 (power (sub 1 'ell) (sub m a))
3528 ($diff (mul (power 'ell
3529 (sub (sub (sub m n) (inv 2)) a))
3530 (power (sub 1 'ell)
3531 (inv -2))
3532 fun)
3533 'ell (- m n)))
3534 'ell n)))
3536 ;; F(-1/2+n, 1+m; 1/2+l; z)
3537 (defun hyp-atanh (a b c arg)
3538 ;; We start with F(-1/2,1;1/2;z) = 1-sqrt(z)*atanh(sqrt(z)). We can
3539 ;; derive the remaining forms by differentiating this enough times.
3541 ;; FIXME: Do we need to assume z > 0? We do that anyway, here.
3542 (let* ((s (gensym (symbol-name '#:hyp-atanh-)))
3543 (n (add a 1//2))
3544 (m (sub b 1))
3545 (ell (sub c 1//2))
3546 (res (sub 1 (mul (power s '((rat simp) 1 2))
3547 (take '(%atanh) (power s '((rat simp) 1 2))))))
3548 (new-a -1//2)
3549 (new-b 1)
3550 (new-c 1//2))
3551 ;; The total number of derivates we compute is n + m + ell. We
3552 ;; should do something to reduce the number of derivatives.
3553 #+nil
3554 (progn
3555 (format t "a ,b ,c = ~a ~a ~a~%" a b c)
3556 (format t "n, m, ell = ~a ~a ~a~%" n m ell)
3557 (format t "init a b c = ~a ~a ~a~%" new-a new-b new-c))
3558 (cond ((alike1 (sub n ell) 0)
3559 ;; n = ell so we can use A&S 15.2.7 or A&S 15.2.8
3560 (cond ((plusp n)
3561 (setf res (as-15.2.7 new-a new-b new-c n s res)))
3563 (setf res (as-15.2.8 new-a new-b new-c (- n) s res))))
3564 (setf new-a (add new-a n))
3565 (setf new-c (add new-c n)))
3567 ;; Adjust ell and then n. (Does order matter?)
3568 (cond ((plusp ell)
3569 (setf res (as-15.2.6 new-a new-b new-c ell s res))
3570 (setf new-c (add new-c ell)))
3572 (setf res (as-15.2.4 new-a new-b new-c (- ell) s res))
3573 (setf new-c (add new-c ell))))
3574 #+nil
3575 (progn
3576 (maxima-display res)
3577 (format t "new a b c = ~a ~a ~a~%" new-a new-b new-c))
3578 (cond ((plusp n)
3579 ;; A&S 15.2.3
3580 (setf res (as-15.2.3 new-a new-b new-c n s res))
3581 (setf new-a (add new-a n)))
3583 ;; A&S 15.2.5
3584 (setf res (as-15.2.5 new-a new-b new-c (- n) s res))
3585 (setf new-a (add new-a n))))))
3586 #+nil
3587 (progn
3588 (format t "new a b c = ~a ~a ~a~%" new-a new-b new-c)
3589 (maxima-display res))
3590 ;; Finally adjust m by swapping the a and b parameters, since the
3591 ;; hypergeometric function is symmetric in a and b.
3592 (cond ((plusp m)
3593 (setf res (as-15.2.3 new-b new-a new-c m s res))
3594 (setf new-b (add new-b m)))
3596 (setf res (as-15.2.5 new-b new-a new-c (- m) s res))
3597 (setf new-b (add new-b m))))
3598 #+nil
3599 (progn
3600 (format t "new a b c = ~a ~a ~a~%" new-a new-b new-c)
3601 (maxima-display res))
3602 ;; Substitute the argument into the expression and simplify the result.
3603 (sratsimp (maxima-substitute arg s res))))