1 ;;; -*- Mode: Lisp; Package: Maxima; Syntax: Common-Lisp; Base: 10 -*- ;;;;
3 ;;; ** (c) Copyright 1976, 1983 Massachusetts Institute of Technology **
6 ;;; These are the main routines for finding the Laplace Transform
7 ;;; of special functions --- written by Yannis Avgoustis
8 ;;; --- modified by Edward Lafferty
9 ;;; Latest mod by jpg 8/21/81
11 ;;; This program uses the programs on ELL;HYP FASL.
15 ;;; Definite integration using the generalized hypergeometric functions
16 ;;; Avgoustis, Ioannis Dimitrios
17 ;;; Thesis. 1977. M.S.--Massachusetts Institute of Technology. Dept.
18 ;;; of Electrical Engineering and Computer Science
19 ;;; http://dspace.mit.edu/handle/1721.1/16269
21 ;;; Avgoustis, I. D., Symbolic Laplace Transforms of Special Functions,
22 ;;; Proceedings of the 1977 MACSYMA Users' Conference, pp 21-41
26 (macsyma-module hypgeo
)
28 (declare-top (special checkcoefsignlist $exponentialize $radexpand $logexpand
31 (defmvar $prefer_d nil
32 "When NIL express a parabolic cylinder function as hypergeometric function.")
34 ;; The properties NOUN and VERB give correct linear display
35 (defprop $specint %specint verb
)
36 (defprop %specint $specint noun
)
38 (defvar *hyp-return-noun-form-p
* t
39 "Return noun form instead of internal Lisp symbols for integrals
40 that we don't support.")
42 (defvar *hyp-return-noun-flag
* nil
43 "Flag to signal that we have no result and to return a noun.")
45 (defvar *debug-hypgeo
* nil
46 "Print debug information if enabled.")
48 ;; The variables *var* and *par* are global to this file only.
49 ;; They are initialized in the routine defexec. The values are never changed.
50 ;; These globals are introduced to avoid passing the values of *par* and *var*
51 ;; through all functions of this file.
54 "Variable of integration of Laplace transform.")
56 "Parameter of Laplace transform.")
58 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
60 ;;; Helper function for this file
62 (defun substl (p1 p2 p3
)
64 (t (maxima-substitute p1 p2 p3
))))
66 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
68 ;;; Test functions for pattern match, which use the globals var and *par*
74 (cond ((equal m
0) nil
)
77 ;;; Test functions which do not depend on globals
79 ;; Test if EXP is 1 or %e.
87 (defun free-not-zero-p (expr x
)
88 (and (not (equal expr
0)) (free expr x
)))
90 (defun free2 (expr x y
)
91 (and (free expr x
) (free expr y
)))
93 (defun has-not-alike1-p (expr x
)
94 (and (not (alike1 expr x
)) (has expr x
)))
96 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
98 ;;; Some shortcuts for special functions.
100 ;; Lommel's little s[u,v](z) function.
101 (defun littleslommel (m n z
)
102 (list '(mqapply simp
) (list '($%s simp array
) m n
) z
))
104 ;; Whittaker's M function
105 (defun mwhit (a i1 i2
)
106 (list '(mqapply simp
) (list '($%m simp array
) i1 i2
) a
))
108 ;; Whittaker's W function
109 (defun wwhit (a i1 i2
)
110 (list '(mqapply simp
) (list '($%w simp array
) i1 i2
) a
))
113 (defun pjac (x n a b
)
114 (list '(mqapply simp
) (list '($%p simp array
) n a b
) x
))
116 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
118 ;;; Two general pattern for the routine lt-sf-log.
120 ;; Recognize c*u^v + a and a=0.
121 (defun m2-arbpow1 (expr var
)
125 (c free
,var
) ; more special to ensure that c is constant
126 ((mexpt) (u has
,var
) (v free
,var
)))
127 ((coeffpp) (a zerp
)))))
129 ;; Recognize c*u^v*(a+b*u)^w+d and d=0. This is a generalization of arbpow1.
130 (defun m2-arbpow2 (expr var
)
134 ((coefftt) (c free
,var
))
135 ((mexpt) (u equal
,var
) (v free
,var
))
137 ((mplus) (a free
,var
) ((coefft) (b free
,var
) (u equal
,var
)))
138 (w free-not-zero-p
,var
)))
139 ((coeffpp) (d zerp
)))))
141 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
143 ;;; The pattern to match special functions in the routine lt-sf-log.
146 (defun m2-asin (expr var
)
149 ((coeffpt) (u nonzerp
)
150 ((%asin
) (w has
,var
)))
151 ((coeffpp) (a equal
0)))))
154 (defun m2-atan (expr var
)
157 ((coeffpt) (u nonzerp
)
158 ((%atan
) (w has
,var
)))
159 ((coeffpp) (a equal
0)))))
161 ;; Recognize bessel_j(v,w)
162 (defun m2-onej (expr var
)
167 ((%bessel_j
) (v free
,var
) (w has
,var
)))
168 ((coeffpp) (a zerp
)))))
170 ;; Recognize bessel_j(v1,w1)*bessel_j(v2,w2)
171 (defun m2-twoj (expr var
)
176 ((%bessel_j
) (v1 free
,var
) (w1 has
,var
))
177 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
178 ((coeffpp) (a zerp
)))))
180 ;; Recognize bessel_y(v1,w1)*bessel_y(v2,w2)
181 (defun m2-twoy (expr var
)
186 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
187 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
188 ((coeffpp) (a zerp
)))))
190 ;; Recognize bessel_k(v1,w1)*bessel_k(v2,w2)
191 (defun m2-twok (expr var
)
196 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
197 ((%bessel_k
) (v2 free
,var
) (w2 has
,var
)))
198 ((coeffpp) (a zerp
)))))
200 ;; Recognize bessel_k(v1,w1)*bessel_y(v2,w2)
201 (defun m2-onekoney (expr var
)
206 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
207 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
208 ((coeffpp) (a zerp
)))))
210 ;; Recognize bessel_j(v,w)^2
211 (defun m2-onej^
2 (expr var
)
217 ((%bessel_j
) (v free
,var
) (w has
,var
))
219 ((coeffpp) (a zerp
)))))
221 ;; Recognize bessel_y(v,w)^2
222 (defun m2-oney^
2 (expr var
)
228 ((%bessel_y
) (v free
,var
) (w has
,var
))
230 ((coeffpp) (a zerp
)))))
232 ;; Recognize bessel_k(v,w)^2
233 (defun m2-onek^
2 (expr var
)
239 ((%bessel_k
) (v free
,var
) (w has
,var
))
241 ((coeffpp) (a zerp
)))))
243 ;; Recognize bessel_i(v,w)
244 (defun m2-onei (expr var
)
249 ((%bessel_i
) (v free
,var
) (w has
,var
)))
250 ((coeffpp) (a zerp
)))))
252 ;; Recognize bessel_i(v1,w1)*bessel_i(v2,w2)
253 (defun m2-twoi (expr var
)
258 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
259 ((%bessel_i
) (v2 free
,var
) (w2 has
,var
)))
260 ((coeffpp) (a zerp
)))))
262 ;; Recognize hankel_1(v1,w1)*hankel_1(v2,w2), product of 2 Hankel 1 functions.
263 (defun m2-two-hankel_1 (expr var
)
268 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
269 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
270 ((coeffpp) (a zerp
)))))
272 ;; Recognize hankel_2(v1,w1)*hankel_2(v2,w2), product of 2 Hankel 2 functions.
273 (defun m2-two-hankel_2 (expr var
)
278 ((%hankel_2
) (v1 free
,var
) (w1 has
,var
))
279 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
280 ((coeffpp) (a zerp
)))))
282 ;; Recognize hankel_1(v1,w1)*hankel_2(v2,w2), product of 2 Hankel functions.
283 (defun m2-hankel_1*hankel_2
(expr var
)
288 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
289 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
290 ((coeffpp) (a zerp
)))))
292 ;; Recognize bessel_y(v1,w1)*bessel_j(v2,w2)
293 (defun m2-oneyonej (expr var
)
298 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
299 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
300 ((coeffpp) (a zerp
)))))
302 ;; Recognize bessel_k(v1,w1)*bessel_j(v2,w2)
303 (defun m2-onekonej (expr var
)
308 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
309 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
310 ((coeffpp) (a zerp
)))))
312 ;; Recognize bessel_y(v1,w1)*hankel_1(v2,w2)
313 (defun m2-bessel_y*hankel_1
(expr var
)
318 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
319 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
320 ((coeffpp) (a zerp
)))))
322 ;; Recognize bessel_y(v1,w1)*hankel_2(v2,w2)
323 (defun m2-bessel_y*hankel_2
(expr var
)
328 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
329 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
330 ((coeffpp) (a zerp
)))))
332 ;; Recognize bessel_k(v1,w1)*hankel_1(v2,w2)
333 (defun m2-bessel_k*hankel_1
(expr var
)
338 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
339 ((%hankel_1
) (v1 free
,var
) (w2 has
,var
)))
340 ((coeffpp) (a zerp
)))))
342 ;; Recognize bessel_k(v1,w1)*hankel_2(v2,w2)
343 (defun m2-bessel_k*hankel_2
(expr var
)
348 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
349 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
350 ((coeffpp) (a zerp
)))))
352 ;; Recognize bessel_i(v1,w1)*bessel_j(v2,w2)
353 (defun m2-oneionej (expr var
)
358 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
359 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
360 ((coeffpp) (a zerp
)))))
362 ;; Recognize bessel_i(v1,w1)*hankel_1(v2,w2)
363 (defun m2-bessel_i*hankel_1
(expr var
)
368 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
369 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
370 ((coeffpp) (a zerp
)))))
372 ;; Recognize bessel_i(v1,w1)*hankel_2(v2,w2)
373 (defun m2-bessel_i*hankel_2
(expr var
)
378 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
379 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
380 ((coeffpp) (a zerp
)))))
382 ;; Recognize hankel_1(v1,w1)*bessel_j(v2,w2)
383 (defun m2-hankel_1*bessel_j
(expr var
)
388 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
389 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
390 ((coeffpp) (a zerp
)))))
392 ;; Recognize hankel_2(v1,w1)*bessel_j(v2,w2)
393 (defun m2-hankel_2*bessel_j
(expr var
)
398 ((%hankel_2
) (v1 free
,var
) (w1 has
,var
))
399 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
400 ((coeffpp) (a zerp
)))))
402 ;; Recognize bessel_i(v1,w1)*bessel_y(v2,w2)
403 (defun m2-oneioney (expr var
)
408 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
409 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
410 ((coeffpp) (a zerp
)))))
412 ;; Recognize bessel_i(v1,w1)*bessel_k(v2,w2)
413 (defun m2-oneionek (expr var
)
418 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
419 ((%bessel_k
) (v2 free
,var
) (w2 has
,var
)))
420 ((coeffpp) (a zerp
)))))
422 ;; Recognize bessel_i(v,w)^2
423 (defun m2-onei^
2 (expr var
)
429 ((%bessel_i
) (v free
,var
) (w has
,var
))
431 ((coeffpp) (a zerp
)))))
433 ;; Recognize hankel_1(v,w)^2
434 (defun m2-hankel_1^
2 (expr var
)
440 ((%hankel_1
) (v free
,var
) (w has
,var
))
442 ((coeffpp) (a zerp
)))))
444 ;; Recognize hankel_2(v,w)^2
445 (defun m2-hankel_2^
2 (expr var
)
451 ((%hankel_2
) (v free
,var
) (w has
,var
))
453 ((coeffpp) (a zerp
)))))
455 ;; Recognize bessel_y(v,w)
456 (defun m2-oney (expr var
)
461 ((%bessel_y
) (v free
,var
) (w has
,var
)))
462 ((coeffpp) (a zerp
)))))
464 ;; Recognize bessel_k(v,w)
465 (defun m2-onek (expr var
)
470 ((%bessel_k
) (v free
,var
) (w has
,var
)))
471 ((coeffpp) (a zerp
)))))
473 ;; Recognize hankel_1(v,w)
474 (defun m2-hankel_1 (expr var
)
479 ((%hankel_1
) (v free
,var
) (w has
,var
)))
480 ((coeffpp) (a zerp
)))))
482 ;; Recognize hankel_2(v,w)
483 (defun m2-hankel_2 (expr var
)
488 ((%hankel_2
) (v free
,var
) (w has
,var
)))
489 ((coeffpp) (a zerp
)))))
492 (defun m2-onelog (expr var
)
497 ((%log
) (w has
,var
)))
498 ((coeffpp) (a zerp
)))))
501 (defun m2-onerf (expr var
)
506 ((%erf
) (w has
,var
)))
507 ((coeffpp) (a zerp
)))))
510 (defun m2-onerfc (expr var
)
515 ((%erfc
) (w has
,var
)))
516 ((coeffpp) (a zerp
)))))
518 ;; Recognize fresnel_s(w)
519 (defun m2-onefresnel_s (expr var
)
524 ((%fresnel_s
) (w has
,var
)))
525 ((coeffpp) (a zerp
)))))
527 ;; Recognize fresnel_c(w)
528 (defun m2-onefresnel_c (expr var
)
533 ((%fresnel_c
) (w has
,var
)))
534 ((coeffpp) (a zerp
)))))
536 ;; Recognize expintegral_e(v,w)
537 (defun m2-oneexpintegral_e (expr var
)
542 ((%expintegral_e
) (v free
,var
) (w has
,var
)))
543 ((coeffpp) (a zerp
)))))
545 ;; Recognize expintegral_ei(w)
546 (defun m2-oneexpintegral_ei (expr var
)
551 ((%expintegral_ei
) (w has
,var
)))
552 ((coeffpp) (a zerp
)))))
554 ;; Recognize expintegral_e1(w)
555 (defun m2-oneexpintegral_e1 (expr var
)
560 ((%expintegral_e1
) (w has
,var
)))
561 ((coeffpp) (a zerp
)))))
563 ;; Recognize expintegral_si(w)
564 (defun m2-oneexpintegral_si (expr var
)
569 ((%expintegral_si
) (w has
,var
)))
570 ((coeffpp) (a zerp
)))))
572 ;; Recognize expintegral_shi(w)
573 (defun m2-oneexpintegral_shi (expr var
)
578 ((%expintegral_shi
) (w has
,var
)))
579 ((coeffpp) (a zerp
)))))
581 ;; Recognize expintegral_ci(w)
582 (defun m2-oneexpintegral_ci (expr var
)
587 ((%expintegral_ci
) (w has
,var
)))
588 ((coeffpp) (a zerp
)))))
590 ;; Recognize expintegral_chi(w)
591 (defun m2-oneexpintegral_chi (expr var
)
596 ((%expintegral_chi
) (w has
,var
)))
597 ((coeffpp) (a zerp
)))))
599 ;; Recognize kelliptic(w), (new would be elliptic_kc)
600 (defun m2-onekelliptic (expr var
)
605 (($kelliptic
) (w has
,var
)))
606 ((coeffpp) (a zerp
)))))
608 ;; Recognize elliptic_kc
609 (defun m2-elliptic_kc (expr var
)
614 ((%elliptic_kc
) (w has
,var
)))
615 ((coeffpp) (a zerp
)))))
617 ;; Recognize %e(w), (new would be elliptic_ec)
618 (defun m2-onee (expr var
)
623 (($%e
) (w has
,var
)))
624 ((coeffpp) (a zerp
)))))
626 ;; Recognize elliptic_ec
627 (defun m2-elliptic_ec (expr var
)
632 ((%elliptic_ec
) (w has
,var
)))
633 ((coeffpp) (a zerp
)))))
635 ;; Recognize gamma_incomplete(w1, w2), Incomplete Gamma function
636 (defun m2-onegammaincomplete (expr var
)
641 ((%gamma_incomplete
) (w1 free
,var
) (w2 has
,var
)))
642 ((coeffpp) (a zerp
)))))
644 ;; Recognize gamma_incomplete_lower(w1,w2), gamma(a)-gamma_incomplete(w1,w2)
645 (defun m2-onegamma-incomplete-lower (expr var
)
650 ((%gamma_incomplete_lower
) (w1 free
,var
) (w2 has
,var
)))
651 ((coeffpp) (a zerp
)))))
653 ;; Recognize Struve H function.
654 (defun m2-struve_h (expr var
)
659 ((%struve_h
) (v free
,var
) (w has
,var
)))
660 ((coeffpp)(a zerp
)))))
662 ;; Recognize Struve L function.
663 (defun m2-struve_l (expr var
)
668 ((%struve_l
) (v free
,var
) (w has
,var
)))
669 ((coeffpp) (a zerp
)))))
671 ;; Recognize Lommel s[v1,v2](w) function.
672 (defun m2-ones (expr var
)
678 (($%s array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
679 ((coeffpp) (a zerp
)))))
681 ;; Recognize S[v1,v2](w), Lommel function
682 (defun m2-oneslommel (expr var
)
688 (($slommel array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
689 ((coeffpp) (a zerp
)))))
691 ;; Recognize parabolic_cylinder_d function
692 (defun m2-parabolic_cylinder_d (expr var
)
697 (($parabolic_cylinder_d
) (v free
,var
) (w has
,var
)))
698 ((coeffpp) (a zerp
)))))
700 ;; Recognize kbatmann(v,w), Batemann function
701 (defun m2-onekbateman (expr var
)
707 (($kbateman array
) (v free
,var
)) (w has
,var
)))
708 ((coeffpp) (a zerp
)))))
710 ;; Recognize %l[v1,v2](w), Generalized Laguerre function
711 (defun m2-onel (expr var
)
717 (($%l array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
718 ((coeffpp) (a zerp
)))))
720 ;; Recognize gen_laguerre(v1,v2,w), Generalized Laguerre function
721 (defun m2-one-gen-laguerre (expr var
)
726 ((%gen_laguerre
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
727 ((coeffpp) (a zerp
)))))
729 ;; Recognize laguerre(v1,w), Laguerre function
730 (defun m2-one-laguerre (expr var
)
735 ((%laguerre
) (v1 free
,var
) (w has
,var
)))
736 ((coeffpp) (a zerp
)))))
738 ;; Recognize %c[v1,v2](w), Gegenbauer function
739 (defun m2-onec (expr var
)
745 (($%c array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
746 ((coeffpp) (a zerp
)))))
748 ;; Recognize %t[v1](w), Chebyshev function of the first kind
749 (defun m2-onet (expr var
)
754 ((mqapply) (($%t array
) (v1 free
,var
)) (w has
,var
)))
755 ((coeffpp) (a zerp
)))))
757 ;; Recognize %u[v1](w), Chebyshev function of the second kind
758 (defun m2-oneu (expr var
)
763 ((mqapply) (($%u array
) (v1 free
,var
)) (w has
,var
)))
764 ((coeffpp) (a zerp
)))))
766 ;; Recognize %p[v1,v2,v3](w), Jacobi function
767 (defun m2-onepjac (expr var
)
774 (v1 free
,var
) (v2 free
,var
) (v3 free
,var
)) (w has
,var
)))
775 ((coeffpp) (a zerp
)))))
777 ;; Recognize jacobi_p function
778 (defun m2-jacobi_p (expr var
)
784 (v1 free
,var
) (v2 free
,var
) (v3 free
,var
) (w has
,var
)))
785 ((coeffpp) (a zerp
)))))
787 ;; Recognize %p[v1,v2](w), Associated Legendre P function
788 (defun m2-hyp-onep (expr var
)
794 (($%p array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
795 ((coeffpp) (a zerp
)))))
797 ;; Recognize assoc_legendre_p function
798 (defun m2-assoc_legendre_p (expr var
)
803 (($assoc_legendre_p
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
804 ((coeffpp) (a zerp
)))))
806 ;; Recognize %p[v1](w), Legendre P function
807 (defun m2-onep0 (expr var
)
812 ((mqapply)(($%p array
) (v1 free
,var
)) (w has
,var
)))
813 ((coeffpp) (a zerp
)))))
815 ;; Recognize %p[v1](w), Legendre P function
816 (defun m2-legendre_p (expr var
)
821 (($legendre_p
) (v free
,var
)) (w has
,var
))
822 ((coeffpp) (a zerp
)))))
824 ;; Recognize hermite(v1,w), Hermite function
825 (defun m2-one-hermite (expr var
)
830 ((%hermite
) (v1 free
,var
) (w has
,var
)))
831 ((coeffpp) (a zerp
)))))
833 ;; Recognize %q[v1,v2](w), Associated Legendre function of the second kind
834 (defun m2-oneq (expr var
)
840 (($%q array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
841 ((coeffpp) (a zerp
)))))
843 ;; Recognize assoc_legendre_q function
844 (defun m2-assoc_legendre_q (expr var
)
849 (($assoc_legendre_q
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
850 ((coeffpp) (a zerp
)))))
852 ;; Recognize %w[v1,v2](w), Whittaker W function.
853 (defun m2-onew (expr var
)
859 (($%w array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
860 ((coeffpp) (a zerp
)))))
862 ;; Recognize whittaker_w function.
863 (defun m2-whittaker_w (expr var
)
868 (($whittaker_w
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
869 ((coeffpp) (a zerp
)))))
871 ;; Recognize %m[v1,v2](w), Whittaker M function
872 (defun m2-onem (expr var
)
878 (($%m array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
879 ((coeffpp) (a zerp
)))))
881 ;; Recognize whittaker_m function.
882 (defun m2-whittaker_m (expr var
)
887 (($whittaker_m
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
888 ((coeffpp) (a zerp
)))))
890 ;; Recognize %f[v1,v2](w1,w2,w3), Hypergeometric function
891 (defun m2-onef (expr var
)
897 (($%f array
) (v1 free
,var
) (v2 free
,var
))
901 ((coeffpp) (a zerp
)))))
903 ;; Recognize hypergeometric function
904 (defun m2-hypergeometric (expr var
)
909 (($hypergeometric
) (w1 free
,var
) (w2 free
,var
) (w3 has
,var
)))
910 ((coeffpp) (a zerp
)))))
912 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
914 ;;; Pattern for the routine hypgeo-exec.
915 ;;; RECOGNIZES L.T.E. "U*%E^(A*X+E*F(X)-P*X+C)+D".
917 (defun m2-ltep (expr var par
)
925 ((coeffpt) (a free2
,var
,par
) (x alike1
,var
))
926 ((coeffpt) (e free2
,var
,par
) (f has
,var
))
927 ((mtimes) -
1 (p alike1
,par
) (x alike1
,var
))
928 ((coeffpp) (c free2
,var
,par
)))))
929 ((coeffpp) (d equal
0)))))
931 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
933 ;;; Pattern for the routine defexec.
934 ;;; This is trying to match EXP to u*%e^(a*x+e*f+c)+d
935 ;;; where a, c, and e are free of x, f is free of p, and d is 0.
937 (defun m2-defltep (expr var
)
945 ((coeffpt) (a free
,var
) (x alike1
,var
))
946 ((coeffpt) (e free
,var
) (f has-not-alike1-p
,var
))
947 ((coeffpp) (c free
,var
)))))
948 ((coeffpp) (d equal
0)))))
950 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
952 ;;; $specint is the Maxima User function
954 (defmfun $specint
(expr var
)
955 (prog ($radexpand checkcoefsignlist
)
956 (setq $radexpand
'$all
)
957 (return (defintegrate expr var
))))
959 (defun defintegrate (expr var
)
960 ;; This used to have $exponentialize enabled for everything, but I
961 ;; don't think we should do that. If various routines want
962 ;; $exponentialize, let them set it themselves. So, for here, we
963 ;; want to expand the form with $exponentialize to convert trig and
964 ;; hyperbolic functions to exponential functions that we can handle.
965 (let ((form (let (($exponentialize t
))
966 ($factor
(resimplify expr
))))) ; At first factor the integrand.
968 ;; Because we call defintegrate recursively, we add code to end the
972 (cond ((and (numberp form
) (zerop form
)) (return-from defintegrate
0))
973 (t (return-from defintegrate
(list '(%specint simp
) form var
)))))
975 ;; We try to find a constant denominator. This is necessary to get results
976 ;; for integrands like u(t)/(a+b+c+...).
978 (let ((den ($denom form
)))
979 (when (and (not (equal 1 den
)) ($freeof var den
))
980 (return-from defintegrate
981 (div (defintegrate (mul den form
) var
) den
))))
983 ;; We search for a sum of Exponential functions which we can integrate.
984 ;; This code finds result for Trigonometric or Hyperbolic functions with
985 ;; a factor t^-1 or t^-2 e.g. t^-1*sin(a*t).
987 (let* ((l (m2-defltep form var
))
988 (s (mul -
1 (cdras 'a l
)))
989 (u ($expand
(cdras 'u l
)))
992 ((setq l1
(m2-sum-with-exp-case1 u var
))
993 ;; c * t^-1 * (%e^(-a*t) - %e^(-b*t)) + d
994 (let ((c (cdras 'c l1
))
995 (a (mul -
1 (cdras 'a l1
)))
996 (b (mul -
1 (cdras 'b l1
)))
998 (add (mul c
(take '(%log
) (div (add s b
) (add s a
))))
999 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1001 ((setq l1
(m2-sum-with-exp-case2 u var
))
1002 ;; c * t^(-3/2) * (%e^(-a*t) - %e^(-b*t)) + d
1003 (let ((c (cdras 'c l1
))
1004 (a (mul -
1 (cdras 'a l1
)))
1005 (b (mul -
1 (cdras 'b l1
)))
1008 (power '$%pi
'((rat simp
) 1 2))
1009 (sub (power (add s b
) '((rat simp
) 1 2))
1010 (power (add s a
) '((rat simp
) 1 2))))
1011 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1013 ((setq l1
(m2-sum-with-exp-case3 u var
))
1014 ;; c * t^-2 * (1 - 2 * %e^(-a*t) + %e^(2*a*t)) + d
1015 (let ((c (cdras 'c l1
))
1016 (a (div (cdras 'a l1
) -
2))
1019 (add (mul (add s a a
) (take '(%log
) (add s a a
)))
1020 (mul s
(take '(%log
) s
))
1021 (mul -
2 (add s a
) (take '(%log
) (add s a
)))))
1022 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1024 ((setq l1
(m2-sum-with-exp-case4 u var
))
1025 ;; c * t^-1 * (1 - 2 * %e^(-a*t) + %e^(2*a*t)) + d
1026 (let ((c (cdras 'c l1
))
1027 (a (div (cdras 'a l1
) (mul 4 '$%i
)))
1033 (mul (sub s
(mul 2 '$%i a
))
1034 (sub s
(mul 2 '$%i a
)))))))
1035 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1037 ((setq l1
(m2-sum-with-exp-case5 u var
))
1038 ;; c * t^-1 * (1 - %e^(2*a*t)) + d
1039 (let ((c (cdras 'c l1
))
1042 (add (mul c
(take '(%log
) (div (sub s a
) s
)))
1043 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1046 ;; At this point we expand the integrand.
1047 (distrdefexecinit ($expand form
) var
))))))
1049 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1051 ;;; Five pattern to find sums of Exponential functions which we can integrate.
1053 ;; Case 1: c * u^-1 * (%e^(-a*u) - %e^(-b*u))
1054 (defun m2-sum-with-exp-case1 (expr var
)
1059 ((mexpt) (u alike1
,var
) -
1)
1062 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1063 ((coeffpp) (z1 zerp
)))))
1065 (c2 equal-times-minus-one c
)
1066 ((mexpt) (u alike1
,var
) -
1)
1069 ((coeffpt) (b nonzerp
) (u alike1
,var
))
1070 ((coeffpp) (z2 zerp
)))))
1071 ((coeffpp) (d true
)))))
1073 ;; Case 2: c * u^(-3/2) * (%e^(-a*u) - %e^(-b*u))
1074 (defun m2-sum-with-exp-case2 (expr var
)
1079 ((mexpt) (u alike1
,var
) ((rat) -
3 2))
1082 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1083 ((coeffpp) (z1 zerp
)))))
1085 (c2 equal-times-minus-one c
)
1086 ((mexpt) (u alike1
,var
) ((rat) -
3 2))
1089 ((coeffpt) (b nonzerp
) (u alike1
,var
))
1090 ((coeffpp) (z2 zerp
)))))
1091 ((coeffpp) (d true
)))))
1093 ;; Case 3: c * u^-2 * (1 - 2 * %e^(-a*u) + %e^(2*a*u))
1094 (defun m2-sum-with-exp-case3 (expr var
)
1099 ((mexpt) (u alike1
,var
) -
2))
1102 ((mexpt) (u alike1
,var
) -
2)
1105 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1106 ((coeffpp) (z1 zerp
)))))
1108 (c3 equal-times-minus-two c
)
1109 ((mexpt) (u alike1
,var
) -
2)
1112 ((coeffpt) (b equal-div-two a
) (u alike1
,var
))
1113 ((coeffpp) (z2 zerp
)))))
1114 ((coeffpp) (d true
)))))
1116 ;; Case 4: c * t^-1 * (1 - 2 * %e^(-a*t) + %e^(2*a*t))
1117 (defun m2-sum-with-exp-case4 (expr var
)
1122 ((mexpt) (u alike1
,var
) -
1))
1125 ((mexpt) (u alike1
,var
) -
1)
1128 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1129 ((coeffpp) (z1 zerp
)))))
1131 (c3 equal-times-minus-two c
)
1132 ((mexpt) (u alike1
,var
) -
1)
1135 ((coeffpt) (b equal-div-two a
) (u alike1
,var
))
1136 ((coeffpp) (z2 zerp
)))))
1137 ((coeffpp) (d true
)))))
1139 ;; Case 5: c* t^-1 * (1 - %e^(2*a*t))
1140 (defun m2-sum-with-exp-case5 (expr var
)
1145 ((mexpt) (u alike1
,var
) -
1))
1147 (c2 equal-times-minus-one c
)
1148 ((mexpt) (u alike1
,var
) -
1)
1151 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1152 ((coeffpp) (z1 zerp
)))))
1153 ((coeffpp) (d true
)))))
1155 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1157 ;;; Test functions for the pattern m2-sum-with-exp-case<n>
1159 (defun equal-times-minus-one (a b
)
1160 (equal a
(mul -
1 b
)))
1162 (defun equal-times-minus-two (a b
)
1163 (equal a
(mul -
2 b
)))
1165 (defun equal-div-two (a b
)
1166 (equal a
(div b
2)))
1168 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1170 ;;; Called by defintegrate.
1171 ;;; Call for every term of a sum defexec and add up the results.
1173 ;; Evaluate the transform of a sum as sum of transforms.
1174 (defun distrdefexecinit (expr var
)
1175 (cond ((equal (caar expr
) 'mplus
)
1176 (distrdefexec (cdr expr
) var
))
1177 (t (defexec expr var
))))
1179 ;; FUN is a list of addends. Compute the transform of each addend and
1181 (defun distrdefexec (expr var
)
1182 (cond ((null expr
) 0)
1183 (t (add (defexec (car expr
) var
)
1184 (distrdefexec (cdr expr
) var
)))))
1186 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1188 ;;; Called after transformation of a integrand to a new representation.
1189 ;;; Evalutate the tansform of a sum as sum of transforms.
1191 (defun sendexec (r a
)
1192 (distrexecinit ($expand
(mul (init r
) a
))))
1194 ;; Compute r*exp(-p*t), where t is the variable of integration and
1195 ;; p is the parameter of the Laplace transform.
1197 (mul r
(power '$%e
(mul -
1 *var
* *par
*))))
1199 (defun distrexecinit (expr)
1200 (cond ((and (consp expr
)
1202 (equal (caar expr
) 'mplus
))
1203 (distrexec (cdr expr
)))
1204 (t (hypgeo-exec expr
))))
1206 (defun distrexec (expr)
1207 (cond ((null expr
) 0)
1208 (t (add (hypgeo-exec (car expr
))
1209 (distrexec (cdr expr
))))))
1211 ;; It dispatches according to the kind of transform it matches.
1212 (defun hypgeo-exec (expr)
1214 (cond ((setq l
(m2-ltep expr
*var
* *par
*))
1215 (setq u
(cdras 'u l
)
1220 (return (ltscale u c a e f
))))
1221 (return (setq *hyp-return-noun-flag
* 'other-trans-to-follow
))))
1223 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1225 ;;; Compute transform of EXP wrt the variable of integration VAR.
1227 (defun defexec (expr var
)
1228 (let* ((*par
* 'psey
) ; Set parameter of Laplace transform
1229 (*var
* var
) ; Set variable of integration
1230 (*hyp-return-noun-flag
* nil
) ; Reset the flag
1232 (l (m2-defltep expr var
))
1233 (s (cdras 'a l
))) ; Get the parameter of the Laplace transform.
1235 ;; If we have not found a parameter, we try to factor the integrand.
1237 (when (and (numberp s
) (equal s
0))
1238 (setq l
(m2-defltep ($factor form
) *var
*))
1239 (setq s
(cdras 'a l
)))
1242 ;; EXP is an expression of the form u*%e^(s*t+e*f+c). So s
1243 ;; is basically the parameter of the Laplace transform.
1244 (let ((result (negtest l s
)))
1245 ;; At this point we construct the noun form if one of the
1246 ;; called routines set the global flag. If the global flag
1247 ;; is not set, the noun form has been already constructed.
1248 (if (and *hyp-return-noun-form-p
* *hyp-return-noun-flag
*)
1249 (list '(%specint simp
) expr
*var
*)
1252 ;; If necessary we construct the noun form.
1253 (if *hyp-return-noun-form-p
*
1254 (list '(%specint simp
) expr
*var
*)
1255 'other-defint-to-follow-defexec
)))))
1257 ;; L is the integrand of the transform, after pattern matching. S is
1258 ;; the parameter (p) of the transform.
1259 (defun negtest (l s
)
1261 (cond ((eq ($asksign
($realpart s
)) '$neg
)
1262 ;; The parameter of transform must have a negative
1263 ;; realpart. Break out the integrand into its various
1265 (setq u
(cdras 'u l
)
1269 (when (equal e
0) (setq f
1))
1270 ;; To compute the transform, we replace A with PSEY for
1271 ;; simplicity. After the transform is computed, replace
1274 ;; FIXME: Sometimes maxima will ask for the sign of PSEY.
1275 ;; But that doesn't occur in the original expression, so
1276 ;; it's very confusing. What should we do?
1278 ;; We know psey must be positive. psey is a substitution
1279 ;; for the paratemter a and we have checked the sign.
1280 ;; So it is the best to add a rule for the sign of psey.
1282 (mfuncall '$assume
`((mgreaterp) ,*par
* 0))
1289 (ltscale u c
0 e f
))
1291 ;; We forget the rule after finishing the calculation.
1292 (mfuncall '$forget
`((mgreaterp) ,*par
* 0))))))
1295 (setq *hyp-return-noun-flag
* 'other-defint-to-follow-negtest
))))
1297 ;; Compute the transform of
1299 ;; U * %E^(-VAR * (*PAR* - PAR0) + E*F + C)
1300 (defun ltscale (u c par0 e f
)
1302 (substl (sub *par
* par0
) *par
* (lt-exec u e f
))))
1304 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1306 ;;; Compute the transform of u*%e^(-p*t+e*f)
1308 (defun lt-exec (u e f
)
1310 (cond ((setq l
(m2-sum u
*var
*))
1311 ;; We have found a summation.
1314 (sendexec 1 (cdras 'u l
))
1319 ((setq l
(m2-unit_step u
*var
*))
1320 ;; We have found the Unit Step function.
1321 (setq u
(cdras 'u l
)
1323 (mul (power '$%e
(mul a
*par
*))
1324 (sendexec (cond (($freeof
*var
* u
) u
)
1325 (t (maxima-substitute (sub *var
* a
) *var
* u
)))
1329 ;; The simple case of u*%e^(-p*t)
1331 ((and (not (equal e
0))
1332 (setq l
(m2-c*t^v u
*var
*)))
1333 ;; We have u*%e^(-p*t+e*f). Try to see if U is of the form
1334 ;; c*t^v. If so, we can handle it here.
1337 ;; The complicated case. Remove the e*f term and move it to u.
1338 (lt-sf-log (mul u
(power '$%e
(mul e f
))))))))
1340 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1342 ;;; Pattern for the routine lt-exec
1344 ;; Recognize c*sum(u,index,low,high)
1345 (defun m2-sum (expr var
)
1350 ((%sum
) (u true
) (i true
) (l true
) (h true
)))
1351 ((coeffpp) (d zerp
)))))
1353 ;; Recognize u(t)*unit_step(x-a)
1354 (defun m2-unit_step (expr var
)
1359 (($unit_step
) ((mplus) (x alike1
,var
) ((coeffpp) (a true
)))))
1360 ((coeffpp) (d zerp
)))))
1363 ;; This is a duplicate of m2-arbpow1. Look if we can use it.
1364 (defun m2-c*t^v
(expr var
)
1367 ((coefftt) (c free
,var
))
1368 ((mexpt) (u alike1
,var
) (v free
,var
)))))
1370 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1372 ;;; Algorithm 1: Laplace transform of c*t^v*exp(-s*t+e*f)
1374 ;;; L contains the pattern for c*t^v.
1375 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1377 (defun lt-exp (l e f
)
1378 (let ((c (cdras 'c l
))
1380 (cond ((m2-t^
2 f
*var
*)
1381 (setq e
(inv (mul -
8 e
)) v
(add v
1))
1382 (f24p146test c v e
))
1383 ((m2-sqroott f
*var
*)
1384 ;; We don't do the transformation at this place. Because we take the
1385 ;; square of e we lost the sign and get wrong results.
1386 ;(setq e (mul* e e (inv 4)) v (add v 1))
1387 (f35p147test c v e
))
1389 (setq e
(mul -
4 e
) v
(add v
1))
1390 (f29p146test c v e
)) ; We have to call with the constant c.
1391 ((and (equal v
0) ; We have to test for v=0 and to call
1393 (f36p147 c e
)) ; with the constant c.
1394 ((and (equal v
0) (m2-e^t f
*var
*))
1395 (f37p147 c
(mul -
1 e
)))
1397 (setq *hyp-return-noun-flag
* 'other-lt-exponential-to-follow
)))))
1399 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1401 ;;; Pattern for the routine lt-exp
1404 (defun m2-t^
2 (expr var
)
1405 (m2 expr
`((mexpt) (u alike1
,var
) 2)))
1407 ;; Recognize sqrt(t)
1408 (defun m2-sqroott (expr var
)
1409 (m2 expr
`((mexpt) (u alike1
,var
) ((rat) 1 2))))
1412 (defun m2-t^-
1 (expr var
)
1413 (m2 expr
`((mexpt) (u alike1
,var
) -
1)))
1416 (defun m2-e^-t
(expr var
)
1417 (m2 expr
`((mexpt) $%e
((mtimes) -
1 (u alike1
,var
)))))
1420 (defun m2-e^t
(expr var
)
1421 (m2 expr
`((mexpt) $%e
(u alike1
,var
))))
1423 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1425 ;;; Algorithm 1.1: Laplace transform of c*t^v*exp(-a*t^2)
1427 ;;; Table of Integral Transforms
1429 ;;; p. 146, formula 24:
1431 ;;; t^(v-1)*exp(-t^2/8/a)
1432 ;;; -> gamma(v)*2^v*a^(v/2)*exp(a*p^2)*D[-v](2*p*sqrt(a))
1434 ;;; Re(a) > 0, Re(v) > 0
1435 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1437 (defun f24p146test (c v a
)
1438 (cond ((and (eq ($asksign a
) '$pos
)
1439 (eq ($asksign v
) '$pos
))
1440 ;; Both a and v must be positive
1443 (setq *hyp-return-noun-flag
* 'fail-on-f24p146test
))))
1445 (defun f24p146 (c v a
)
1450 (power '$%e
(mul a
*par
* *par
*))
1451 (dtford (mul 2 *par
* (power a
'((rat simp
) 1 2)))
1454 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1456 ;;; Algorithm 1.2: Laplace transform of c*t^v*exp(-a*sqrt(t))
1458 ;;; Table of Integral Transforms
1460 ;;; p. 147, formula 35:
1462 ;;; (2*t)^(v-1)*exp(-2*sqrt(a)*sqrt(t))
1463 ;;; -> gamma(2*v)*p^(-v)*exp(a/p/2)*D[-2*v](sqrt(2*a/p))
1465 ;;; Re(v) > 0, Re(p) > 0
1466 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1468 ;; Check if conditions for f35p147 hold
1469 (defun f35p147test (c v a
)
1470 (cond ((eq ($asksign
(add v
1)) '$pos
)
1471 ;; v must be positive
1474 ;; Set a global flag. When *hyp-return-noun-form-p* is T the noun
1475 ;; form will be constructed in the routine DEFEXEC.
1476 (setq *hyp-return-noun-flag
* 'fail-on-f35p147test
))))
1478 (defun f35p147 (c v a
)
1479 ;; We have not done the calculation v->v+1 and a-> a^2/4
1480 ;; and subsitute here accordingly.
1481 (let ((v (add v
1)))
1483 (take '(%gamma
) (add v v
))
1484 (power 2 (sub 1 v
)) ; Is this supposed to be here?
1485 (power *par
* (mul -
1 v
))
1486 (power '$%e
(mul a a
'((rat simp
) 1 8) (inv *par
*)))
1487 ;; We need an additional factor -1 to get the expected results.
1488 ;; What is the mathematically reason?
1489 (dtford (mul -
1 a
(inv (power (mul 2 *par
*) '((rat simp
) 1 2))))
1492 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1494 ;; Express a parabolic cylinder function as either a parabolic
1495 ;; cylinder function or as hypergeometric function.
1497 ;; Tables of Integral Transforms, p. 386 gives this definition:
1499 ;; D[v](z) = 2^(v/2+1/4)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
1502 (let ((inv4 (inv 4)))
1503 (cond ((or $prefer_d
1504 (whittindtest (add (div v
2) inv4
) inv4
))
1505 ;; At this time the Parabolic Cylinder D function is not implemented
1506 ;; as a simplifying function. We call nevertheless the simplifer
1507 ;; to simplify the arguments. When we implement the function
1508 ;; The symbol has to be changed to the noun form.
1509 (take '($parabolic_cylinder_d
) v z
))
1510 (t (simpdtf z v
)))))
1512 (defun whittindtest (i1 i2
)
1513 (or (maxima-integerp (add i2 i2
))
1514 (neginp (sub (sub '((rat simp
) 1 2) i2
) i1
))
1515 (neginp (sub (add '((rat simp
) 1 2) i2
) i1
))))
1517 ;; Return T if a is a non-positive integer.
1518 ;; (Do we really want maxima-integerp or hyp-integerp here?)
1520 (cond ((maxima-integerp a
)
1522 (eq ($sign a
) '$neg
)))))
1524 ;; Express parabolic cylinder function as a hypergeometric function.
1528 ;; U(a,x) = D[-a-1/2](x)
1530 ;; and A&S 19.12.3 gives
1532 ;; U(a,+/-x) = sqrt(%pi)*2^(-1/4-a/2)*exp(-x^2/4)/gamma(3/4+a/2)
1533 ;; *M(a/2+1/4,1/2,x^2/2)
1534 ;; -/+ sqrt(%pi)*2^(1/4-a/2)*x*exp(-x^2/4)/gamma(1/4+a/2)
1535 ;; *M(a/2+3/4,3/2,x^2/2)
1539 ;; D[v](z) = U(-v-1/2,z)
1540 ;; = sqrt(%pi)*2^(v/2+1/2)*x*exp(-x^2/4)
1541 ;; *M(1/2-v/2,3/2,x^2/2)/gamma(-v/2)
1542 ;; + sqrt(%pi)*2^(v/2)*exp(-x^2/4)/gamma(1/2-v/2)
1543 ;; *M(-v/2,1/2,x^2/2)
1545 (defun simpdtf (z v
)
1546 (let ((inv2 '((rat simp
) 1 2))
1547 (pow (power '$%e
(mul z z
'((rat simp
) -
1 4)))))
1548 (add (mul (power 2 (div (sub v
1) 2))
1550 -
2 (power '$%pi inv2
) ; gamma(-1/2)
1551 (inv (take '(%gamma
) (mul v -
1 inv2
)))
1553 (hgfsimp-exec (list (sub inv2
(div v
2)))
1554 (list '((rat simp
) 3 2))
1556 (mul (power 2 (div v
2))
1557 (power '$%pi inv2
) ; gamma(1/2)
1559 (inv (take '(%gamma
) (sub inv2
(mul v inv2
))))
1560 (hgfsimp-exec (list (mul v -
1 inv2
))
1564 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1566 ;;; Algorithm 1.3: Laplace transform of t^v*exp(1/t)
1568 ;;; Table of Integral Transforms
1570 ;;; p. 146, formula 29:
1572 ;;; t^(v-1)*exp(-a/t/4)
1573 ;;; -> 2*(a/p/4)^(v/2)*bessel_k(v, sqrt(a)*sqrt(p))
1576 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1578 ;; Check if conditions for f29p146test hold
1579 (defun f29p146test (c v a
)
1580 (cond ((eq ($asksign a
) '$pos
)
1583 (setq *hyp-return-noun-flag
* 'fail-on-f29p146test
))))
1585 (defun f29p146 (c v a
)
1587 (power (mul a
'((rat simp
) 1 4) (inv *par
*))
1591 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1593 ;; bessel_k(v, sqrt(a)*sqrt(p)) in terms of bessel_k or in terms of
1594 ;; hypergeometric functions.
1596 ;; Choose bessel_k if the order v is an integer. (Why?)
1599 (let ((z (power (mul a
*par
*) '((rat simp
) 1 2))))
1600 (cond ((maxima-integerp v
)
1601 (take '(%bessel_k
) v z
))
1605 ;; Express the Bessel K function in terms of hypergeometric functions.
1607 ;; K[v](z) = %pi/2*(bessel_i(-v,z)-bessel(i,z))/sin(v*%pi)
1611 ;; bessel_i(v,z) = (z/2)^v/gamma(v+1)*0F1(;v+1;z^2/4)
1613 (defun simpktf (z v
)
1614 (let ((dz2 (div z
2)))
1617 (inv (take '(%sin
) (mul v
'$%pi
)))
1618 (sub (mul (power dz2
(mul -
1 v
))
1619 (inv (take '(%gamma
) (sub 1 v
)))
1622 (mul z z
'((rat simp
) 1 4))))
1624 (inv (take '(%gamma
) (add v
1)))
1627 (mul z z
'((rat simp
) 1 4))))))))
1629 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1631 ;;; Algorithm 1.4: Laplace transform of exp(exp(-t))
1633 ;;; Table of Integral Transforms
1635 ;;; p. 147, formula 36:
1638 ;;; -> a^(-p)*gamma(p,a)
1639 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1641 (defun f36p147 (c a
)
1642 (let ((-a (mul -
1 a
)))
1644 (power -a
(mul -
1 *par
*))
1645 `((%gamma_incomplete_lower simp
) ,*par
* ,-a
))))
1647 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1649 ;;; Algorithm 1.5: Laplace transform of exp(exp(t))
1651 ;;; Table of Integral Transforms
1653 ;;; p. 147, formula 36:
1656 ;;; -> a^(-p)*gamma_incomplete(-p,a)
1657 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1659 (defun f37p147 (c a
)
1662 (take '(%gamma_incomplete
) (mul -
1 *par
*) a
)))
1664 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1666 ;;; Algorithm 2: Laplace transform of u(t)*%e^(-p*t).
1668 ;;; u contains one or more special functions. Dispatches according to the
1669 ;;; special functions involved in the Laplace transformable expression.
1671 ;;; We have three general types of integrands:
1673 ;;; 1. Call a function to return immediately the Laplace transform,
1674 ;;; e.g. call lt-arbpow, lt-arbpow2, lt-log, whittest to return the
1675 ;;; Laplace transform.
1676 ;;; 2. Call lt-ltp directly or via an "expert function on Laplace transform",
1677 ;;; transform the special function to a representation in terms of one
1678 ;;; hypergeometric function and do the integration
1679 ;;; e.g. for a direct call of lt-ltp asin, atan or via lt2j for
1680 ;;; an integrand with two bessel function.
1681 ;;; 3. Call fractest, fractest1, ... which transform the involved special
1682 ;;; function to a new representation. Send the transformed expression with
1683 ;;; the routine sendexec to the integrator and try to integrate the new
1684 ;;; representation, e.g. gamma_incomplete is first transformed to a new
1687 ;;; The ordering of the calls to match a pattern is important.
1688 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1690 (defun lt-sf-log (u)
1691 (prog (l index1 index11 index2 index21 arg1 arg2 rest
)
1693 ;; Laplace transform of asin(w)
1694 (cond ((setq l
(m2-asin u
*var
*))
1695 (setq arg1
(cdras 'w l
)
1697 (return (lt-ltp 'asin rest arg1 nil
))))
1699 ;; Laplace transform of atan(w)
1700 (cond ((setq l
(m2-atan u
*var
*))
1701 (setq arg1
(cdras 'w l
)
1703 (return (lt-ltp 'atan rest arg1 nil
))))
1705 ;; Laplace transform of two Bessel J functions
1706 (cond ((setq l
(m2-twoj u
*var
*))
1707 (setq index1
(cdras 'v1 l
)
1708 index2
(cdras 'v2 l
)
1712 (return (lt2j rest arg1 arg2 index1 index2
))))
1714 ;; Laplace transform of two hankel_1 functions
1715 (cond ((setq l
(m2-two-hankel_1 u
*var
*))
1716 (setq index1
(cdras 'v1 l
)
1717 index2
(cdras 'v2 l
)
1721 ;; Call the code for the symbol %h.
1722 (return (fractest rest arg1 arg2 index1
1 index2
1 '2htjory
))))
1724 ;; Laplace transform of two hankel_2 functions
1725 (cond ((setq l
(m2-two-hankel_2 u
*var
*))
1726 (setq index1
(cdras 'v1 l
)
1727 index2
(cdras 'v2 l
)
1731 ;; Call the code for the symbol %h.
1732 (return (fractest rest arg1 arg2 index1
2 index2
2 '2htjory
))))
1734 ;; Laplace transform of hankel_1 * hankel_2
1735 (cond ((setq l
(m2-hankel_1*hankel_2 u
*var
*))
1736 (setq index1
(cdras 'v1 l
)
1737 index2
(cdras 'v2 l
)
1741 ;; Call the code for the symbol %h.
1742 (return (fractest rest arg1 arg2 index1
1 index2
2 '2htjory
))))
1744 ;; Laplace transform of two Bessel Y functions
1745 (cond ((setq l
(m2-twoy u
*var
*))
1746 (setq index1
(cdras 'v1 l
)
1747 index2
(cdras 'v2 l
)
1751 (return (fractest rest arg1 arg2 index1 nil index2 nil
'2ytj
))))
1753 ;; Laplace transform of two Bessel K functions
1754 (cond ((setq l
(m2-twok u
*var
*))
1755 (setq index1
(cdras 'v1 l
)
1756 index2
(cdras 'v2 l
)
1760 (return (fractest rest arg1 arg2 index1 nil index2 nil
'2kti
))))
1762 ;; Laplace transform of Bessel K and Bessel Y functions
1763 (cond ((setq l
(m2-onekoney u
*var
*))
1764 (setq index1
(cdras 'v1 l
)
1765 index2
(cdras 'v2 l
)
1769 (return (fractest rest arg1 arg2 index1 nil index2 nil
'ktiytj
))))
1771 ;; Laplace transform of Bessel I and Bessel J functions
1772 (cond ((setq l
(m2-oneionej u
*var
*))
1773 (setq index1
(cdras 'v1 l
)
1774 index2
(cdras 'v2 l
)
1775 index21
(cdras 'v21 l
)
1776 arg1
(mul '$%i
(cdras 'w1 l
))
1778 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1779 (return (lt2j rest arg1 arg2 index1 index2
))))
1781 ;; Laplace transform of Bessel I and Hankel 1 functions
1782 (cond ((setq l
(m2-bessel_i*hankel_1 u
*var
*))
1783 (setq index1
(cdras 'v1 l
)
1784 index2
(cdras 'v2 l
)
1785 arg1
(mul '$%i
(cdras 'w1 l
))
1787 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1788 (return (fractest1 rest arg1 arg2 index1 index2
1 'besshtjory
))))
1790 ;; Laplace transform of Bessel I and Hankel 2 functions
1791 (cond ((setq l
(m2-bessel_i*hankel_2 u
*var
*))
1792 (setq index1
(cdras 'v1 l
)
1793 index2
(cdras 'v2 l
)
1794 arg1
(mul '$%i
(cdras 'w1 l
))
1796 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1797 (return (fractest1 rest arg1 arg2 index1 index2
2 'besshtjory
))))
1799 ;; Laplace transform of Bessel Y and Bessel J functions
1800 (cond ((setq l
(m2-oneyonej u
*var
*))
1801 (setq index1
(cdras 'v1 l
)
1802 index2
(cdras 'v2 l
)
1806 (return (fractest1 rest arg2 arg1 index2 index1 nil
'bessytj
))))
1808 ;; Laplace transform of Bessel K and Bessel J functions
1809 (cond ((setq l
(m2-onekonej u
*var
*))
1810 (setq index1
(cdras 'v1 l
)
1811 index2
(cdras 'v2 l
)
1815 (return (fractest1 rest arg2 arg1 index2 index1 nil
'besskti
))))
1817 ;; Laplace transform of Hankel 1 and Bessel J functions
1818 (cond ((setq l
(m2-hankel_1*bessel_j u
*var
*))
1819 (setq index1
(cdras 'v1 l
)
1820 index2
(cdras 'v2 l
)
1824 (return (fractest1 rest arg2 arg1 index2 index1
1 'besshtjory
))))
1826 ;; Laplace transform of Hankel 2 and Bessel J functions
1827 (cond ((setq l
(m2-hankel_2*bessel_j u
*var
*))
1828 (setq index1
(cdras 'v1 l
)
1829 index2
(cdras 'v2 l
)
1833 (return (fractest1 rest arg2 arg1 index2 index1
2 'besshtjory
))))
1835 ;; Laplace transform of Bessel Y and Hankel 1 functions
1836 (cond ((setq l
(m2-bessel_y*hankel_1 u
*var
*))
1837 (setq index1
(cdras 'v1 l
)
1838 index2
(cdras 'v2 l
)
1842 (return (fractest1 rest arg2 arg1 index2 index1
1 'htjoryytj
))))
1844 ;; Laplace transform of Bessel Y and Hankel 2 functions
1845 (cond ((setq l
(m2-bessel_y*hankel_2 u
*var
*))
1846 (setq index1
(cdras 'v1 l
)
1847 index2
(cdras 'v2 l
)
1851 (return (fractest1 rest arg2 arg1 index2 index1
2 'htjoryytj
))))
1853 ;; Laplace transform of Bessel K and Hankel 1 functions
1854 (cond ((setq l
(m2-bessel_k*hankel_1 u
*var
*))
1855 (setq index1
(cdras 'v1 l
)
1856 index2
(cdras 'v2 l
)
1860 (return (fractest1 rest arg2 arg1 index2 index1
1 'htjorykti
))))
1862 ;; Laplace transform of Bessel K and Hankel 2 functions
1863 (cond ((setq l
(m2-bessel_k*hankel_2 u
*var
*))
1864 (setq index1
(cdras 'v1 l
)
1865 index2
(cdras 'v2 l
)
1869 (return (fractest1 rest arg2 arg1 index2 index1
2 'htjorykti
))))
1871 ;; Laplace transform of Bessel I and Bessel Y functions
1872 (cond ((setq l
(m2-oneioney u
*var
*))
1873 (setq index1
(cdras 'v1 l
)
1874 index2
(cdras 'v2 l
)
1875 arg1
(mul '$%i
(cdras 'w1 l
))
1877 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1878 (return (fractest1 rest arg1 arg2 index1 index2 nil
'bessytj
))))
1880 ;; Laplace transform of Bessel I and Bessel K functions
1881 (cond ((setq l
(m2-oneionek u
*var
*))
1882 (setq index1
(cdras 'v1 l
)
1883 index2
(cdras 'v2 l
)
1884 arg1
(mul '$%i
(cdras 'w1 l
))
1886 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1887 (return (fractest1 rest arg1 arg2 index1 index2 nil
'besskti
))))
1889 ;; Laplace transform of Struve H function
1890 (cond ((setq l
(m2-struve_h u
*var
*))
1891 (setq index1
(cdras 'v l
)
1894 (return (lt1hstruve rest arg1 index1
))))
1896 ;; Laplace transform of Struve L function
1897 (cond ((setq l
(m2-struve_l u
*var
*))
1898 (setq index1
(cdras 'v l
)
1901 (return (lt1lstruve rest arg1 index1
))))
1903 ;; Laplace transform of little Lommel s function
1904 (cond ((setq l
(m2-ones u
*var
*))
1905 (setq index1
(cdras 'v1 l
)
1906 index2
(cdras 'v2 l
)
1909 (return (lt1s rest arg1 index1 index2
))))
1911 ;; Laplace transform of Lommel S function
1912 (cond ((setq l
(m2-oneslommel u
*var
*))
1913 (setq index1
(cdras 'v1 l
)
1914 index2
(cdras 'v2 l
)
1917 (return (fractest2 rest arg1 index1 index2
'slommel
))))
1919 ;; Laplace transform of Bessel Y function
1920 (cond ((setq l
(m2-oney u
*var
*))
1921 (setq index1
(cdras 'v l
)
1924 (return (lt1yref rest arg1 index1
))))
1926 ;; Laplace transform of Bessel K function
1927 (cond ((setq l
(m2-onek u
*var
*))
1928 (setq index1
(cdras 'v l
)
1931 (cond ((zerop1 index1
)
1932 ;; Special case for a zero index
1933 (return (lt-bessel_k0 rest arg1
)))
1935 (return (fractest2 rest arg1 index1 nil
'kti
))))))
1937 ;; Laplace transform of Parabolic Cylinder function
1938 (cond ((setq l
(m2-parabolic_cylinder_d u
*var
*))
1939 (setq index1
(cdras 'v l
)
1942 (return (fractest2 rest arg1 index1 nil
'd
))))
1944 ;; Laplace transform of Incomplete Gamma function
1945 (cond ((setq l
(m2-onegammaincomplete u
*var
*))
1946 (setq arg1
(cdras 'w1 l
)
1949 (return (fractest2 rest arg1 arg2 nil
'gamma_incomplete
))))
1951 ;; Laplace transform of Batemann function
1952 (cond ((setq l
(m2-onekbateman u
*var
*))
1953 (setq index1
(cdras 'v l
)
1956 (return (fractest2 rest arg1 index1 nil
'kbateman
))))
1958 ;; Laplace transform of Bessel J function
1959 (cond ((setq l
(m2-onej u
*var
*))
1960 (setq index1
(cdras 'v l
)
1963 (return (lt1j rest arg1 index1
))))
1965 ;; Laplace transform of lower incomplete Gamma function
1966 (cond ((setq l
(m2-onegamma-incomplete-lower u
*var
*))
1967 (setq arg1
(cdras 'w1 l
)
1970 (return (lt1gamma-incomplete-lower rest arg1 arg2
))))
1972 ;; Laplace transform of Hankel 1 function
1973 (cond ((setq l
(m2-hankel_1 u
*var
*))
1974 (setq index1
(cdras 'v l
)
1977 (return (fractest2 rest arg1 index1
1 'htjory
))))
1979 ;; Laplace transform of Hankel 2 function
1980 (cond ((setq l
(m2-hankel_2 u
*var
*))
1981 (setq index1
(cdras 'v l
)
1984 (return (fractest2 rest arg1 index1
2 'htjory
))))
1986 ;; Laplace transform of Whittaker M function
1987 (cond ((setq l
(m2-onem u
*var
*))
1988 (setq index1
(cdras 'v1 l
)
1989 index11
(cdras 'v2 l
)
1992 (return (lt1m rest arg1 index1 index11
))))
1994 ;; Laplace transform of Whittaker M function
1995 (cond ((setq l
(m2-whittaker_m u
*var
*))
1996 (setq index1
(cdras 'v1 l
)
1997 index2
(cdras 'v2 l
)
2000 (return (lt1m rest arg1 index1 index2
))))
2002 ;; Laplace transform of the Generalized Laguerre function, %l[v1,v2](w)
2003 (cond ((setq l
(m2-onel u
*var
*))
2004 (setq index1
(cdras 'v1 l
)
2005 index11
(cdras 'v2 l
)
2008 (return (integertest rest arg1 index1 index11
'l
))))
2010 ;; Laplace transform for the Generalized Laguerre function
2011 ;; We call the routine for %l[v1,v2](w).
2012 (cond ((setq l
(m2-one-gen-laguerre u
*var
*))
2013 (setq index1
(cdras 'v1 l
)
2014 index2
(cdras 'v2 l
)
2017 (return (integertest rest arg1 index1 index2
'l
))))
2019 ;; Laplace transform for the Laguerre function
2020 ;; We call the routine for %l[v1,0](w).
2021 (cond ((setq l
(m2-one-laguerre u
*var
*))
2022 (setq index1
(cdras 'v1 l
)
2025 (return (integertest rest arg1 index1
0 'l
))))
2027 ;; Laplace transform of Gegenbauer function
2028 (cond ((setq l
(m2-onec u
*var
*))
2029 (setq index1
(cdras 'v1 l
)
2030 index11
(cdras 'v2 l
)
2033 (return (integertest rest arg1 index1 index11
'c
))))
2035 ;; Laplace transform of Chebyshev function of the first kind
2036 (cond ((setq l
(m2-onet u
*var
*))
2037 (setq index1
(cdras 'v1 l
)
2040 (return (integertest rest arg1 index1 nil
't
))))
2042 ;; Laplace transform of Chebyshev function of the second kind
2043 (cond ((setq l
(m2-oneu u
*var
*))
2044 (setq index1
(cdras 'v1 l
)
2047 (return (integertest rest arg1 index1 nil
'u
))))
2049 ;; Laplace transform for the Hermite function, hermite(index1,arg1)
2050 (cond ((setq l
(m2-one-hermite u
*var
*))
2051 (setq index1
(cdras 'v1 l
)
2055 (cond ((and (maxima-integerp index1
)
2058 ;; When index1 is an even or odd integer, we transform
2059 ;; directly to a hypergeometric function. For this case we
2060 ;; get a Laplace transform when the arg is the
2061 ;; square root of the variable.
2062 (sendexec rest
(hermite-to-hypergeometric index1 arg1
)))
2064 (integertest rest arg1 index1 nil
'he
))))))
2066 ;; Laplace transform of %p[v1,v2](w), Associated Legendre P function
2067 (cond ((setq l
(m2-hyp-onep u
*var
*))
2068 (setq index1
(cdras 'v1 l
)
2069 index11
(cdras 'v2 l
)
2072 (return (lt1p rest arg1 index1 index11
))))
2074 ;; Laplace transform of Associated Legendre P function
2075 (cond ((setq l
(m2-assoc_legendre_p u
*var
*))
2076 (setq index1
(cdras 'v1 l
)
2077 index2
(cdras 'v2 l
)
2080 (return (lt1p rest arg1 index1 index2
))))
2082 ;; Laplace transform of %p[v1,v2,v3](w), Jacobi function
2083 (cond ((setq l
(m2-onepjac u
*var
*))
2084 (setq index1
(cdras 'v1 l
)
2085 index2
(cdras 'v2 l
)
2086 index21
(cdras 'v3 l
)
2089 (return (pjactest rest arg1 index1 index2 index21
))))
2091 ;; Laplace transform of Jacobi P function
2092 (cond ((setq l
(m2-jacobi_p u
*var
*))
2093 (setq index1
(cdras 'v1 l
)
2094 index2
(cdras 'v2 l
)
2095 index21
(cdras 'v3 l
)
2098 (return (pjactest rest arg1 index1 index2 index21
))))
2100 ;; Laplace transform of Associated Legendre function of the second kind
2101 (cond ((setq l
(m2-oneq u
*var
*))
2102 (setq index1
(cdras 'v1 l
)
2103 index11
(cdras 'v2 l
)
2106 (return (lt1q rest arg1 index1 index11
))))
2108 ;; Laplace transform of Associated Legendre function of the second kind
2109 (cond ((setq l
(m2-assoc_legendre_q u
*var
*))
2110 (setq index1
(cdras 'v1 l
)
2111 index2
(cdras 'v2 l
)
2114 (return (lt1q rest arg1 index1 index2
))))
2116 ;; Laplace transform of %p[v1](w), Legendre P function
2117 (cond ((setq l
(m2-onep0 u
*var
*))
2118 (setq index1
(cdras 'v1 l
)
2122 (return (lt1p rest arg1 index1 index11
))))
2124 ;; Laplace transform of Legendre P function
2125 (cond ((setq l
(m2-legendre_p u
*var
*))
2126 (setq index1
(cdras 'v1 l
)
2129 (return (lt1p rest arg1 index1
0))))
2131 ;; Laplace transform of Whittaker W function
2132 (cond ((setq l
(m2-onew u
*var
*))
2133 (setq index1
(cdras 'v1 l
)
2134 index11
(cdras 'v2 l
)
2137 (return (whittest rest arg1 index1 index11
))))
2139 ;; Laplace transform of Whittaker W function
2140 (cond ((setq l
(m2-whittaker_w u
*var
*))
2141 (setq index1
(cdras 'v1 l
)
2142 index2
(cdras 'v2 l
)
2145 (return (whittest rest arg1 index1 index2
))))
2147 ;; Laplace transform of square of Bessel J function
2148 (cond ((setq l
(m2-onej^
2 u
*var
*))
2149 (setq index1
(cdras 'v l
)
2152 (return (lt1j^
2 rest arg1 index1
))))
2154 ;; Laplace transform of square of Hankel 1 function
2155 (cond ((setq l
(m2-hankel_1^
2 u
*var
*))
2156 (setq index1
(cdras 'v l
)
2159 (return (fractest rest arg1 arg1 index1
1 index1
1 '2htjory
))))
2161 ;; Laplace transform of square of Hankel 2 function
2162 (cond ((setq l
(m2-hankel_2^
2 u
*var
*))
2163 (setq index1
(cdras 'v l
)
2166 (return (fractest rest arg1 arg1 index1
2 index1
2 '2htjory
))))
2168 ;; Laplace transform of square of Bessel Y function
2169 (cond ((setq l
(m2-oney^
2 u
*var
*))
2170 (setq index1
(cdras 'v l
)
2173 (return (fractest rest arg1 arg1 index1 nil index1 nil
'2ytj
))))
2175 ;; Laplace transform of square of Bessel K function
2176 (cond ((setq l
(m2-onek^
2 u
*var
*))
2177 (setq index1
(cdras 'v l
)
2180 (return (fractest rest arg1 arg1 index1 nil index1 nil
'2kti
))))
2182 ;; Laplace transform of two Bessel I functions
2183 (cond ((setq l
(m2-twoi u
*var
*))
2184 (setq index1
(cdras 'v1 l
)
2185 index2
(cdras 'v2 l
)
2186 arg1
(mul '$%i
(cdras 'w1 l
))
2187 arg2
(mul '$%i
(cdras 'w2 l
))
2188 rest
(mul (power '$%i
(neg index1
))
2189 (power '$%i
(neg index1
))
2191 (return (lt2j rest arg1 arg2 index1 index2
))))
2193 ;; Laplace transform of Bessel I. We use I[v](w)=%i^n*J[n](%i*w).
2194 (cond ((setq l
(m2-onei u
*var
*))
2195 (setq index1
(cdras 'v l
)
2196 arg1
(mul '$%i
(cdras 'w l
))
2197 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
2198 (return (lt1j rest arg1 index1
))))
2200 ;; Laplace transform of square of Bessel I function
2201 (cond ((setq l
(m2-onei^
2 u
*var
*))
2202 (setq index1
(cdras 'v l
)
2203 arg1
(mul '$%i
(cdras 'w l
))
2204 rest
(mul (power '$%i
(neg index1
))
2205 (power '$%i
(neg index1
))
2207 (return (lt1j^
2 rest arg1 index1
))))
2209 ;; Laplace transform of Erf function
2210 (cond ((setq l
(m2-onerf u
*var
*))
2211 (setq arg1
(cdras 'w l
)
2213 (return (lt1erf rest arg1
))))
2215 ;; Laplace transform of the logarithmic function.
2216 ;; We add an algorithm for the Laplace transform and call the routine
2217 ;; lt-log. The old code is still present, but isn't called.
2218 (cond ((setq l
(m2-onelog u
*var
*))
2219 (setq arg1
(cdras 'w l
)
2221 (return (lt-log rest arg1
))))
2223 ;; Laplace transform of Erfc function
2224 (cond ((setq l
(m2-onerfc u
*var
*))
2225 (setq arg1
(cdras 'w l
)
2227 (return (fractest2 rest arg1 nil nil
'erfc
))))
2229 ;; Laplace transform of expintegral_ei.
2230 ;; Maxima uses the build in transformation to the gamma_incomplete
2231 ;; function and simplifies the log functions of the transformation. We do
2232 ;; not use the dispatch mechanism of fractest2, but call sendexec directly
2233 ;; with the transformed function.
2234 (cond ((setq l
(m2-oneexpintegral_ei u
*var
*))
2235 (setq arg1
(cdras 'w l
)
2237 (let (($expintrep
'%gamma_incomplete
)
2239 (return (sratsimp (sendexec rest
(ftake '%expintegral_ei arg1
)))))))
2241 ;; Laplace transform of expintegral_e1
2242 (cond ((setq l
(m2-oneexpintegral_e1 u
*var
*))
2243 (setq arg1
(cdras 'w l
)
2245 (let (($expintrep
'%gamma_incomplete
)
2247 (return (sratsimp (sendexec rest
(ftake '%expintegral_e1 arg1
)))))))
2249 ;; Laplace transform of expintegral_e
2250 (cond ((setq l
(m2-oneexpintegral_e u
*var
*))
2251 (setq arg1
(cdras 'v l
)
2254 (let (($expintrep
'%gamma_incomplete
)
2256 (return (sratsimp (sendexec rest
(ftake '%expintegral_e arg1 arg2
)))))))
2258 ;; Laplace transform of expintegral_si
2259 (cond ((setq l
(m2-oneexpintegral_si u
*var
*))
2260 (setq arg1
(cdras 'w l
)
2262 ;; We transform to the hypergeometric representation.
2264 (sendexec rest
(expintegral_si-to-hypergeometric arg1
)))))
2266 ;; Laplace transform of expintegral_shi
2267 (cond ((setq l
(m2-oneexpintegral_shi u
*var
*))
2268 (setq arg1
(cdras 'w l
)
2270 ;; We transform to the hypergeometric representation.
2272 (sendexec rest
(expintegral_shi-to-hypergeometric arg1
)))))
2274 ;; Laplace transform of expintegral_ci
2275 (cond ((setq l
(m2-oneexpintegral_ci u
*var
*))
2276 (setq arg1
(cdras 'w l
)
2278 ;; We transform to the hypergeometric representation.
2279 ;; Because we have Logarithmic terms in the transformation,
2280 ;; we switch on the flag $logexpand and do a ratsimp.
2281 (let (($logexpand
'$super
))
2284 (sendexec rest
(expintegral_ci-to-hypergeometric arg1
)))))))
2286 ;; Laplace transform of expintegral_chi
2287 (cond ((setq l
(m2-oneexpintegral_chi u
*var
*))
2288 (setq arg1
(cdras 'w l
)
2290 ;; We transform to the hypergeometric representation.
2291 ;; Because we have Logarithmic terms in the transformation,
2292 ;; we switch on the flag $logexpand and do a ratsimp.
2293 (let (($logexpand
'$super
))
2296 (sendexec rest
(expintegral_chi-to-hypergeometric arg1
)))))))
2298 ;; Laplace transform of Complete elliptic integral of the first kind
2299 (cond ((setq l
(m2-onekelliptic u
*var
*))
2300 (setq arg1
(cdras 'w l
)
2302 (return (lt1kelliptic rest arg1
))))
2304 ;; Laplace transform of Complete elliptic integral of the first kind
2305 (cond ((setq l
(m2-elliptic_kc u
*var
*))
2306 (setq arg1
(cdras 'w l
)
2308 (return (lt1kelliptic rest arg1
))))
2310 ;; Laplace transform of Complete elliptic integral of the second kind
2311 (cond ((setq l
(m2-onee u
*var
*))
2312 (setq arg1
(cdras 'w l
)
2314 (return (lt1e rest arg1
))))
2316 ;; Laplace transform of Complete elliptic integral of the second kind
2317 (cond ((setq l
(m2-elliptic_ec u
*var
*))
2318 (setq arg1
(cdras 'w l
)
2320 (return (lt1e rest arg1
))))
2322 ;; Laplace transform of %f[v1,v2](w1,w2,w3), Hypergeometric function
2323 ;; We support the Laplace transform of the build in symbol %f. We do
2324 ;; not use the mechanism of defining an "Expert on Laplace transform",
2325 ;; the expert function does a call to lt-ltp. We do this call directly.
2326 (cond ((setq l
(m2-onef u
*var
*))
2327 (setq rest
(cdras 'u l
)
2329 index1
(list (cdras 'w1 l
) (cdras 'w2 l
)))
2330 (return (lt-ltp 'f rest arg1 index1
))))
2332 ;; Laplace transform of Hypergeometric function
2333 (cond ((setq l
(m2-hypergeometric u
*var
*))
2334 (setq rest
(cdras 'u l
)
2336 index1
(list (cdras 'w1 l
) (cdras 'w2 l
)))
2337 (return (lt-ltp 'f rest arg1 index1
))))
2339 ;; Laplace transform of c * t^v * (a+t)^w
2340 ;; It is possible to combine arbpow2 and arbpow.
2341 (cond ((setq l
(m2-arbpow2 u
*var
*))
2342 (setq rest
(cdras 'c l
)
2346 index2
(cdras 'w l
))
2347 (return (lt-arbpow2 rest arg1 arg2 index1 index2
))))
2349 ;; Laplace transform of c * t^v
2350 (cond ((setq l
(m2-arbpow1 u
*var
*))
2351 (setq arg1
(cdras 'u l
)
2353 index1
(cdras 'v l
))
2354 (return (mul arg2
(lt-arbpow arg1 index1
)))))
2356 ;; We have specialized the pattern for arbpow1. Now a lot of integrals
2357 ;; will fail correctly and we have to return a noun form.
2358 (return (setq *hyp-return-noun-flag
* 'other-j-cases-next
))))
2360 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2362 ;;; Algorithm 2.1: Laplace transform of c*t^u
2364 ;;; Table of Integral Transforms
2366 ;;; p. 137, formula 1:
2369 ;;; -> gamma(u+1)*p^(-u-1)
2372 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2374 (defun lt-arbpow (expr pow
)
2375 (cond ((or (eq expr
*var
*) (equal pow
0))
2378 (setq *hyp-return-noun-flag
* 'lt-arbow-failed
))))
2380 ;; Check if conditions for f1p137 hold
2381 (defun f1p137test (pow)
2382 (cond ((eq ($asksign
($realpart
(add pow
1))) '$pos
)
2385 (setq *hyp-return-noun-flag
* 'fail-in-arbpow
))))
2388 (mul (take '(%gamma
) (add pow
1))
2389 (power *par
* (sub (mul -
1 pow
) 1))))
2391 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2393 ;;; Algorithm 2.2: Laplace transform of c*t^v*(1+t)^w
2395 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2397 (defun lt-arbpow2 (c a b pow1 pow2
)
2398 (if (eq ($asksign
(add pow1
1)) '$pos
)
2401 ;; The Laplace transform is an Incomplete Gamma function.
2403 (power a
(add pow2
1))
2405 (power (mul *par
* a
(inv b
)) (mul -
1 (add pow2
1)))
2406 (power '$%e
(mul *par
* a
(inv b
)))
2407 (take '(%gamma_incomplete
) (add pow2
1) (mul *par
* a
(inv b
)))))
2408 ((not (maxima-integerp (add pow1 pow2
2)))
2409 ;; The general result is a Hypergeometric U function U(a,b,z) which can
2410 ;; be represented by two Hypergeometic 1F1 functions for the special
2411 ;; case that the index b is not an integer value.
2413 (power a
(add pow1 pow2
1))
2414 (inv (power b
(add pow1
1)))
2415 (take '(%gamma
) (add pow1 pow2
1))
2416 (power (mul *par
* a
(inv b
)) (mul -
1 (add pow1 pow2
1)))
2417 (hgfsimp-exec (list (mul -
1 pow2
))
2418 (list (mul -
1 (add pow1 pow2
)))
2419 (mul *par
* a
(inv b
))))
2421 (power a
(add pow1 pow2
1))
2422 (inv (power b
(add pow1
1)))
2423 (take '(%gamma
) (add pow1
1))
2424 (take '(%gamma
) (mul -
1 (add pow1 pow2
1)))
2425 (inv (take '(%gamma
) (mul -
1 pow2
)))
2426 (hgfsimp-exec (list (add pow1
1))
2427 (list (add pow1 pow2
2))
2428 (mul *par
* a
(inv b
))))))
2430 ;; The most general case is a result with the Hypergeometric U function.
2432 (power a
(add pow1 pow2
1))
2433 (inv (power b
(add pow1
1)))
2434 (take '(%gamma
) (add pow1
1))
2435 (take '($hypergeometric_u
)
2438 (mul *par
* a
(inv b
))))))
2439 (setq *hyp-return-noun-flag
* 'lt-arbpow2-failed
)))
2441 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2443 ;;; Algorithm 2.3: Laplace transform of the Logarithmic function
2445 ;;; c*t^(v-1)*log(a*t)
2446 ;;; -> c*gamma(v)*s^(-v)*(psi[0](v)-log(s/a))
2448 ;;; This is the formula for an expression with log(t) scaled like 1/a*F(s/a).
2450 ;;; For the following cases we have to add further algorithm:
2451 ;;; log(1+a*x), log(x+a), log(x)^2.
2452 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2454 (defun lt-log (rest arg
)
2455 (let* ((l (m2-c*t^v rest
*var
*))
2457 (v (add (cdras 'v l
) 1))) ; because v -> v-1
2459 ((and l
(eq ($asksign v
) '$pos
))
2460 (let* ((l1 (m2-a*t arg
*var
*))
2465 (inv (power *par
* v
))
2466 (sub (take '(mqapply) (list '($psi array
) 0) v
)
2467 (take '(%log
) (div *par
* a
)))))
2469 (setq *hyp-return-noun-flag
* 'lt-log-failed
)))))
2471 (setq *hyp-return-noun-flag
* 'lt-log-failed
)))))
2473 ;; Pattern for lt-log.
2474 ;; Extract the argument of a function: a*t+c for c=0.
2475 (defun m2-a*t
(expr var
)
2478 ((mtimes) (u alike1
,var
) (a free
,var
))
2479 ((coeffpp) (c equal
0)))))
2481 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2482 ;;; Algorithm 2.4: Laplace transfom of the Whittaker function
2484 ;;; Test for Whittaker W function. Simplify this if possible, or
2485 ;;; convert to Whittaker M function.
2487 ;;; We have r * %w[i1,i2](a)
2489 ;;; Formula 16, p. 217
2491 ;;; t^(v-1)*%w[k,u](a*t)
2492 ;;; -> gamma(u+v+1/2)*gamma(v-u+1/2)*a^(u+1/2)/
2493 ;;; (gamma(v-k+1)*(p+a/2)^(u+v+1/2)
2494 ;;; *2f1(u+v+1/2,u-k+1/2;v-k+1;(p-a/2)/(p+a/2))
2496 ;;; For Re(v +/- mu) > -1/2
2497 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2499 (defun whittest (r a i1 i2
)
2501 (cond ((f16p217test r a i1 i2
))
2504 (and (maxima-integerp (setq n
(sub (sub '((rat simp
) 1 2) i2
) i1
)))
2505 (member ($sign n
) '($zero $neg $nz
))))
2507 (and (maxima-integerp (setq m
(sub (add '((rat simp
) 1 2) i2
) i1
)))
2508 (member ($sign m
) '($zero $neg $nz
)))))
2509 ;; 1/2-u-k and 1/2+u-k are not zero or a negative integer
2510 ;; Transform to Whittaker M and try again.
2511 (distrexecinit ($expand
(mul (init r
) (wtm a i1 i2
)))))
2513 ;; Both conditions fails, return a noun form.
2514 (setq *hyp-return-noun-flag
* 'whittest-failed
)))))
2516 (defun f16p217test (r a i1 i2
)
2517 ;; We have r*%w[i1,i2](a)
2518 (let ((l (m2-c*t^v r
*var
*)))
2519 ;; Make sure r is of the form c*t^v
2521 (let* ((v (add (cdras 'v l
) 1))
2523 ;; Check that v + i2 + 1/2 > 0 and v - i2 + 1/2 > 0.
2524 (when (and (eq ($asksign
(add (add v i2
) '((rat simp
) 1 2))) '$pos
)
2525 (eq ($asksign
(add (sub v i2
) '((rat simp
) 1 2))) '$pos
))
2526 ;; Ok, we satisfy the conditions. Now extract the arg.
2527 ;; The transformation is only valid for an argument a*t. We have
2528 ;; to special the pattern to make sure that we satisfy the condition.
2529 (let ((l (m2-a*t a
*var
*)))
2531 (let ((a (cdras 'a l
)))
2532 ;; We're ready now to compute the transform.
2534 (power a
(add i2
'((rat simp
) 1 2)))
2535 (take '(%gamma
) (add (add v i2
) '((rat simp
) 1 2)))
2536 (take '(%gamma
) (add (sub v i2
) '((rat simp
) 1 2)))
2537 (inv (mul (take '(%gamma
) (add (sub v i1
) 1))
2538 (power (add *par
* (div a
2))
2539 (add (add i2 v
) '((rat simp
) 1 2)))))
2540 (hgfsimp-exec (list (add (add i2 v
'((rat simp
) 1 2)))
2541 (add (sub i2 i1
) '((rat simp
) 1 2)))
2542 (list (add (sub v i1
) 1))
2543 (div (sub *par
* (div a
2))
2544 (add *par
* (div a
2)))))))))))))
2546 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2547 ;;; Algorithm 2.5: Laplace transfom of bessel_k(0,a*t)
2549 ;;; The general algorithm handles the Bessel K function for an order |v|<1.
2550 ;;; but does not include the special case v=0. Return the Laplace transform:
2552 ;;; bessel_k(0,a*t) --> acosh(s/a)/sqrt(s^2-a^2)
2554 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2556 (defun lt-bessel_k0 (rest arg
)
2557 (let* ((l (m2-c*t^v rest
*var
*))
2560 (l (m2-a*t arg
*var
*))
2562 (cond ((and l
(zerop1 v
))
2564 (take '(%acosh
) (div *par
* a
))
2565 (inv (power (sub (mul *par
* *par
*) (mul a a
))
2566 '((rat simp
) 1 2)))))
2568 (setq *hyp-return-noun-flag
* 'lt-bessel_k-failed
)))))
2570 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2572 ;;; DISPATCH FUNCTIONS TO CHANGE THE REPRESENTATION OF SPECIAL FUNCTIONS
2574 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2576 ;; Laplace transform of a product of Bessel functions. A1, A2 are
2577 ;; the args of the two functions. I1, I2 are the indices of each
2578 ;; function. I11, I21 are secondary indices of each function, if any.
2579 ;; FLG is a symbol indicating how we should handle the special
2580 ;; functions (and also indicates what the special functions are.)
2582 ;; I11 and I21 are for the Hankel functions.
2584 (defun fractest (r a1 a2 i1 i11 i2 i21 flg
)
2585 (cond ((or (and (listp i1
) (equal (caar i1
) 'rat
)
2586 (listp i2
) (equal (caar i2
) 'rat
))
2588 ;; We have to Bessel or Hankel functions. Both indizes have to be
2589 ;; rational numbers or we have two Hankel functions.
2591 (cond ((eq flg
'2ytj
)
2595 (mul (htjory i1 i11 a1
)
2596 (htjory i2 i21 a2
)))
2604 (setq *hyp-return-noun-flag
* 'product-of-y-with-nofract-indices
))))
2606 ;; Laplace transform of a product of Bessel functions. A1, A2 are
2607 ;; the args of the two functions. I1, I2 are the indices of each
2608 ;; function. I is a secondary index to one function, if any.
2609 ;; FLG is a symbol indicating how we should handle the special
2610 ;; functions (and also indicates what the special functions are.)
2612 ;; I is for the kind of Hankel function.
2614 (defun fractest1 (r a1 a2 i1 i2 i flg
)
2615 (cond ((or (and (listp i2
)
2616 (equal (caar i2
) 'rat
))
2617 (eq flg
'besshtjory
))
2618 ;; We have two Bessel or Hankel functions. The second index has to
2619 ;; be a rational number or one of the functions is a Hankel function
2620 ;; and the second function is Bessel J or Bessel I
2622 (cond ((eq flg
'bessytj
)
2623 (mul (take '(%bessel_j
) i1 a1
)
2625 ((eq flg
'besshtjory
)
2626 (mul (take '(%bessel_j
) i1 a1
)
2628 ((eq flg
'htjoryytj
)
2629 (mul (htjory i1 i a1
)
2632 (mul (take '(%bessel_j
) i1 a1
)
2634 ((eq flg
'htjorykti
)
2635 (mul (htjory i1 i a1
)
2638 (setq *hyp-return-noun-flag
* 'product-of-i-y-of-nofract-index
))))
2640 ;; Laplace transform of a single special function. A is the arg of
2641 ;; the special function. I1, I11 are the indices of the function. FLG
2642 ;; is a symbol indicating how we should handle the special functions
2643 ;; (and also indicates what the special functions are.)
2645 ;; I11 is the kind of Hankel function
2647 (defun fractest2 (r a1 i1 i11 flg
)
2648 (cond ((or (and (listp i1
)
2649 (equal (caar i1
) 'rat
))
2652 (eq flg
'gamma_incomplete
)
2657 ;; The index is a rational number or flg has the value of one of the
2658 ;; above special functions.
2660 (cond ((eq flg
'ytj
)
2668 ((eq flg
'gamma_incomplete
)
2669 (gamma_incomplete-to-gamma-incomplete-lower a1 i1
))
2675 (slommeltjandy i1 i11 a1
)))))
2677 (setq *hyp-return-noun-flag
* 'y-of-nofract-index
))))
2679 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2681 (defun integertest (r arg i1 i2 flg
)
2682 (cond ((maxima-integerp i1
)
2683 (dispatchpoltrans r arg i1 i2 flg
))
2685 (setq *hyp-return-noun-flag
* 'index-should-be-an-integer-in-polys
))))
2687 (defun dispatchpoltrans (r x i1 i2 flg
)
2689 (cond ((eq flg
'l
)(ltw x i1 i2
))
2690 ((eq flg
'he
)(hetd x i1
))
2691 ((eq flg
'c
)(ctpjac x i1 i2
))
2692 ((eq flg
't
)(ttpjac x i1
))
2693 ((eq flg
'u
)(utpjac x i1
)))))
2695 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2696 ;; (-1)^n*n!*laguerre(n,a,x) = U(-n,a+1,x)
2698 ;; W[k,u](z) = exp(-z/2)*z^(u+1/2)*U(1/2+u-k,1+2*u,z)
2702 ;; laguerre(n,a,x) = (-1)^n*U(-n,a+1,x)/n!
2704 ;; U(-n,a+1,x) = exp(z/2)*z^(-a/2-1/2)*W[1/2+a/2+n,a/2](z)
2708 ;; laguerre(n,a,x) = (-1)^n/n!*exp(z/2)*z^(-a/2-1/2)*M[1/2+a/2+n,a/2](z)
2710 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2713 (let ((diva2 (div a
2)))
2714 (mul (take '(%gamma
) (add n a
1))
2715 (inv (take '(%gamma
) (add a
1)))
2716 (inv (take '(%gamma
) (add n
1)))
2717 (power x
(sub '((rat simp
) -
1 2) diva2
))
2718 (power '$%e
(div x
2))
2719 (list '(mqapply simp
)
2720 (list '($%m simp array
)
2721 (add '((rat simp
) 1 2) diva2 n
) diva2
) x
))))
2723 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2724 ;; Hermite He function as a parabolic cylinder function
2726 ;; Tables of Integral Transforms
2730 ;; D[n](z) = (-1)^n*exp(z^2/4)*diff(exp(-z^2/2),z,n);
2734 ;; He[n](x) = (-1)^n*exp(x^2/2)*diff(exp(-x^2/2),x,n)
2735 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2738 (mul (power '$%e
(mul x x
'((rat simp
) 1 4)))
2739 ;; At this time the Parabolic Cylinder D function is not implemented
2740 ;; as a simplifying function. We call nevertheless the simplifer.
2741 (take '($parabolic_cylinder_d
) n x
)))
2743 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2745 ;; Transform Gegenbauer function to Jacobi P function
2746 ;; ultraspherical(n,v,x) = gamma(2*v+n)*gamma(v+1/2)/gamma(2*v)/gamma(v+n+1/2)
2747 ;; *jacobi_p(n,v-1/2,v-1/2,x)
2748 (defun ctpjac (x n v
)
2749 (let ((inv2 '((rat simp
) 1 2)))
2750 (mul (take '(%gamma
) (add v v n
))
2751 (inv (take '(%gamma
) (add v v
)))
2752 (take '(%gamma
) (add inv2 v
))
2753 (inv (take '(%gamma
) (add v inv2 n
)))
2754 (pjac x n
(sub v inv2
) (sub v inv2
)))))
2756 ;; Transform Chebyshev T function to Jacobi P function
2757 ;; chebyshev_t(n,x) = gamma(n+1)*sqrt(%pi)/gamma(n+1/2)*jacobi_p(n,-1/2,-1/2,x)
2759 (let ((inv2 '((rat simp
) 1 2)))
2760 (mul (take '(%gamma
) n
1)
2761 (power '$%pi inv2
) ; gamma(1/2)
2762 (inv (take '(%gamma
) (add inv2 n
)))
2763 (pjac x n
(mul -
1 inv2
) (mul -
1 inv2
)))))
2765 ;; Transform Chebyshev U function to Jacobi P function
2766 ;; chebyshev_u(n,x) = gamma(n+2)*sqrt(%pi)/2/gamma(3/2+n)*jacobi_p(n,1/2,1/2,x)
2768 (let ((inv2 '((rat simp
) 1 2)))
2769 (mul (take '(%gamma
) (add n
2))
2771 (power '$%pi inv2
) ; gamma(1/2)
2772 (inv (take '(%gamma
) (add inv2 n
1)))
2773 (pjac x n inv2 inv2
))))
2775 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2777 (defun pjactest (rest arg index1 index2 index3
)
2778 (cond ((maxima-integerp index1
)
2782 (list index1 index2 index3
)))
2784 (setq *hyp-return-noun-flag
* 'ind-should-be-an-integer-in-polys
))))
2786 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2788 ;;; Laplace transform of a single Bessel Y function.
2790 ;;; REST is the multiplier, ARG1 is the arg, and INDEX1 is the order of
2791 ;;; the Bessel Y function.
2792 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2794 (defun lt1yref (rest arg1 index1
)
2795 ;; If the index is an integer, use LT1Y. Otherwise, convert Bessel
2796 ;; Y to Bessel J and compute the transform of that. We do this
2797 ;; because converting Y to J for an integer index doesn't work so
2798 ;; well without taking limits.
2799 (cond ((maxima-integerp index1
)
2800 ;; Do not call lt1y but lty directly.
2801 ;; lt1y calls lt-ltp with the flag 'oney. lt-ltp checks this flag
2802 ;; and calls lty. So we can do it at this place and the algorithm is
2804 (lty rest arg1 index1
))
2805 (t (fractest2 rest arg1 index1 nil
'ytj
))))
2807 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2809 ;;; TRANSFORMATIONS TO CHANGE THE REPRESENTATION OF SPECIAL FUNCTIONS
2811 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2813 ;; erfc in terms of D, parabolic cylinder function
2815 ;; Tables of Integral Transforms
2818 ;; erfc(x) = (%pi*x)^(-1/2)*exp(-x^2/2)*W[-1/4,1/4](x^2)
2821 ;; D[v](z) = 2^(v/2+1/2)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
2825 ;; erfc(x) = %pi^(-1/2)*2^(1/4)*exp(-x^2/2)*D[-1](x*sqrt(2))
2828 (let ((inv2 '((rat simp
) 1 2)))
2829 (mul (power 2 inv2
) ; Should this be 2^(1/4)?
2830 (inv (power '$%pi inv2
))
2831 (power '$%e
(mul -
1 inv2 x x
))
2832 ;; At this time the Parabolic Cylinder D function is not implemented
2833 ;; as a simplifying function. We call nevertheless the simplifer
2834 ;; to simplify the arguments. When we implement the function
2835 ;; The symbol has to be changed to the noun form.
2836 (take '($parabolic_cylinder_d
) -
1 (mul (power 2 inv2
) x
)))))
2838 ;; Lommel S function in terms of Bessel J and Bessel Y.
2841 ;; S[u,v](z) = s[u,v](z) + {2^(u-1)*gamma((u-v+1)/2)*gamma((u+v+1)/2)}
2842 ;; * {sin[(u-v)*%pi/2]*bessel_j(v,z)
2843 ;; - cos[(u-v)*%pi/2]*bessel_y(v,z)
2845 (defun slommeltjandy (m n z
)
2846 (let ((arg (mul '((rat simp
) 1 2) '$%pi
(sub m n
))))
2847 (add (littleslommel m n z
)
2848 (mul (power 2 (sub m
1))
2849 (take '(%gamma
) (div (sub (add m
1) n
) 2))
2850 (take '(%gamma
) (div (add m n
1) 2))
2851 (sub (mul (take '(%sin
) arg
)
2852 (take '(%bessel_j
) n z
))
2853 (mul (take '(%cos
) arg
)
2854 (take '(%bessel_y
) n z
)))))))
2856 ;; Whittaker W function in terms of Whittaker M function
2860 ;; W[k,u](z) = gamma(-2*u)/gamma(1/2-u-k)*M[k,u](z)
2861 ;; + gamma(2*u)/gamma(1/2+u-k)*M[k,-u](z)
2863 (defun wtm (a i1 i2
)
2864 (add (mul (take '(%gamma
) (mul -
2 i2
))
2866 (inv (take '(%gamma
) (sub (sub '((rat simp
) 1 2) i2
) i1
))))
2867 (mul (take '(%gamma
) (add i2 i2
))
2868 (mwhit a i1
(mul -
1 i2
))
2869 (inv (take '(%gamma
) (sub (add '((rat simp
) 1 2) i2
) i1
))))))
2871 ;; Incomplete gamma function in terms of Whittaker W function
2873 ;; Tables of Integral Transforms, p. 387
2875 ;; gamma_incomplete(a,x) = x^((a-1)/2)*exp(-x/2)*W[(a-1)/2,a/2](x)
2877 (defun gammaincompletetw (a x
)
2878 (mul (power x
(div (sub a
1) 2))
2879 (power '$%e
(div x -
2))
2880 (wwhit x
(div (sub a
1) 2)(div a
2))))
2882 ;;; Incomplete Gamma function in terms of lower incomplete Gamma function
2884 ;;; Only for a=0 we use the general representation as a Whittaker W function:
2886 ;;; gamma_incomplete(a,x) = x^((a-1)/2)*exp(-x/2)*W[(a-1)/2,a/2](x)
2888 ;;; In all other cases we transform to a lower incomplete Gamma function:
2890 ;;; gamma_incomplete(a,x) = gamma(a)- gamma_incomplete_lower(a,x)
2892 ;;; The lower incomplete Gamma function will be further transformed to a Hypergeometric 1F1
2893 ;;; representation. With this change we get more simple and correct results for
2894 ;;; the Laplace transform of the Incomplete Gamma function.
2896 (defun gamma_incomplete-to-gamma-incomplete-lower (a x
)
2897 (if (or (eq ($sign a
) '$zero
)
2898 (and (integerp a
) (< a
0)))
2899 ;; The representation as a Whittaker W function for a=0 or a negative
2900 ;; integer (The gamma function is not defined for this case.)
2901 (mul (power x
(div (sub a
1) 2))
2902 (power '$%e
(div x -
2))
2903 (wwhit x
(div (sub a
1) 2) (div a
2)))
2904 ;; In all other cases the representation as a lower incomplete Gamma function
2905 (sub (take '(%gamma
) a
)
2906 (list '(%gamma_incomplete_lower simp
) a x
))))
2908 ;; Bessel Y in terms of Bessel J
2912 ;; bessel_y(v,z) = bessel_j(v,z)*cot(v*%pi)-bessel_j(-v,z)/sin(v*%pi)
2915 (sub (mul (take '(%bessel_j
) i a
)
2916 (take '(%cot
) (mul i
'$%pi
)))
2917 (mul (take '(%bessel_j
) (mul -
1 i
) a
)
2918 (inv (take '(%sin
) (mul i
'$%pi
))))))
2920 ;; Parabolic cylinder function in terms of Whittaker W function.
2922 ;; See Table of Integral Transforms, p.386:
2924 ;; D[v](z) = 2^(v/2+1/4)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
2927 (mul (power 2 (add (div i
2) '((rat simp
) 1 4)))
2928 (power a
'((rat simp
) -
1 2))
2929 (wwhit (mul a a
'((rat simp
) 1 2))
2930 (add (div i
2) '((rat simp
) 1 4))
2931 '((rat simp
) 1 4))))
2933 ;; Bateman's function in terms of Whittaker W function
2935 ;; See Table of Integral Transforms, p.386:
2937 ;; k[2*v](z) = 1/gamma(v+1)*W[v,1/2](2*z)
2939 (defun kbatemantw (v a
)
2940 (div (wwhit (add a a
) (div v
2) '((rat simp
) 1 2))
2941 (take '(%gamma
) (add (div v
2) 1))))
2943 ;; Bessel K in terms of Bessel I
2947 ;; bessel_k(v,z) = %pi/2*(bessel_i(-v,z)-bessel_i(v,z))/sin(v*%pi)
2952 (inv (take '(%sin
) (mul i
'$%pi
)))
2953 (sub (take '(%bessel_i
) (mul -
1 i
) a
)
2954 (take '(%bessel_i
) i a
))))
2956 ;; Express Hankel function in terms of Bessel J or Y function.
2960 ;; H[v,1](z) = %i*csc(v*%pi)*(exp(-v*%pi*%i)*bessel_j(v,z) - bessel_j(-v,z))
2963 ;; H[v,2](z) = %i*csc(v*%pi)*(bessel_j(-v,z) - exp(-v*%pi*%i)*bessel_j(v,z))
2965 ;; Both formula are not valid for v an integer.
2966 ;; For this case use the definitions of the Hankel functions:
2967 ;; H[v,1](z) = bessel_j(v,z) + %i* bessel_y(v,z)
2968 ;; H[v,2](z) = bessel_j(v,z) - %i* bessel_y(v,z)
2970 ;; All this can be implemented more simple.
2971 ;; We do not need the transformation to bessel_j for rational numbers,
2972 ;; because the correct transformation for bessel_y is already implemented.
2973 ;; It is enough to use the definitions for the Hankel functions.
2975 (defun htjory (v sort z
)
2976 ;; V is the order, SORT is the kind of Hankel function (1 or 2), Z
2978 (cond ((and (consp v
)
2980 (equal (caar v
) 'rat
))
2981 ;; If the order is a rational number of some sort,
2983 ;; (bessel_j(-v,z) - bessel_j(v,z)*exp(-v*%pi*%i))/(%i*sin(v*%pi*%i))
2984 (div (numjory v sort z
'j
)
2985 (mul '$%i
(take '(%sin
) (mul v
'$%pi
)))))
2987 ;; Transform hankel_1(v,z) to bessel_j(v,z)+%i*bessel_y(v,z)
2988 (add (take '(%bessel_j
) v z
)
2989 (mul '$%i
(take '(%bessel_y
) v z
))))
2991 ;; Transform hankel_2(v,z) to bessel_j(v,z)-%i*bessel_y(v,z)
2992 (sub (take '(%bessel_j
) v z
)
2993 (mul '$%i
(take '(%bessel_y
) v z
))))
2995 ;; We should never reach this point of code.
2996 ;; Problem: The user input for the symbol %h[v,sort](t) is not checked.
2997 ;; Therefore the user can generate this error as long as we do not cut
2998 ;; out the support for the Laplace transform of the symbol %h.
2999 (merror "htjory: Called with wrong sort of Hankel functions."))))
3001 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3003 ;;; Three helper functions only used by htjory
3005 ;; Bessel J or Y, depending on if FLG is 'J or not.
3006 (defun desjy (v z flg
)
3008 (take '(%bessel_j
) v z
))
3010 (take '(%bessel_y
) v z
))))
3012 (defun numjory (v sort z flg
)
3013 (cond ((equal sort
1)
3014 ;; bessel(-v, z) - exp(-v*%pi*%i)*bessel(v, z)
3016 ;; Where bessel is bessel_j if FLG is 'j. Otherwise, bessel
3019 ;; bessel_y(-v, z) - exp(-v*%pi*%i)*bessel_y(v, z)
3020 (sub (desjy (mul -
1 v
) z flg
)
3021 (mul (power '$%e
(mul -
1 v
'$%pi
'$%i
))
3024 ;; exp(-v*%pi*%i)*bessel(v,z) - bessel(-v,z), where bessel is
3025 ;; bessel_j or bessel_y, depending on if FLG is 'j or not.
3026 (sub (mul (power '$%e
(mul v
'$%pi
'$%i
))
3028 (desmjy (mul -
1 v
) z flg
)))))
3030 (defun desmjy (v z flg
)
3033 (take '(%bessel_j
) v z
))
3036 (mul -
1 (take '(%bessel_y
) v z
)))))
3038 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3040 ;;; TRANSFORM TO HYPERGEOMETRIC FUNCTION WITHOUT USING THE ROUTINE REF
3042 ;;; This functions are called in the routine lt-sf-log to get the
3043 ;;; representation in terms of a hypergeometric function.
3044 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3046 ;; The algorithm of the implemented Hermite function %he does not work for
3047 ;; the known Laplace transforms. For an even or odd integer order, we
3048 ;; can represent the Hermite function by the Hypergeometric function 1F1.
3049 ;; With this representations we get the expected Laplace transforms.
3050 (defun hermite-to-hypergeometric (order arg
)
3052 ((and (maxima-integerp order
)
3054 ;; Transform to 1F1 for order an even integer
3055 (mul (power 2 order
)
3056 (power '$%pi
(div 1 2))
3057 (inv (take '(%gamma
) (div (sub 1 order
) 2)))
3058 (list '(mqapply simp
)
3059 (list '($%f array simp
) 1 1)
3060 (list '(mlist) (div order -
2))
3061 (list '(mlist) '((rat simp
) 1 2))
3064 ((and (maxima-integerp order
)
3066 ;; Transform to 1F1 for order an odd integer
3069 (power '$%pi
'((rat simp
) 1 2))
3070 (inv (take '(%gamma
) (div order -
2)))
3071 (list '(mqapply simp
)
3072 (list '($%f simp array
) 1 1)
3073 (list '(mlist simp
) (div (sub 1 order
) 2))
3074 (list '(mlist simp
) '((rat simp
) 3 2))
3077 ;; The general case, transform to 2F0
3078 ;; For this case we have no Laplace transform.
3079 (mul (power (mul 2 arg
) order
)
3080 (list '(mqapply simp
)
3081 (list '($%f array simp
) 2 0)
3082 (list '(mlist simp
) (div order
2) (div (sub 1 order
) 2))
3083 (list '(mlist simp
))
3084 (div -
1 (mul arg arg
)))))))
3086 ;;; Hypergeometric representation of the Exponential Integral Si
3087 ;;; Si(z) = z*1F2([1/2],[3/2,3/2],-z^2/4)
3088 (defun expintegral_si-to-hypergeometric (arg)
3090 (list '(mqapply simp
)
3091 (list '($%f array simp
) 1 2)
3092 (list '(mlist simp
) '((rat simp
) 1 2))
3093 (list '(mlist simp
) '((rat simp
) 3 2) '((rat simp
) 3 2))
3094 (div (mul -
1 arg arg
) 4))))
3096 ;;; Hypergeometric representation of the Exponential Integral Shi
3097 ;;; Shi(z) = z*1F2([1/2],[3/2,3/2],z^2/4)
3098 (defun expintegral_shi-to-hypergeometric (arg)
3100 (list '(mqapply simp
)
3101 (list '($%f simp array
) 1 2)
3102 (list '(mlist simp
) '((rat simp
) 1 2))
3103 (list '(mlist simp
) '((rat simp
) 3 2) '((rat simp
) 3 2))
3104 (div (mul arg arg
) 4))))
3106 ;;; Hypergeometric representation of the Exponential Integral Ci
3107 ;;; Ci(z) = -z^2/4*2F3([1,1],[2,2,3/2],-z^2/4)+log(z)+%gamma
3108 (defun expintegral_ci-to-hypergeometric (arg)
3109 (add (mul (div (mul -
1 arg arg
) 4)
3110 (list '(mqapply simp
)
3111 (list '($%f simp array
) 2 3)
3112 (list '(mlist simp
) 1 1)
3113 (list '(mlist simp
) 2 2 '((rat simp
) 3 2))
3114 (div (mul -
1 arg arg
) 4)))
3118 ;;; Hypergeometric representation of the Exponential Integral Chi
3119 ;;; Chi(z) = z^2/4*2F3([1,1],[2,2,3/2],z^2/4)+log(z)+%gamma
3120 (defun expintegral_chi-to-hypergeometric (arg)
3121 (add (mul (div (mul arg arg
) 4)
3122 (list '(mqapply simp
)
3123 (list '($%f array simp
) 2 3)
3124 (list '(mlist simp
) 1 1)
3125 (list '(mlist simp
) 2 2 '((rat simp
) 3 2))
3126 (div (mul arg arg
) 4)))
3130 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3132 ;;; EXPERTS ON LAPLACE TRANSFORMS
3134 ;;; LT<foo> functions are various experts on Laplace transforms of the
3135 ;;; function <foo>. The expression being transformed is
3136 ;;; r*<foo>(args). The first arg of each expert is r, The remaining
3137 ;;; args are the arg(s) and/or parameters.
3138 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3140 ;; Laplace transform of one Bessel J
3141 (defun lt1j (rest arg index
)
3142 (lt-ltp 'onej rest arg index
))
3144 ;; Laplace transform of two Bessel J functions.
3145 ;; The argument of each must be the same, but the orders may be different.
3146 (defun lt2j (rest arg1 arg2 index1 index2
)
3147 (cond ((not (equal arg1 arg2
))
3148 (setq *hyp-return-noun-flag
* 'product-of-bessel-with-different-args
))
3150 ;; Call lt-ltp two transform and integrate two Bessel J functions.
3151 (lt-ltp 'twoj rest arg1
(list 'list index1 index2
)))))
3153 ;; Laplace transform of a square of a Bessel J function
3154 (defun lt1j^
2 (rest arg index
)
3155 (cond ((alike1 index
'((rat) -
1 2))
3156 ;; Special case: Laplace transform of bessel_j(v,arg)^2, v = -1/2.
3157 ;; For this case the algorithm for the product of two Bessel functions
3158 ;; does not work. Call the integrator with the hypergeometric
3159 ;; representation: 2/%pi/arg*1/2*(1+0F1([], [1/2], -arg*arg)).
3160 (sendexec (mul (div 2 '$%pi
)
3163 (add '((rat simp
) 1 2)
3164 (mul '((rat simp
) 1 2)
3165 (list '(mqapply simp
)
3166 (list '($%f simp array
) 1 0)
3167 (list '(mlist simp
) )
3168 (list '(mlist simp
) '((rat simp
) 1 2))
3169 (mul -
1 (mul arg arg
)))))))
3171 (lt-ltp 'twoj rest arg
(list 'list index index
)))))
3173 ;; Laplace transform of Incomplete Gamma function
3174 (defun lt1gamma-incomplete-lower (rest arg1 arg2
)
3175 (lt-ltp 'gamma-incomplete-lower rest arg2 arg1
))
3177 ;; Laplace transform of Whittaker M function
3178 (defun lt1m (r a i1 i2
)
3179 (lt-ltp 'onem r a
(list i1 i2
)))
3181 ;; Laplace transform of Jacobi function
3182 (defun lt1p (r a i1 i2
)
3183 (lt-ltp 'hyp-onep r a
(list i1 i2
)))
3185 ;; Laplace transform of Associated Legendre function of the second kind
3186 (defun lt1q (r a i1 i2
)
3187 (lt-ltp 'oneq r a
(list i1 i2
)))
3189 ;; Laplace transform of Erf function
3190 (defun lt1erf (rest arg
)
3191 (lt-ltp 'onerf rest arg nil
))
3193 ;; Laplace transform of Log function
3194 (defun lt1log (rest arg
)
3195 (lt-ltp 'onelog rest arg nil
))
3197 ;; Laplace transform of Complete elliptic integral of the first kind
3198 (defun lt1kelliptic (rest arg
)
3199 (lt-ltp 'onekelliptic rest arg nil
))
3201 ;; Laplace transform of Complete elliptic integral of the second kind
3202 (defun lt1e (rest arg
)
3203 (lt-ltp 'onee rest arg nil
))
3205 ;; Laplace transform of Struve H function
3206 (defun lt1hstruve (rest arg1 index1
)
3207 (lt-ltp 'hs rest arg1 index1
))
3209 ;; Laplace transform of Struve L function
3210 (defun lt1lstruve (rest arg1 index1
)
3211 (lt-ltp 'hl rest arg1 index1
))
3213 ;; Laplace transform of Lommel s function
3214 (defun lt1s (rest arg1 index1 index2
)
3215 (lt-ltp 's rest arg1
(list index1 index2
)))
3217 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3219 ;;; TRANSFORM TO HYPERGEOMETRIC FUNCTION AND DO THE INTEGRATION
3221 ;;; FLG = special function we're transforming
3222 ;;; REST = other stuff
3223 ;;; ARG = arg of special function
3224 ;;; INDEX = index of special function.
3226 ;;; So we're transforming REST*FLG(INDEX, ARG).
3227 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3229 (defun lt-ltp (flg rest arg index
)
3231 (when (and (listp index
)
3232 (eq (car index
) 'list
))
3233 (setq index
(list (cadr index
) (caddr index
))))
3237 (car (setq l1
(ref flg index arg
)))))
3239 ;; Convert the special function to a hypgergeometric
3240 ;; function. L1 is the special function converted to the
3241 ;; hypergeometric function. d*x^m*%e^a*x looks for that
3242 ;; factor in the expanded form.
3245 ;; We currently don't know how to handle this yet.
3246 ;; We add the return of a noun form.
3247 (setq *hyp-return-noun-flag
* 'other-ca-later
)))))
3249 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3250 ;; Dispatch function to convert the given function to a hypergeometric
3253 ;; The first arg is a symbol naming the function; the last is the
3254 ;; argument to the function. The second arg is the index (or list of
3255 ;; indices) to the function. Not used if the function doesn't have
3258 ;; The result is a list of 2 elements: The first element is a
3259 ;; multiplier; the second, the hypergeometric function itself.
3260 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3262 (defun ref (flg index arg
)
3264 (onej (j1tf index arg
))
3265 (twoj (j2tf (car index
) (cadr index
) arg
))
3266 (hs (hstf index arg
))
3267 (hl (lstf index arg
))
3268 (s (stf (car index
) (cadr index
) arg
))
3270 (onelog (logtf arg
))
3271 (onekelliptic (kelliptictf arg
)) ; elliptic_kc
3272 (onee (etf arg
)) ; elliptic_ec
3273 (onem (mtf (car index
) (cadr index
) arg
))
3274 (hyp-onep (ptf (car index
) (cadr index
) arg
))
3275 (oneq (qtf (car index
) (cadr index
) arg
))
3276 (gamma-incomplete-lower (gamma-incomplete-lower-tf index arg
))
3277 (onepjac (pjactf (car index
) (cadr index
) (caddr index
) arg
))
3281 ;; Transform %f to internal representation FPQ
3282 (list 1 (ref-fpq (rest (car index
)) (rest (cadr index
)) arg
)))))
3284 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3286 ;;; TRANSFORM FUNCTION IN TERMS OF HYPERGEOMETRIC FUNCTION
3288 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3290 ;; Create a hypergeometric form that we recognize. The function is
3291 ;; %f[n,m](p; q; arg). We represent this as a list of the form
3292 ;; (fpq (<length p> <length q>) <p> <q> <arg>)
3294 (defun ref-fpq (p q arg
)
3295 (list 'fpq
(list (length p
) (length q
))
3298 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3300 ;; Whittaker M function in terms of hypergeometric function
3304 ;; M[k,u](z) = exp(-z/2)*z^(1/2+u)*M(1/2+u-k,1+2*u,z)
3306 (defun mtf (i1 i2 arg
)
3307 (list (mul (power arg
(add i2
'((rat simp
) 1 2)))
3308 (power '$%e
(div arg -
2)))
3309 (ref-fpq (list (add '((rat simp
) 1 2) i2
(mul -
1 i1
)))
3310 (list (add i2 i2
1))
3313 ;; Jacobi P in terms of hypergeometric function
3317 ;; F(-n,a+1+b+n; a+1; x) = n!/poch(a+1,n)*jacobi_p(n,a,b,1-2*x)
3319 ;; jacobi_p(n,a,b,x) = poch(a+1,n)/n!*F(-n,a+1+b+n; a+1; (1-x)/2)
3320 ;; = gamma(a+n+1)/gamma(a+1)/n!*F(-n,a+1+b+n; a+1; (1-x)/2)
3322 ;; We have a problem:
3323 ;; We transform the argument x -> (1-x)/2. But an argument with a constant
3324 ;; part is not integrable with the implemented algorithm for the hypergeometric
3325 ;; function. It might be possible to get a result for an argument y=1-2*x
3326 ;; with x=t^-2. But for this case the routine lt-ltp fails. The routine
3327 ;; recognize a constant term in the argument, but does not take into account
3328 ;; that the constant term might vanish, when we transform to a hypergeometric
3330 ;; Because of this the Laplace transform for the following functions does not
3331 ;; work too: Legendre P, Chebyshev T, Chebyshev U, and Gegenbauer.
3333 (defun pjactf (n a b x
)
3334 (list (mul (take '(%gamma
) (add n a
1))
3335 (inv (take '(%gamma
) (add a
1)))
3336 (inv (take '(%gamma
) (add n
1))))
3337 (ref-fpq (list (mul -
1 n
) (add n a b
1))
3339 (sub '((rat simp
) 1 2) (div x
2)))))
3341 ;; ArcSin in terms of hypergeometric function
3345 ;; F(1/2,1/2; 3/2; z^2) = asin(z)/z
3347 ;; asin(z) = z*F(1/2,1/2; 3/2; z^2)
3350 (let ((inv2 '((rat simp
) 1 2)))
3352 (ref-fpq (list inv2 inv2
)
3353 (list '((rat simp
) 3 2))
3356 ;; ArcTan in terms of hypergeometric function
3360 ;; F(1/2,1; 3/2; -z^2) = atan(z)/z
3362 ;; atan(z) = z*F(1/2,1; 3/2; -z^2)
3366 (ref-fpq (list '((rat simp
) 1 2) 1)
3367 (list '((rat simp
) 3 2))
3370 ;; Associated Legendre function P in terms of hypergeometric function
3374 ;; assoc_legendre_p(v,u,z) = ((z+1)/(z-2))^(u/2)/gamma(1-u)*F(-v,v+1;1-u,(1-z)/2)
3376 ;; FIXME: What about the branch cut? 8.1.2 is for z not on the real
3377 ;; line with -1 < z < 1.
3380 (list (mul (inv (take '(%gamma
) (sub 1 m
)))
3381 (power (div (add z
1)
3384 (ref-fpq (list (mul -
1 n
) (add n
1))
3386 (sub '((rat simp
) 1 2) (div z
2)))))
3388 ;; Associated Legendre function Q in terms of hypergeometric function
3392 ;; assoc_legendre_q(v,u,z)
3393 ;; = exp(%i*u*%pi)*2^(-v-1)*sqrt(%pi) *
3394 ;; gamma(v+u+1)/gamma(v+3/2)*z^(-v-u-1)*(z^2-1)^(u/2) *
3395 ;; F(1+v/2+u/2, 1/2+v/2+u/2; v+3/2; 1/z^2)
3397 ;; FIXME: What about the branch cut?
3400 (list (mul (power '$%e
(mul m
'$%pi
'$%i
))
3401 (power '$%pi
'((rat simp
) 1 2))
3402 (take '(%gamma
) (add m n
1))
3403 (power 2 (sub -
1 n
))
3404 (inv (take '(%gamma
) (add n
'((rat simp
) 3 2))))
3405 (power z
(mul -
1 (add m n
1)))
3406 (power (sub (mul z z
) 1) (div m
2)))
3407 (ref-fpq (list (div (add m n
1) 2) (div (add m n
2) 2))
3408 (list (add n
'((rat simp
) 3 2)))
3411 ;; Incomplete lower Gamma in terms of hypergeometric function
3415 ;; M(a,a+1,-x) = a*x^(-a)*gamma_incomplete_lower(a,x)
3417 ;; gamma_incomplete_lower(a,x) = x^a/a*M(a,a+1,-x)
3419 (defun gamma-incomplete-lower-tf (a x
)
3420 (list (mul (inv a
) (power x a
))
3425 ;; Complete elliptic K in terms of hypergeometric function
3429 ;; K(k) = %pi/2*2F1(1/2,1/2; 1; k^2)
3431 (defun kelliptictf (k)
3432 (let ((inv2 '((rat simp
) 1 2)))
3433 (list (mul inv2
'$%pi
)
3434 (ref-fpq (list inv2 inv2
)
3438 ;; Complete elliptic E in terms of hypergeometric function
3442 ;; E(k) = %pi/2*2F1(-1/2,1/2;1;k^2)
3445 (let ((inv2 '((rat simp
) 1 2)))
3446 (list (mul inv2
'$%pi
)
3447 (ref-fpq (list (mul -
1 inv2
) inv2
)
3451 ;; erf in terms of hypgeometric function.
3455 ;; erf(z) = 2*z/sqrt(%pi)*M(1/2,3/2,-z^2)
3456 ;; = 2*z/sqrt(%pi)*exp(-z^2)*M(1,3/2,z^2)
3459 (list (mul 2 arg
(power '$%pi
'((rat simp
) -
1 2)))
3460 (ref-fpq (list '((rat simp
) 1 2))
3461 (list '((rat simp
) 3 2))
3464 ;; log in terms of hypergeometric function
3466 ;; We know from A&S 15.1.3 that
3468 ;; F(1,1;2;z) = -log(1-z)/z.
3470 ;; So log(z) = (z-1)*F(1,1;2;1-z)
3478 ;; Bessel J function expressed as a hypergeometric function.
3482 ;; bessel_j(v,z) = (z/2)^v*sum (-z^2/4)^k/k!/gamma(v+k+1)
3485 ;; = (z/2)^v/gamma(v+1)*sum 1/poch(v+1,k)*(-z^2/4)^k/k!
3487 ;; = (z/2)^v/gamma(v+1) * 0F1(; v+1; -z^2/4)
3490 (list (mul (inv (power 2 v
))
3492 (inv (take '(%gamma
) (add v
1))))
3495 (mul '((rat simp
) -
1 4) (power z
2)))))
3497 ;; Product of 2 Bessel J functions in terms of hypergeometric function
3499 ;; See Y. L. Luke, formula 39, page 216:
3501 ;; bessel_j(u,z)*bessel_j(v,z)
3502 ;; = (z/2)^(u+v)/gamma(u+1)/gamma(v+1) *
3503 ;; 2F3((u+v+1)/2, (u+v+2)/2; u+1, v+1, u+v+1; -z^2)
3505 (defun j2tf (n m arg
)
3506 (list (mul (inv (take '(%gamma
) (add n
1)))
3507 (inv (take '(%gamma
) (add m
1)))
3508 (inv (power 2 (add n m
)))
3509 (power arg
(add n m
)))
3510 (ref-fpq (list (add '((rat simp
) 1 2) (div n
2) (div m
2))
3511 (add 1 (div n
2) (div m
2)))
3512 (list (add 1 n
) (add 1 m
) (add 1 n m
))
3513 (mul -
1 (power arg
2)))))
3515 ;; Struve H function in terms of hypergeometric function.
3517 ;; A&S 12.1.2 gives the following series for the Struve H function:
3520 ;; H[v](z) = (z/2)^(v+1)*sum (-1)^k*(z/2)^(2*k)/gamma(k+3/2)/gamma(k+v+3/2)
3523 ;; We can write this in the form
3525 ;; H[v](z) = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2)
3528 ;; * sum n!/poch(3/2,n)/poch(v+3/2,n)*(-z^2/4)^n/n!
3531 ;; = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2) * 1F2(1;3/2,v+3/2;(-z^2/4))
3533 ;; See also A&S 12.1.21.
3536 (let ((d32 (div 3 2)))
3537 (list (mul (power (div z
2)(add v
1))
3538 (inv (take '(%gamma
) d32
))
3539 (inv (take '(%gamma
) (add v d32
))))
3541 (list d32
(add v d32
))
3542 (mul '((rat simp
) -
1 4) z z
)))))
3544 ;; Struve L function in terms of hypergeometric function
3548 ;; L[v](z) = -%i*exp(-v*%i*%pi/2)*H[v](%i*z)
3550 ;; This function computes exactly this way. (But why is %i written as
3551 ;; exp(%i*%pi/2) instead of just %i)
3553 ;; A&S 12.2.1 gives the series expansion as
3556 ;; L[v](z) = (z/2)^(v+1)*sum (z/2)^(2*k)/gamma(k+3/2)/gamma(k+v+3/2)
3559 ;; It's quite easy to derive
3561 ;; L[v](z) = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2) * 1F2(1;3/2,v+3/2;(z^2/4))
3564 (let ((d32 (div 3 2)))
3565 (list (mul (power (div z
2) (add v
1))
3566 (inv (take '(%gamma
) d32
))
3567 (inv (take '(%gamma
) (add v d32
))))
3569 (list d32
(add v d32
))
3570 (mul '((rat simp
) 1 4) z z
)))))
3572 ;; Lommel s function in terms of hypergeometric function
3574 ;; See Y. L. Luke, p 217, formula 1
3576 ;; s(u,v,z) = z^(u+1)/(u-v+1)/(u+v+1)*1F2(1; (u-v+3)/2, (u+v+3)/2; -z^2/4)
3579 (list (mul (power z
(add m
1))
3580 (inv (sub (add m
1) n
))
3583 (list (div (sub (add m
3) n
) 2)
3584 (div (add m n
3) 2))
3585 (mul '((rat simp
) -
1 4) z z
))))
3587 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3589 ;;; Algorithm 3: Laplace transform of a hypergeometric function
3591 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3593 (defun %$etest
(l l1
)
3595 (setq q
(cdras 'q l
))
3596 (cond ((equal q
1)(setq a
0)(go loop
)))
3597 (setq a
(cdras 'a l
))
3599 (return (substl (sub *par
* a
)
3601 (execf19 l
(cadr l1
))))))
3603 (defun execf19 (l1 l2
)
3605 (setq ans
(execargmatch (car (cddddr l2
))))
3606 (cond ((eq (car ans
) 'dionimo
)
3607 (return (dionarghyp l1 l2
(cadr ans
)))))
3608 ;; We specialized the pattern for the argument. Now the test fails
3609 ;; correctly and we have to add the return of a noun form.
3611 (setq *hyp-return-noun-flag
* 'next-for-other-args
))))
3613 ;; Executive for recognizing the sort of argument to the
3614 ;; hypergeometric function. We look to see if the arg is of the form
3615 ;; a*x^m + c. Return a list of 'dionimo (what does that mean?) and
3618 (defun execargmatch (arg)
3620 (cond ((setq l1
(m2-a*x^m
+c
($factor arg
) *var
*))
3621 (return (list 'dionimo l1
))))
3622 (cond ((setq l1
(m2-a*x^m
+c
($expand arg
) *var
*))
3623 (return (list 'dionimo l1
))))
3624 ;; The return value has to be a list.
3625 (return (list 'other-case-args-to-follow
))))
3627 ;; We have hypergeometric function whose arg looks like a*x^m+c. L1
3628 ;; matches the d*x^m... part, L2 is the hypergeometric function and
3629 ;; arg is the match for a*x^m+c.
3631 (defun dionarghyp (l1 l2 arg
)
3633 (setq a
(cdras 'a arg
)
3637 ((and (integerp m
) ; The power of the argument has to be
3638 (plusp m
) ; a positive integer.
3639 (equal 0 c
)) ; No const term to the argument.
3640 (return (f19cond a m l1 l2
)))
3643 (setq *hyp-return-noun-flag
* 'prop4-and-other-cases-to-follow
))))))
3645 (defun f19cond (a m l1 l2
)
3653 ;; At this point, we have the function d*x^s*%f[p,q](l1, l2, (a*t)^m).
3654 ;; Check to see if Formula 19, p 220 applies.
3655 (cond ((and (not (eq ($asksign
(sub (add p m -
1) q
)) '$pos
))
3656 (eq ($asksign
(add s
1)) '$pos
))
3658 (f19p220-simp (add s
1)
3663 (return (setq *hyp-return-noun-flag
*
3664 'failed-on-f19cond-multiply-the-other-cases-with-d
))))
3666 ;; Table of Laplace transforms, p 220, formula 19:
3668 ;; If m + k <= n + 1, and Re(s) > 0, the Laplace transform of
3670 ;; t^(s-1)*%f[m,n]([a1,...,am],[p1,...,pn],(c*t)^k)
3673 ;; gamma(s)/p^s*%f[m+k,n]([a1,...,am,s/k,(s+1)/k,...,(s+k-1)/k],[p1,...,pm],(k*c/p)^k)
3675 ;; with Re(p) > 0 if m + k <= n, Re(p+k*c*exp(2*%pi*%i*r/k)) > 0 for r
3676 ;; = 0, 1,...,k-1, if m + k = n + 1.
3678 ;; The args below are s, [a's], [p's], c^k, k.
3680 (defun f19p220-simp (s l1 l2 cf k
)
3681 (mul (take '(%gamma
) s
)
3682 (inv (power *par
* s
))
3683 (hgfsimp-exec (append l1
(addarglist s k
))
3687 (power (inv *par
*) k
)))))
3689 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3691 ;; Return a list of s/k, (s+1)/k, ..., (s+|k|-1)/k
3692 (defun addarglist (s k
)
3693 (let ((abs-k (abs k
))
3696 (push (div (add s n
) k
) res
))
3699 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3701 ;;; Pattern for the Laplace transform of the hypergeometric function
3703 ;; Match d*x^m*%e^(a*x). If we match, Q is the e^(a*x) part, A is a,
3704 ;; M is M, and D is d.
3705 (defun m2-d*x^m
*%e^a
*x
(expr var par
)
3708 ((coefftt) (d free2
,var
,par
))
3709 ((mexpt) (x alike1
,var
) (m free2
,var
,par
))
3713 ((coefftt) (a free2
,var
,par
))
3714 (x alike1
,var
))))))
3717 (defun m2-f+c
(expr var
)
3720 ((coeffpt) (f has
,var
))
3721 ((coeffpp) (c free
,var
)))))
3724 ;; The pattern was too general. We match also a*t^2+b*t. But that's not correct.
3725 (defun m2-a*x^m
+c
(expr var
)
3728 ((coefft) ; more special (not coeffpt)
3730 ((mexpt) (x alike1
,var
) (m free-not-zero-p
,var
)))
3731 ((coeffpp) (c free
,var
)))))
3733 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3735 ;;; Algorithm 4: SPECIAL HANDLING OF Bessel Y for an integer order
3737 ;;; This is called for one Bessel Y function, when the order is an integer.
3738 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3740 (defun lty (rest arg index
)
3742 (cond ((setq l
(m2-d*x^m
*%e^a
*x rest
*var
* *par
*))
3743 (return (execfy l arg index
))))
3744 (return (setq *hyp-return-noun-flag
* 'fail-in-lty
))))
3746 (defun execfy (l arg index
)
3748 (setq ans
(execargmatch arg
))
3749 (cond ((eq (car ans
) 'dionimo
)
3750 (return (dionarghyp-y l index
(cadr ans
)))))
3751 (return (setq *hyp-return-noun-flag
* 'fail-in-execfy
))))
3753 (defun dionarghyp-y (l index arg
)
3755 (setq a
(cdras 'a arg
)
3758 (cond ((and (zerp c
) (equal m
1))
3759 (let ((ans (f2p105v2cond a l index
)))
3760 (unless (symbolp ans
)
3762 (cond ((and (zerp c
) (equal m
'((rat simp
) 1 2)))
3763 (let ((ans (f50cond a l index
)))
3764 (unless (symbolp ans
)
3766 (return (setq *hyp-return-noun-flag
* 'fail-in-dionarghyp-y
))))
3768 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3770 ;;; Algorithm 4.1: Laplace transform of t^n*bessel_y(v,a*t)
3771 ;;; v is an integer and n>=v
3773 ;;; Table of Integral Transforms
3775 ;;; Volume 2, p 105, formula 2 is a formula for the Y-transform of
3777 ;;; f(x) = x^(u-3/2)*exp(-a*x)
3779 ;;; where the Y-transform is defined by
3781 ;;; integrate(f(x)*bessel_y(v,x*y)*sqrt(x*y), x, 0, inf)
3785 ;;; -2/%pi*gamma(u+v)*sqrt(y)*(y^2+a^2)^(-u/2)
3786 ;;; *assoc_legendre_q(u-1,-v,a/sqrt(y^2+a^2))
3788 ;;; with a > 0, Re u > |Re v|.
3790 ;;; In particular, with a slight change of notation, we have
3792 ;;; integrate(x^(u-1)*exp(-p*x)*bessel_y(v,a*x)*sqrt(a), x, 0, inf)
3794 ;;; which is the Laplace transform of x^(u-1/2)*bessel_y(v,x).
3796 ;;; Thus, the Laplace transform is
3798 ;;; -2/%pi*gamma(u+v)*sqrt(a)*(a^2+p^2)^(-u/2)
3799 ;;; *assoc_legendre_q(u-1,-v,p/sqrt(a^2+p^2))
3800 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3802 (defun f2p105v2cond (a l index
)
3804 (setq d
(cdras 'd l
) ; contains constant part of integrand
3807 (cond ((eq ($asksign
($realpart
(sub m index
))) '$pos
)
3808 (return (mul d
(f2p105v2cond-simp m index a
)))))
3809 (return (setq *hyp-return-noun-flag
* 'fail-in-f2p105v2cond
))))
3811 (defun f2p105v2cond-simp (m v a
)
3814 (take '(%gamma
) (add m v
))
3815 (power (add (mul a a
) (mul *par
* *par
*))
3816 (mul -
1 '((rat simp
) 1 2) m
))
3817 ;; Call Associated Legendre Q function, which simplifies accordingly.
3818 ;; We have to do a Maxima function call, because $assoc_legendre_q is
3819 ;; not in Maxima core and has to be autoloaded.
3820 (mfuncall '$assoc_legendre_q
3824 (power (add (mul a a
) (mul *par
* *par
*))
3825 '((rat simp
) -
1 2))))))
3827 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3829 ;;; Algorithm 4.2: Laplace transform of t^n*bessel_y(v, a*sqrt(t))
3831 ;;; Table of Integral Transforms
3833 ;;; p. 188, formula 50:
3835 ;;; t^(u-1/2)*bessel_y(2*v,2*sqrt(a)*sqrt(t))
3836 ;;; -> a^(-1/2)*p^(-u)*exp(-a/2/p)
3837 ;;; * [tan((u-v)*%pi)*gamma(u+v+1/2)/gamma(2*v+1)*M[u,v](a/p)
3838 ;;; -sec((u-v)*%pi)*W[u,v](a/p)]
3839 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3841 (defun f50cond (a l v
)
3843 (setq d
(cdras 'd l
)
3844 m
(add (cdras 'm l
) '((rat simp
) 1 2))
3847 ((and (eq ($asksign
($realpart
(add m v
'((rat simp
) 1 2)))) '$pos
)
3848 (eq ($asksign
($realpart
(sub (add m
'((rat simp
) 1 2)) v
))) '$pos
)
3849 (not (maxima-integerp (mul (sub (add m m
) (add v v
1))
3850 '((rat simp
) 1 2)))))
3851 (setq a
(mul a a
'((rat simp
) 1 4)))
3852 (return (f50p188-simp d m v a
))))
3853 (return (setq *hyp-return-noun-flag
* 'fail-in-f50cond
))))
3855 (defun f50p188-simp (d u v a
)
3857 (power a
'((rat simp
) -
1 2))
3858 (power *par
* (mul -
1 u
))
3859 (power '$%e
(div a
(mul -
2 *par
*)))
3860 (sub (mul (take '(%tan
) (mul '$%pi
(sub u v
)))
3861 (take '(%gamma
) (add u v
'((rat simp
) 1 2)))
3862 (inv (take '(%gamma
) (add v v
1)))
3863 (mwhit (div a
*par
*) u v
))
3864 (mul (take '(%sec
) (mul '$%pi
(sub u v
)))
3865 (wwhit (div a
*par
*) u v
)))))
3867 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;