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 $expintrep
))
30 (defmvar $prefer_d nil
31 "When NIL express a parabolic cylinder function as hypergeometric function.")
33 ;; The properties NOUN and VERB give correct linear display
34 (defprop $specint %specint verb
)
35 (defprop %specint $specint noun
)
37 (defvar *hyp-return-noun-form-p
* t
38 "Return noun form instead of internal Lisp symbols for integrals
39 that we don't support.")
41 (defvar *hyp-return-noun-flag
* nil
42 "Flag to signal that we have no result and to return a noun.")
44 (defvar *debug-hypgeo
* nil
45 "Print debug information if enabled.")
47 ;; The variables *var* and *par* are global to this file only.
48 ;; They are initialized in the routine defexec. The values are never changed.
49 ;; These globals are introduced to avoid passing the values of *par* and *var*
50 ;; through all functions of this file.
53 "Variable of integration of Laplace transform.")
55 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
57 ;;; Helper function for this file
59 (defun substl (p1 p2 p3
)
61 (t (maxima-substitute p1 p2 p3
))))
63 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
65 ;;; Test functions for pattern match, which use the globals var and *par*
71 (cond ((equal m
0) nil
)
74 ;;; Test functions which do not depend on globals
76 ;; Test if EXP is 1 or %e.
84 (defun free-not-zero-p (expr x
)
85 (and (not (equal expr
0)) (free expr x
)))
87 (defun free2 (expr x y
)
88 (and (free expr x
) (free expr y
)))
90 (defun has-not-alike1-p (expr x
)
91 (and (not (alike1 expr x
)) (has expr x
)))
93 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
95 ;;; Some shortcuts for special functions.
97 ;; Lommel's little s[u,v](z) function.
98 (defun littleslommel (m n z
)
99 (list '(mqapply simp
) (list '($%s simp array
) m n
) z
))
101 ;; Whittaker's M function
102 (defun mwhit (a i1 i2
)
103 (list '(mqapply simp
) (list '($%m simp array
) i1 i2
) a
))
105 ;; Whittaker's W function
106 (defun wwhit (a i1 i2
)
107 (list '(mqapply simp
) (list '($%w simp array
) i1 i2
) a
))
110 (defun pjac (x n a b
)
111 (list '(mqapply simp
) (list '($%p simp array
) n a b
) x
))
113 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
115 ;;; Two general pattern for the routine lt-sf-log.
117 ;; Recognize c*u^v + a and a=0.
118 (defun m2-arbpow1 (expr var
)
122 (c free
,var
) ; more special to ensure that c is constant
123 ((mexpt) (u has
,var
) (v free
,var
)))
124 ((coeffpp) (a zerp
)))))
126 ;; Recognize c*u^v*(a+b*u)^w+d and d=0. This is a generalization of arbpow1.
127 (defun m2-arbpow2 (expr var
)
131 ((coefftt) (c free
,var
))
132 ((mexpt) (u equal
,var
) (v free
,var
))
134 ((mplus) (a free
,var
) ((coefft) (b free
,var
) (u equal
,var
)))
135 (w free-not-zero-p
,var
)))
136 ((coeffpp) (d zerp
)))))
138 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
140 ;;; The pattern to match special functions in the routine lt-sf-log.
143 (defun m2-asin (expr var
)
146 ((coeffpt) (u nonzerp
)
147 ((%asin
) (w has
,var
)))
148 ((coeffpp) (a equal
0)))))
151 (defun m2-atan (expr var
)
154 ((coeffpt) (u nonzerp
)
155 ((%atan
) (w has
,var
)))
156 ((coeffpp) (a equal
0)))))
158 ;; Recognize bessel_j(v,w)
159 (defun m2-onej (expr var
)
164 ((%bessel_j
) (v free
,var
) (w has
,var
)))
165 ((coeffpp) (a zerp
)))))
167 ;; Recognize bessel_j(v1,w1)*bessel_j(v2,w2)
168 (defun m2-twoj (expr var
)
173 ((%bessel_j
) (v1 free
,var
) (w1 has
,var
))
174 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
175 ((coeffpp) (a zerp
)))))
177 ;; Recognize bessel_y(v1,w1)*bessel_y(v2,w2)
178 (defun m2-twoy (expr var
)
183 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
184 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
185 ((coeffpp) (a zerp
)))))
187 ;; Recognize bessel_k(v1,w1)*bessel_k(v2,w2)
188 (defun m2-twok (expr var
)
193 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
194 ((%bessel_k
) (v2 free
,var
) (w2 has
,var
)))
195 ((coeffpp) (a zerp
)))))
197 ;; Recognize bessel_k(v1,w1)*bessel_y(v2,w2)
198 (defun m2-onekoney (expr var
)
203 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
204 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
205 ((coeffpp) (a zerp
)))))
207 ;; Recognize bessel_j(v,w)^2
208 (defun m2-onej^
2 (expr var
)
214 ((%bessel_j
) (v free
,var
) (w has
,var
))
216 ((coeffpp) (a zerp
)))))
218 ;; Recognize bessel_y(v,w)^2
219 (defun m2-oney^
2 (expr var
)
225 ((%bessel_y
) (v free
,var
) (w has
,var
))
227 ((coeffpp) (a zerp
)))))
229 ;; Recognize bessel_k(v,w)^2
230 (defun m2-onek^
2 (expr var
)
236 ((%bessel_k
) (v free
,var
) (w has
,var
))
238 ((coeffpp) (a zerp
)))))
240 ;; Recognize bessel_i(v,w)
241 (defun m2-onei (expr var
)
246 ((%bessel_i
) (v free
,var
) (w has
,var
)))
247 ((coeffpp) (a zerp
)))))
249 ;; Recognize bessel_i(v1,w1)*bessel_i(v2,w2)
250 (defun m2-twoi (expr var
)
255 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
256 ((%bessel_i
) (v2 free
,var
) (w2 has
,var
)))
257 ((coeffpp) (a zerp
)))))
259 ;; Recognize hankel_1(v1,w1)*hankel_1(v2,w2), product of 2 Hankel 1 functions.
260 (defun m2-two-hankel_1 (expr var
)
265 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
266 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
267 ((coeffpp) (a zerp
)))))
269 ;; Recognize hankel_2(v1,w1)*hankel_2(v2,w2), product of 2 Hankel 2 functions.
270 (defun m2-two-hankel_2 (expr var
)
275 ((%hankel_2
) (v1 free
,var
) (w1 has
,var
))
276 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
277 ((coeffpp) (a zerp
)))))
279 ;; Recognize hankel_1(v1,w1)*hankel_2(v2,w2), product of 2 Hankel functions.
280 (defun m2-hankel_1*hankel_2
(expr var
)
285 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
286 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
287 ((coeffpp) (a zerp
)))))
289 ;; Recognize bessel_y(v1,w1)*bessel_j(v2,w2)
290 (defun m2-oneyonej (expr var
)
295 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
296 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
297 ((coeffpp) (a zerp
)))))
299 ;; Recognize bessel_k(v1,w1)*bessel_j(v2,w2)
300 (defun m2-onekonej (expr var
)
305 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
306 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
307 ((coeffpp) (a zerp
)))))
309 ;; Recognize bessel_y(v1,w1)*hankel_1(v2,w2)
310 (defun m2-bessel_y*hankel_1
(expr var
)
315 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
316 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
317 ((coeffpp) (a zerp
)))))
319 ;; Recognize bessel_y(v1,w1)*hankel_2(v2,w2)
320 (defun m2-bessel_y*hankel_2
(expr var
)
325 ((%bessel_y
) (v1 free
,var
) (w1 has
,var
))
326 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
327 ((coeffpp) (a zerp
)))))
329 ;; Recognize bessel_k(v1,w1)*hankel_1(v2,w2)
330 (defun m2-bessel_k*hankel_1
(expr var
)
335 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
336 ((%hankel_1
) (v1 free
,var
) (w2 has
,var
)))
337 ((coeffpp) (a zerp
)))))
339 ;; Recognize bessel_k(v1,w1)*hankel_2(v2,w2)
340 (defun m2-bessel_k*hankel_2
(expr var
)
345 ((%bessel_k
) (v1 free
,var
) (w1 has
,var
))
346 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
347 ((coeffpp) (a zerp
)))))
349 ;; Recognize bessel_i(v1,w1)*bessel_j(v2,w2)
350 (defun m2-oneionej (expr var
)
355 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
356 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
357 ((coeffpp) (a zerp
)))))
359 ;; Recognize bessel_i(v1,w1)*hankel_1(v2,w2)
360 (defun m2-bessel_i*hankel_1
(expr var
)
365 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
366 ((%hankel_1
) (v2 free
,var
) (w2 has
,var
)))
367 ((coeffpp) (a zerp
)))))
369 ;; Recognize bessel_i(v1,w1)*hankel_2(v2,w2)
370 (defun m2-bessel_i*hankel_2
(expr var
)
375 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
376 ((%hankel_2
) (v2 free
,var
) (w2 has
,var
)))
377 ((coeffpp) (a zerp
)))))
379 ;; Recognize hankel_1(v1,w1)*bessel_j(v2,w2)
380 (defun m2-hankel_1*bessel_j
(expr var
)
385 ((%hankel_1
) (v1 free
,var
) (w1 has
,var
))
386 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
387 ((coeffpp) (a zerp
)))))
389 ;; Recognize hankel_2(v1,w1)*bessel_j(v2,w2)
390 (defun m2-hankel_2*bessel_j
(expr var
)
395 ((%hankel_2
) (v1 free
,var
) (w1 has
,var
))
396 ((%bessel_j
) (v2 free
,var
) (w2 has
,var
)))
397 ((coeffpp) (a zerp
)))))
399 ;; Recognize bessel_i(v1,w1)*bessel_y(v2,w2)
400 (defun m2-oneioney (expr var
)
405 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
406 ((%bessel_y
) (v2 free
,var
) (w2 has
,var
)))
407 ((coeffpp) (a zerp
)))))
409 ;; Recognize bessel_i(v1,w1)*bessel_k(v2,w2)
410 (defun m2-oneionek (expr var
)
415 ((%bessel_i
) (v1 free
,var
) (w1 has
,var
))
416 ((%bessel_k
) (v2 free
,var
) (w2 has
,var
)))
417 ((coeffpp) (a zerp
)))))
419 ;; Recognize bessel_i(v,w)^2
420 (defun m2-onei^
2 (expr var
)
426 ((%bessel_i
) (v free
,var
) (w has
,var
))
428 ((coeffpp) (a zerp
)))))
430 ;; Recognize hankel_1(v,w)^2
431 (defun m2-hankel_1^
2 (expr var
)
437 ((%hankel_1
) (v free
,var
) (w has
,var
))
439 ((coeffpp) (a zerp
)))))
441 ;; Recognize hankel_2(v,w)^2
442 (defun m2-hankel_2^
2 (expr var
)
448 ((%hankel_2
) (v free
,var
) (w has
,var
))
450 ((coeffpp) (a zerp
)))))
452 ;; Recognize bessel_y(v,w)
453 (defun m2-oney (expr var
)
458 ((%bessel_y
) (v free
,var
) (w has
,var
)))
459 ((coeffpp) (a zerp
)))))
461 ;; Recognize bessel_k(v,w)
462 (defun m2-onek (expr var
)
467 ((%bessel_k
) (v free
,var
) (w has
,var
)))
468 ((coeffpp) (a zerp
)))))
470 ;; Recognize hankel_1(v,w)
471 (defun m2-hankel_1 (expr var
)
476 ((%hankel_1
) (v free
,var
) (w has
,var
)))
477 ((coeffpp) (a zerp
)))))
479 ;; Recognize hankel_2(v,w)
480 (defun m2-hankel_2 (expr var
)
485 ((%hankel_2
) (v free
,var
) (w has
,var
)))
486 ((coeffpp) (a zerp
)))))
489 (defun m2-onelog (expr var
)
494 ((%log
) (w has
,var
)))
495 ((coeffpp) (a zerp
)))))
498 (defun m2-onerf (expr var
)
503 ((%erf
) (w has
,var
)))
504 ((coeffpp) (a zerp
)))))
507 (defun m2-onerfc (expr var
)
512 ((%erfc
) (w has
,var
)))
513 ((coeffpp) (a zerp
)))))
515 ;; Recognize fresnel_s(w)
516 (defun m2-onefresnel_s (expr var
)
521 ((%fresnel_s
) (w has
,var
)))
522 ((coeffpp) (a zerp
)))))
524 ;; Recognize fresnel_c(w)
525 (defun m2-onefresnel_c (expr var
)
530 ((%fresnel_c
) (w has
,var
)))
531 ((coeffpp) (a zerp
)))))
533 ;; Recognize expintegral_e(v,w)
534 (defun m2-oneexpintegral_e (expr var
)
539 ((%expintegral_e
) (v free
,var
) (w has
,var
)))
540 ((coeffpp) (a zerp
)))))
542 ;; Recognize expintegral_ei(w)
543 (defun m2-oneexpintegral_ei (expr var
)
548 ((%expintegral_ei
) (w has
,var
)))
549 ((coeffpp) (a zerp
)))))
551 ;; Recognize expintegral_e1(w)
552 (defun m2-oneexpintegral_e1 (expr var
)
557 ((%expintegral_e1
) (w has
,var
)))
558 ((coeffpp) (a zerp
)))))
560 ;; Recognize expintegral_si(w)
561 (defun m2-oneexpintegral_si (expr var
)
566 ((%expintegral_si
) (w has
,var
)))
567 ((coeffpp) (a zerp
)))))
569 ;; Recognize expintegral_shi(w)
570 (defun m2-oneexpintegral_shi (expr var
)
575 ((%expintegral_shi
) (w has
,var
)))
576 ((coeffpp) (a zerp
)))))
578 ;; Recognize expintegral_ci(w)
579 (defun m2-oneexpintegral_ci (expr var
)
584 ((%expintegral_ci
) (w has
,var
)))
585 ((coeffpp) (a zerp
)))))
587 ;; Recognize expintegral_chi(w)
588 (defun m2-oneexpintegral_chi (expr var
)
593 ((%expintegral_chi
) (w has
,var
)))
594 ((coeffpp) (a zerp
)))))
596 ;; Recognize kelliptic(w), (new would be elliptic_kc)
597 (defun m2-onekelliptic (expr var
)
602 (($kelliptic
) (w has
,var
)))
603 ((coeffpp) (a zerp
)))))
605 ;; Recognize elliptic_kc
606 (defun m2-elliptic_kc (expr var
)
611 ((%elliptic_kc
) (w has
,var
)))
612 ((coeffpp) (a zerp
)))))
614 ;; Recognize %e(w), (new would be elliptic_ec)
615 (defun m2-onee (expr var
)
620 (($%e
) (w has
,var
)))
621 ((coeffpp) (a zerp
)))))
623 ;; Recognize elliptic_ec
624 (defun m2-elliptic_ec (expr var
)
629 ((%elliptic_ec
) (w has
,var
)))
630 ((coeffpp) (a zerp
)))))
632 ;; Recognize gamma_incomplete(w1, w2), Incomplete Gamma function
633 (defun m2-onegammaincomplete (expr var
)
638 ((%gamma_incomplete
) (w1 free
,var
) (w2 has
,var
)))
639 ((coeffpp) (a zerp
)))))
641 ;; Recognize gamma_incomplete_lower(w1,w2), gamma(a)-gamma_incomplete(w1,w2)
642 (defun m2-onegamma-incomplete-lower (expr var
)
647 ((%gamma_incomplete_lower
) (w1 free
,var
) (w2 has
,var
)))
648 ((coeffpp) (a zerp
)))))
650 ;; Recognize Struve H function.
651 (defun m2-struve_h (expr var
)
656 ((%struve_h
) (v free
,var
) (w has
,var
)))
657 ((coeffpp)(a zerp
)))))
659 ;; Recognize Struve L function.
660 (defun m2-struve_l (expr var
)
665 ((%struve_l
) (v free
,var
) (w has
,var
)))
666 ((coeffpp) (a zerp
)))))
668 ;; Recognize Lommel s[v1,v2](w) function.
669 (defun m2-ones (expr var
)
675 (($%s array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
676 ((coeffpp) (a zerp
)))))
678 ;; Recognize S[v1,v2](w), Lommel function
679 (defun m2-oneslommel (expr var
)
685 (($slommel array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
686 ((coeffpp) (a zerp
)))))
688 ;; Recognize parabolic_cylinder_d function
689 (defun m2-parabolic_cylinder_d (expr var
)
694 (($parabolic_cylinder_d
) (v free
,var
) (w has
,var
)))
695 ((coeffpp) (a zerp
)))))
697 ;; Recognize kbatmann(v,w), Batemann function
698 (defun m2-onekbateman (expr var
)
704 (($kbateman array
) (v free
,var
)) (w has
,var
)))
705 ((coeffpp) (a zerp
)))))
707 ;; Recognize %l[v1,v2](w), Generalized Laguerre function
708 (defun m2-onel (expr var
)
714 (($%l array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
715 ((coeffpp) (a zerp
)))))
717 ;; Recognize gen_laguerre(v1,v2,w), Generalized Laguerre function
718 (defun m2-one-gen-laguerre (expr var
)
723 ((%gen_laguerre
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
724 ((coeffpp) (a zerp
)))))
726 ;; Recognize laguerre(v1,w), Laguerre function
727 (defun m2-one-laguerre (expr var
)
732 ((%laguerre
) (v1 free
,var
) (w has
,var
)))
733 ((coeffpp) (a zerp
)))))
735 ;; Recognize %c[v1,v2](w), Gegenbauer function
736 (defun m2-onec (expr var
)
742 (($%c array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
743 ((coeffpp) (a zerp
)))))
745 ;; Recognize %t[v1](w), Chebyshev function of the first kind
746 (defun m2-onet (expr var
)
751 ((mqapply) (($%t array
) (v1 free
,var
)) (w has
,var
)))
752 ((coeffpp) (a zerp
)))))
754 ;; Recognize %u[v1](w), Chebyshev function of the second kind
755 (defun m2-oneu (expr var
)
760 ((mqapply) (($%u array
) (v1 free
,var
)) (w has
,var
)))
761 ((coeffpp) (a zerp
)))))
763 ;; Recognize %p[v1,v2,v3](w), Jacobi function
764 (defun m2-onepjac (expr var
)
771 (v1 free
,var
) (v2 free
,var
) (v3 free
,var
)) (w has
,var
)))
772 ((coeffpp) (a zerp
)))))
774 ;; Recognize jacobi_p function
775 (defun m2-jacobi_p (expr var
)
781 (v1 free
,var
) (v2 free
,var
) (v3 free
,var
) (w has
,var
)))
782 ((coeffpp) (a zerp
)))))
784 ;; Recognize %p[v1,v2](w), Associated Legendre P function
785 (defun m2-hyp-onep (expr var
)
791 (($%p array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
792 ((coeffpp) (a zerp
)))))
794 ;; Recognize assoc_legendre_p function
795 (defun m2-assoc_legendre_p (expr var
)
800 (($assoc_legendre_p
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
801 ((coeffpp) (a zerp
)))))
803 ;; Recognize %p[v1](w), Legendre P function
804 (defun m2-onep0 (expr var
)
809 ((mqapply)(($%p array
) (v1 free
,var
)) (w has
,var
)))
810 ((coeffpp) (a zerp
)))))
812 ;; Recognize %p[v1](w), Legendre P function
813 (defun m2-legendre_p (expr var
)
818 (($legendre_p
) (v free
,var
)) (w has
,var
))
819 ((coeffpp) (a zerp
)))))
821 ;; Recognize hermite(v1,w), Hermite function
822 (defun m2-one-hermite (expr var
)
827 ((%hermite
) (v1 free
,var
) (w has
,var
)))
828 ((coeffpp) (a zerp
)))))
830 ;; Recognize %q[v1,v2](w), Associated Legendre function of the second kind
831 (defun m2-oneq (expr var
)
837 (($%q array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
838 ((coeffpp) (a zerp
)))))
840 ;; Recognize assoc_legendre_q function
841 (defun m2-assoc_legendre_q (expr var
)
846 (($assoc_legendre_q
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
847 ((coeffpp) (a zerp
)))))
849 ;; Recognize %w[v1,v2](w), Whittaker W function.
850 (defun m2-onew (expr var
)
856 (($%w array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
857 ((coeffpp) (a zerp
)))))
859 ;; Recognize whittaker_w function.
860 (defun m2-whittaker_w (expr var
)
865 (($whittaker_w
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
866 ((coeffpp) (a zerp
)))))
868 ;; Recognize %m[v1,v2](w), Whittaker M function
869 (defun m2-onem (expr var
)
875 (($%m array
) (v1 free
,var
) (v2 free
,var
)) (w has
,var
)))
876 ((coeffpp) (a zerp
)))))
878 ;; Recognize whittaker_m function.
879 (defun m2-whittaker_m (expr var
)
884 (($whittaker_m
) (v1 free
,var
) (v2 free
,var
) (w has
,var
)))
885 ((coeffpp) (a zerp
)))))
887 ;; Recognize %f[v1,v2](w1,w2,w3), Hypergeometric function
888 (defun m2-onef (expr var
)
894 (($%f array
) (v1 free
,var
) (v2 free
,var
))
898 ((coeffpp) (a zerp
)))))
900 ;; Recognize hypergeometric function
901 (defun m2-hypergeometric (expr var
)
906 ((%hypergeometric
) (w1 free
,var
) (w2 free
,var
) (w3 has
,var
)))
907 ((coeffpp) (a zerp
)))))
909 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
911 ;;; Pattern for the routine hypgeo-exec.
912 ;;; RECOGNIZES L.T.E. "U*%E^(A*X+E*F(X)-P*X+C)+D".
914 (defun m2-ltep (expr var par
)
922 ((coeffpt) (a free2
,var
,par
) (x alike1
,var
))
923 ((coeffpt) (e free2
,var
,par
) (f has
,var
))
924 ((mtimes) -
1 (p alike1
,par
) (x alike1
,var
))
925 ((coeffpp) (c free2
,var
,par
)))))
926 ((coeffpp) (d equal
0)))))
928 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
930 ;;; Pattern for the routine defexec.
931 ;;; This is trying to match EXP to u*%e^(a*x+e*f+c)+d
932 ;;; where a, c, and e are free of x, f is free of p, and d is 0.
934 (defun m2-defltep (expr var
)
942 ((coeffpt) (a free
,var
) (x alike1
,var
))
943 ((coeffpt) (e free
,var
) (f has-not-alike1-p
,var
))
944 ((coeffpp) (c free
,var
)))))
945 ((coeffpp) (d equal
0)))))
947 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
949 ;;; $specint is the Maxima User function
951 (defmfun $specint
(expr var
)
952 (prog ($radexpand
*checkcoefsignlist
*)
953 (setq $radexpand
'$all
)
954 (return (defintegrate expr var
))))
956 (defun defintegrate (expr var
)
957 ;; This used to have $exponentialize enabled for everything, but I
958 ;; don't think we should do that. If various routines want
959 ;; $exponentialize, let them set it themselves. So, for here, we
960 ;; want to expand the form with $exponentialize to convert trig and
961 ;; hyperbolic functions to exponential functions that we can handle.
962 (let ((form (let (($exponentialize t
))
963 ($factor
(resimplify expr
))))) ; At first factor the integrand.
965 ;; Because we call defintegrate recursively, we add code to end the
969 (cond ((and (numberp form
) (zerop form
)) (return-from defintegrate
0))
970 (t (return-from defintegrate
(list '(%specint simp
) form var
)))))
972 ;; We try to find a constant denominator. This is necessary to get results
973 ;; for integrands like u(t)/(a+b+c+...).
975 (let ((den ($denom form
)))
976 (when (and (not (equal 1 den
)) ($freeof var den
))
977 (return-from defintegrate
978 (div (defintegrate (mul den form
) var
) den
))))
980 ;; We search for a sum of Exponential functions which we can integrate.
981 ;; This code finds result for Trigonometric or Hyperbolic functions with
982 ;; a factor t^-1 or t^-2 e.g. t^-1*sin(a*t).
984 (let* ((l (m2-defltep form var
))
985 (s (mul -
1 (cdras 'a l
)))
986 (u ($expand
(cdras 'u l
)))
989 ((setq l1
(m2-sum-with-exp-case1 u var
))
990 ;; c * t^-1 * (%e^(-a*t) - %e^(-b*t)) + d
991 (let ((c (cdras 'c l1
))
992 (a (mul -
1 (cdras 'a l1
)))
993 (b (mul -
1 (cdras 'b l1
)))
995 (add (mul c
(take '(%log
) (div (add s b
) (add s a
))))
996 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
998 ((setq l1
(m2-sum-with-exp-case2 u var
))
999 ;; c * t^(-3/2) * (%e^(-a*t) - %e^(-b*t)) + d
1000 (let ((c (cdras 'c l1
))
1001 (a (mul -
1 (cdras 'a l1
)))
1002 (b (mul -
1 (cdras 'b l1
)))
1005 (power '$%pi
'((rat simp
) 1 2))
1006 (sub (power (add s b
) '((rat simp
) 1 2))
1007 (power (add s a
) '((rat simp
) 1 2))))
1008 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1010 ((setq l1
(m2-sum-with-exp-case3 u var
))
1011 ;; c * t^-2 * (1 - 2 * %e^(-a*t) + %e^(2*a*t)) + d
1012 (let ((c (cdras 'c l1
))
1013 (a (div (cdras 'a l1
) -
2))
1016 (add (mul (add s a a
) (take '(%log
) (add s a a
)))
1017 (mul s
(take '(%log
) s
))
1018 (mul -
2 (add s a
) (take '(%log
) (add s a
)))))
1019 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1021 ((setq l1
(m2-sum-with-exp-case4 u var
))
1022 ;; c * t^-1 * (1 - 2 * %e^(-a*t) + %e^(2*a*t)) + d
1023 (let ((c (cdras 'c l1
))
1024 (a (div (cdras 'a l1
) (mul 4 '$%i
)))
1030 (mul (sub s
(mul 2 '$%i a
))
1031 (sub s
(mul 2 '$%i a
)))))))
1032 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1034 ((setq l1
(m2-sum-with-exp-case5 u var
))
1035 ;; c * t^-1 * (1 - %e^(2*a*t)) + d
1036 (let ((c (cdras 'c l1
))
1039 (add (mul c
(take '(%log
) (div (sub s a
) s
)))
1040 (defintegrate (mul d
(power '$%e
(mul -
1 s var
))) var
))))
1043 ;; At this point we expand the integrand.
1044 (distrdefexecinit ($expand form
) var
))))))
1046 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1048 ;;; Five pattern to find sums of Exponential functions which we can integrate.
1050 ;; Case 1: c * u^-1 * (%e^(-a*u) - %e^(-b*u))
1051 (defun m2-sum-with-exp-case1 (expr var
)
1056 ((mexpt) (u alike1
,var
) -
1)
1059 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1060 ((coeffpp) (z1 zerp
)))))
1062 (c2 equal-times-minus-one c
)
1063 ((mexpt) (u alike1
,var
) -
1)
1066 ((coeffpt) (b nonzerp
) (u alike1
,var
))
1067 ((coeffpp) (z2 zerp
)))))
1068 ((coeffpp) (d true
)))))
1070 ;; Case 2: c * u^(-3/2) * (%e^(-a*u) - %e^(-b*u))
1071 (defun m2-sum-with-exp-case2 (expr var
)
1076 ((mexpt) (u alike1
,var
) ((rat) -
3 2))
1079 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1080 ((coeffpp) (z1 zerp
)))))
1082 (c2 equal-times-minus-one c
)
1083 ((mexpt) (u alike1
,var
) ((rat) -
3 2))
1086 ((coeffpt) (b nonzerp
) (u alike1
,var
))
1087 ((coeffpp) (z2 zerp
)))))
1088 ((coeffpp) (d true
)))))
1090 ;; Case 3: c * u^-2 * (1 - 2 * %e^(-a*u) + %e^(2*a*u))
1091 (defun m2-sum-with-exp-case3 (expr var
)
1096 ((mexpt) (u alike1
,var
) -
2))
1099 ((mexpt) (u alike1
,var
) -
2)
1102 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1103 ((coeffpp) (z1 zerp
)))))
1105 (c3 equal-times-minus-two c
)
1106 ((mexpt) (u alike1
,var
) -
2)
1109 ((coeffpt) (b equal-div-two a
) (u alike1
,var
))
1110 ((coeffpp) (z2 zerp
)))))
1111 ((coeffpp) (d true
)))))
1113 ;; Case 4: c * t^-1 * (1 - 2 * %e^(-a*t) + %e^(2*a*t))
1114 (defun m2-sum-with-exp-case4 (expr var
)
1119 ((mexpt) (u alike1
,var
) -
1))
1122 ((mexpt) (u alike1
,var
) -
1)
1125 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1126 ((coeffpp) (z1 zerp
)))))
1128 (c3 equal-times-minus-two c
)
1129 ((mexpt) (u alike1
,var
) -
1)
1132 ((coeffpt) (b equal-div-two a
) (u alike1
,var
))
1133 ((coeffpp) (z2 zerp
)))))
1134 ((coeffpp) (d true
)))))
1136 ;; Case 5: c* t^-1 * (1 - %e^(2*a*t))
1137 (defun m2-sum-with-exp-case5 (expr var
)
1142 ((mexpt) (u alike1
,var
) -
1))
1144 (c2 equal-times-minus-one c
)
1145 ((mexpt) (u alike1
,var
) -
1)
1148 ((coeffpt) (a nonzerp
) (u alike1
,var
))
1149 ((coeffpp) (z1 zerp
)))))
1150 ((coeffpp) (d true
)))))
1152 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1154 ;;; Test functions for the pattern m2-sum-with-exp-case<n>
1156 (defun equal-times-minus-one (a b
)
1157 (equal a
(mul -
1 b
)))
1159 (defun equal-times-minus-two (a b
)
1160 (equal a
(mul -
2 b
)))
1162 (defun equal-div-two (a b
)
1163 (equal a
(div b
2)))
1165 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1167 ;;; Called by defintegrate.
1168 ;;; Call for every term of a sum defexec and add up the results.
1170 ;; Evaluate the transform of a sum as sum of transforms.
1171 (defun distrdefexecinit (expr var
)
1172 (cond ((equal (caar expr
) 'mplus
)
1173 (distrdefexec (cdr expr
) var
))
1174 (t (defexec expr var
))))
1176 ;; FUN is a list of addends. Compute the transform of each addend and
1178 (defun distrdefexec (expr var
)
1179 (cond ((null expr
) 0)
1180 (t (add (defexec (car expr
) var
)
1181 (distrdefexec (cdr expr
) var
)))))
1183 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1185 ;;; Called after transformation of a integrand to a new representation.
1186 ;;; Evalutate the tansform of a sum as sum of transforms.
1188 (defun sendexec (r a
)
1189 (distrexecinit ($expand
(mul (init r
) a
))))
1191 ;; Compute r*exp(-p*t), where t is the variable of integration and
1192 ;; p is the parameter of the Laplace transform.
1194 (mul r
(power '$%e
(mul -
1 *var
* *par
*))))
1196 (defun distrexecinit (expr)
1197 (cond ((and (consp expr
)
1199 (equal (caar expr
) 'mplus
))
1200 (distrexec (cdr expr
)))
1201 (t (hypgeo-exec expr
))))
1203 (defun distrexec (expr)
1204 (cond ((null expr
) 0)
1205 (t (add (hypgeo-exec (car expr
))
1206 (distrexec (cdr expr
))))))
1208 ;; It dispatches according to the kind of transform it matches.
1209 (defun hypgeo-exec (expr)
1211 (cond ((setq l
(m2-ltep expr
*var
* *par
*))
1212 (setq u
(cdras 'u l
)
1217 (return (ltscale u c a e f
))))
1218 (return (setq *hyp-return-noun-flag
* 'other-trans-to-follow
))))
1220 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1222 ;;; Compute transform of EXP wrt the variable of integration VAR.
1224 (defun defexec (expr var
)
1225 (let* ((*par
* 'psey
) ; Set parameter of Laplace transform
1226 (*var
* var
) ; Set variable of integration
1227 (*hyp-return-noun-flag
* nil
) ; Reset the flag
1229 (l (m2-defltep expr var
))
1230 (s (cdras 'a l
))) ; Get the parameter of the Laplace transform.
1232 ;; If we have not found a parameter, we try to factor the integrand.
1234 (when (and (numberp s
) (equal s
0))
1235 (setq l
(m2-defltep ($factor form
) *var
*))
1236 (setq s
(cdras 'a l
)))
1239 ;; EXP is an expression of the form u*%e^(s*t+e*f+c). So s
1240 ;; is basically the parameter of the Laplace transform.
1241 (let ((result (negtest l s
)))
1242 ;; At this point we construct the noun form if one of the
1243 ;; called routines set the global flag. If the global flag
1244 ;; is not set, the noun form has been already constructed.
1245 (if (and *hyp-return-noun-form-p
* *hyp-return-noun-flag
*)
1246 (list '(%specint simp
) expr
*var
*)
1249 ;; If necessary we construct the noun form.
1250 (if *hyp-return-noun-form-p
*
1251 (list '(%specint simp
) expr
*var
*)
1252 'other-defint-to-follow-defexec
)))))
1254 ;; L is the integrand of the transform, after pattern matching. S is
1255 ;; the parameter (p) of the transform.
1256 (defun negtest (l s
)
1258 (cond ((eq ($asksign
($realpart s
)) '$neg
)
1259 ;; The parameter of transform must have a negative
1260 ;; realpart. Break out the integrand into its various
1262 (setq u
(cdras 'u l
)
1266 (when (equal e
0) (setq f
1))
1267 ;; To compute the transform, we replace A with PSEY for
1268 ;; simplicity. After the transform is computed, replace
1271 ;; FIXME: Sometimes maxima will ask for the sign of PSEY.
1272 ;; But that doesn't occur in the original expression, so
1273 ;; it's very confusing. What should we do?
1275 ;; We know psey must be positive. psey is a substitution
1276 ;; for the paratemter a and we have checked the sign.
1277 ;; So it is the best to add a rule for the sign of psey.
1279 (mfuncall '$assume
`((mgreaterp) ,*par
* 0))
1286 (ltscale u c
0 e f
))
1288 ;; We forget the rule after finishing the calculation.
1289 (mfuncall '$forget
`((mgreaterp) ,*par
* 0))))))
1292 (setq *hyp-return-noun-flag
* 'other-defint-to-follow-negtest
))))
1294 ;; Compute the transform of
1296 ;; U * %E^(-VAR * (*PAR* - PAR0) + E*F + C)
1297 (defun ltscale (u c par0 e f
)
1299 (substl (sub *par
* par0
) *par
* (lt-exec u e f
))))
1301 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1303 ;;; Compute the transform of u*%e^(-p*t+e*f)
1305 (defun lt-exec (u e f
)
1307 (cond ((setq l
(m2-sum u
*var
*))
1308 ;; We have found a summation.
1311 (sendexec 1 (cdras 'u l
))
1316 ((setq l
(m2-unit_step u
*var
*))
1317 ;; We have found the Unit Step function.
1318 (setq u
(cdras 'u l
)
1320 (mul (power '$%e
(mul a
*par
*))
1321 (sendexec (cond (($freeof
*var
* u
) u
)
1322 (t (maxima-substitute (sub *var
* a
) *var
* u
)))
1326 ;; The simple case of u*%e^(-p*t)
1328 ((and (not (equal e
0))
1329 (setq l
(m2-c*t^v u
*var
*)))
1330 ;; We have u*%e^(-p*t+e*f). Try to see if U is of the form
1331 ;; c*t^v. If so, we can handle it here.
1334 ;; The complicated case. Remove the e*f term and move it to u.
1335 (lt-sf-log (mul u
(power '$%e
(mul e f
))))))))
1337 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1339 ;;; Pattern for the routine lt-exec
1341 ;; Recognize c*sum(u,index,low,high)
1342 (defun m2-sum (expr var
)
1347 ((%sum
) (u true
) (i true
) (l true
) (h true
)))
1348 ((coeffpp) (d zerp
)))))
1350 ;; Recognize u(t)*unit_step(x-a)
1351 (defun m2-unit_step (expr var
)
1356 (($unit_step
) ((mplus) (x alike1
,var
) ((coeffpp) (a true
)))))
1357 ((coeffpp) (d zerp
)))))
1360 ;; This is a duplicate of m2-arbpow1. Look if we can use it.
1361 (defun m2-c*t^v
(expr var
)
1364 ((coefftt) (c free
,var
))
1365 ((mexpt) (u alike1
,var
) (v free
,var
)))))
1367 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1369 ;;; Algorithm 1: Laplace transform of c*t^v*exp(-s*t+e*f)
1371 ;;; L contains the pattern for c*t^v.
1372 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1374 (defun lt-exp (l e f
)
1375 (let ((c (cdras 'c l
))
1377 (cond ((m2-t^
2 f
*var
*)
1378 (setq e
(inv (mul -
8 e
)) v
(add v
1))
1379 (f24p146test c v e
))
1380 ((m2-sqroott f
*var
*)
1381 ;; We don't do the transformation at this place. Because we take the
1382 ;; square of e we lost the sign and get wrong results.
1383 ;(setq e (mul* e e (inv 4)) v (add v 1))
1384 (f35p147test c v e
))
1386 (setq e
(mul -
4 e
) v
(add v
1))
1387 (f29p146test c v e
)) ; We have to call with the constant c.
1388 ((and (equal v
0) ; We have to test for v=0 and to call
1390 (f36p147 c e
)) ; with the constant c.
1391 ((and (equal v
0) (m2-e^t f
*var
*))
1392 (f37p147 c
(mul -
1 e
)))
1394 (setq *hyp-return-noun-flag
* 'other-lt-exponential-to-follow
)))))
1396 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1398 ;;; Pattern for the routine lt-exp
1401 (defun m2-t^
2 (expr var
)
1402 (m2 expr
`((mexpt) (u alike1
,var
) 2)))
1404 ;; Recognize sqrt(t)
1405 (defun m2-sqroott (expr var
)
1406 (m2 expr
`((mexpt) (u alike1
,var
) ((rat) 1 2))))
1409 (defun m2-t^-
1 (expr var
)
1410 (m2 expr
`((mexpt) (u alike1
,var
) -
1)))
1413 (defun m2-e^-t
(expr var
)
1414 (m2 expr
`((mexpt) $%e
((mtimes) -
1 (u alike1
,var
)))))
1417 (defun m2-e^t
(expr var
)
1418 (m2 expr
`((mexpt) $%e
(u alike1
,var
))))
1420 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1422 ;;; Algorithm 1.1: Laplace transform of c*t^v*exp(-a*t^2)
1424 ;;; Table of Integral Transforms
1426 ;;; p. 146, formula 24:
1428 ;;; t^(v-1)*exp(-t^2/8/a)
1429 ;;; -> gamma(v)*2^v*a^(v/2)*exp(a*p^2)*D[-v](2*p*sqrt(a))
1431 ;;; Re(a) > 0, Re(v) > 0
1432 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1434 (defun f24p146test (c v a
)
1435 (cond ((and (eq ($asksign a
) '$pos
)
1436 (eq ($asksign v
) '$pos
))
1437 ;; Both a and v must be positive
1440 (setq *hyp-return-noun-flag
* 'fail-on-f24p146test
))))
1442 (defun f24p146 (c v a
)
1447 (power '$%e
(mul a
*par
* *par
*))
1448 (dtford (mul 2 *par
* (power a
'((rat simp
) 1 2)))
1451 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1453 ;;; Algorithm 1.2: Laplace transform of c*t^v*exp(-a*sqrt(t))
1455 ;;; Table of Integral Transforms
1457 ;;; p. 147, formula 35:
1459 ;;; (2*t)^(v-1)*exp(-2*sqrt(a)*sqrt(t))
1460 ;;; -> gamma(2*v)*p^(-v)*exp(a/p/2)*D[-2*v](sqrt(2*a/p))
1462 ;;; Re(v) > 0, Re(p) > 0
1463 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1465 ;; Check if conditions for f35p147 hold
1466 (defun f35p147test (c v a
)
1467 (cond ((eq ($asksign
(add v
1)) '$pos
)
1468 ;; v must be positive
1471 ;; Set a global flag. When *hyp-return-noun-form-p* is T the noun
1472 ;; form will be constructed in the routine DEFEXEC.
1473 (setq *hyp-return-noun-flag
* 'fail-on-f35p147test
))))
1475 (defun f35p147 (c v a
)
1476 ;; We have not done the calculation v->v+1 and a-> a^2/4
1477 ;; and subsitute here accordingly.
1478 (let ((v (add v
1)))
1480 (take '(%gamma
) (add v v
))
1481 (power 2 (sub 1 v
)) ; Is this supposed to be here?
1482 (power *par
* (mul -
1 v
))
1483 (power '$%e
(mul a a
'((rat simp
) 1 8) (inv *par
*)))
1484 ;; We need an additional factor -1 to get the expected results.
1485 ;; What is the mathematically reason?
1486 (dtford (mul -
1 a
(inv (power (mul 2 *par
*) '((rat simp
) 1 2))))
1489 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1491 ;; Express a parabolic cylinder function as either a parabolic
1492 ;; cylinder function or as hypergeometric function.
1494 ;; Tables of Integral Transforms, p. 386 gives this definition:
1496 ;; D[v](z) = 2^(v/2+1/4)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
1499 (let ((inv4 (inv 4)))
1500 (cond ((or $prefer_d
1501 (whittindtest (add (div v
2) inv4
) inv4
))
1502 ;; At this time the Parabolic Cylinder D function is not implemented
1503 ;; as a simplifying function. We call nevertheless the simplifer
1504 ;; to simplify the arguments. When we implement the function
1505 ;; The symbol has to be changed to the noun form.
1506 (take '($parabolic_cylinder_d
) v z
))
1507 (t (simpdtf z v
)))))
1509 (defun whittindtest (i1 i2
)
1510 (or (maxima-integerp (add i2 i2
))
1511 (neginp (sub (sub '((rat simp
) 1 2) i2
) i1
))
1512 (neginp (sub (add '((rat simp
) 1 2) i2
) i1
))))
1514 ;; Return T if a is a non-positive integer.
1515 ;; (Do we really want maxima-integerp or hyp-integerp here?)
1517 (cond ((maxima-integerp a
)
1519 (eq ($sign a
) '$neg
)))))
1521 ;; Express parabolic cylinder function as a hypergeometric function.
1525 ;; U(a,x) = D[-a-1/2](x)
1527 ;; and A&S 19.12.3 gives
1529 ;; U(a,+/-x) = sqrt(%pi)*2^(-1/4-a/2)*exp(-x^2/4)/gamma(3/4+a/2)
1530 ;; *M(a/2+1/4,1/2,x^2/2)
1531 ;; -/+ sqrt(%pi)*2^(1/4-a/2)*x*exp(-x^2/4)/gamma(1/4+a/2)
1532 ;; *M(a/2+3/4,3/2,x^2/2)
1536 ;; D[v](z) = U(-v-1/2,z)
1537 ;; = sqrt(%pi)*2^(v/2+1/2)*x*exp(-x^2/4)
1538 ;; *M(1/2-v/2,3/2,x^2/2)/gamma(-v/2)
1539 ;; + sqrt(%pi)*2^(v/2)*exp(-x^2/4)/gamma(1/2-v/2)
1540 ;; *M(-v/2,1/2,x^2/2)
1542 (defun simpdtf (z v
)
1543 (let ((inv2 '((rat simp
) 1 2))
1544 (pow (power '$%e
(mul z z
'((rat simp
) -
1 4)))))
1545 (add (mul (power 2 (div (sub v
1) 2))
1547 -
2 (power '$%pi inv2
) ; gamma(-1/2)
1548 (inv (take '(%gamma
) (mul v -
1 inv2
)))
1550 (hgfsimp-exec (list (sub inv2
(div v
2)))
1551 (list '((rat simp
) 3 2))
1553 (mul (power 2 (div v
2))
1554 (power '$%pi inv2
) ; gamma(1/2)
1556 (inv (take '(%gamma
) (sub inv2
(mul v inv2
))))
1557 (hgfsimp-exec (list (mul v -
1 inv2
))
1561 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1563 ;;; Algorithm 1.3: Laplace transform of t^v*exp(1/t)
1565 ;;; Table of Integral Transforms
1567 ;;; p. 146, formula 29:
1569 ;;; t^(v-1)*exp(-a/t/4)
1570 ;;; -> 2*(a/p/4)^(v/2)*bessel_k(v, sqrt(a)*sqrt(p))
1573 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1575 ;; Check if conditions for f29p146test hold
1576 (defun f29p146test (c v a
)
1577 (cond ((eq ($asksign a
) '$pos
)
1580 (setq *hyp-return-noun-flag
* 'fail-on-f29p146test
))))
1582 (defun f29p146 (c v a
)
1584 (power (mul a
'((rat simp
) 1 4) (inv *par
*))
1588 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1590 ;; bessel_k(v, sqrt(a)*sqrt(p)) in terms of bessel_k or in terms of
1591 ;; hypergeometric functions.
1593 ;; Choose bessel_k if the order v is an integer. (Why?)
1596 (let ((z (power (mul a
*par
*) '((rat simp
) 1 2))))
1597 (cond ((maxima-integerp v
)
1598 (take '(%bessel_k
) v z
))
1602 ;; Express the Bessel K function in terms of hypergeometric functions.
1604 ;; K[v](z) = %pi/2*(bessel_i(-v,z)-bessel(i,z))/sin(v*%pi)
1608 ;; bessel_i(v,z) = (z/2)^v/gamma(v+1)*0F1(;v+1;z^2/4)
1610 (defun simpktf (z v
)
1611 (let ((dz2 (div z
2)))
1614 (inv (take '(%sin
) (mul v
'$%pi
)))
1615 (sub (mul (power dz2
(mul -
1 v
))
1616 (inv (take '(%gamma
) (sub 1 v
)))
1619 (mul z z
'((rat simp
) 1 4))))
1621 (inv (take '(%gamma
) (add v
1)))
1624 (mul z z
'((rat simp
) 1 4))))))))
1626 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1628 ;;; Algorithm 1.4: Laplace transform of exp(exp(-t))
1630 ;;; Table of Integral Transforms
1632 ;;; p. 147, formula 36:
1635 ;;; -> a^(-p)*gamma(p,a)
1636 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1638 (defun f36p147 (c a
)
1639 (let ((-a (mul -
1 a
)))
1641 (power -a
(mul -
1 *par
*))
1642 `((%gamma_incomplete_lower simp
) ,*par
* ,-a
))))
1644 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1646 ;;; Algorithm 1.5: Laplace transform of exp(exp(t))
1648 ;;; Table of Integral Transforms
1650 ;;; p. 147, formula 36:
1653 ;;; -> a^(-p)*gamma_incomplete(-p,a)
1654 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1656 (defun f37p147 (c a
)
1659 (take '(%gamma_incomplete
) (mul -
1 *par
*) a
)))
1661 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1663 ;;; Algorithm 2: Laplace transform of u(t)*%e^(-p*t).
1665 ;;; u contains one or more special functions. Dispatches according to the
1666 ;;; special functions involved in the Laplace transformable expression.
1668 ;;; We have three general types of integrands:
1670 ;;; 1. Call a function to return immediately the Laplace transform,
1671 ;;; e.g. call lt-arbpow, lt-arbpow2, lt-log, whittest to return the
1672 ;;; Laplace transform.
1673 ;;; 2. Call lt-ltp directly or via an "expert function on Laplace transform",
1674 ;;; transform the special function to a representation in terms of one
1675 ;;; hypergeometric function and do the integration
1676 ;;; e.g. for a direct call of lt-ltp asin, atan or via lt2j for
1677 ;;; an integrand with two bessel function.
1678 ;;; 3. Call fractest, fractest1, ... which transform the involved special
1679 ;;; function to a new representation. Send the transformed expression with
1680 ;;; the routine sendexec to the integrator and try to integrate the new
1681 ;;; representation, e.g. gamma_incomplete is first transformed to a new
1684 ;;; The ordering of the calls to match a pattern is important.
1685 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
1687 (defun lt-sf-log (u)
1688 (prog (l index1 index11 index2 index21 arg1 arg2 rest
)
1690 ;; Laplace transform of asin(w)
1691 (cond ((setq l
(m2-asin u
*var
*))
1692 (setq arg1
(cdras 'w l
)
1694 (return (lt-ltp 'asin rest arg1 nil
))))
1696 ;; Laplace transform of atan(w)
1697 (cond ((setq l
(m2-atan u
*var
*))
1698 (setq arg1
(cdras 'w l
)
1700 (return (lt-ltp 'atan rest arg1 nil
))))
1702 ;; Laplace transform of two Bessel J functions
1703 (cond ((setq l
(m2-twoj u
*var
*))
1704 (setq index1
(cdras 'v1 l
)
1705 index2
(cdras 'v2 l
)
1709 (return (lt2j rest arg1 arg2 index1 index2
))))
1711 ;; Laplace transform of two hankel_1 functions
1712 (cond ((setq l
(m2-two-hankel_1 u
*var
*))
1713 (setq index1
(cdras 'v1 l
)
1714 index2
(cdras 'v2 l
)
1718 ;; Call the code for the symbol %h.
1719 (return (fractest rest arg1 arg2 index1
1 index2
1 '2htjory
))))
1721 ;; Laplace transform of two hankel_2 functions
1722 (cond ((setq l
(m2-two-hankel_2 u
*var
*))
1723 (setq index1
(cdras 'v1 l
)
1724 index2
(cdras 'v2 l
)
1728 ;; Call the code for the symbol %h.
1729 (return (fractest rest arg1 arg2 index1
2 index2
2 '2htjory
))))
1731 ;; Laplace transform of hankel_1 * hankel_2
1732 (cond ((setq l
(m2-hankel_1*hankel_2 u
*var
*))
1733 (setq index1
(cdras 'v1 l
)
1734 index2
(cdras 'v2 l
)
1738 ;; Call the code for the symbol %h.
1739 (return (fractest rest arg1 arg2 index1
1 index2
2 '2htjory
))))
1741 ;; Laplace transform of two Bessel Y functions
1742 (cond ((setq l
(m2-twoy u
*var
*))
1743 (setq index1
(cdras 'v1 l
)
1744 index2
(cdras 'v2 l
)
1748 (return (fractest rest arg1 arg2 index1 nil index2 nil
'2ytj
))))
1750 ;; Laplace transform of two Bessel K functions
1751 (cond ((setq l
(m2-twok u
*var
*))
1752 (setq index1
(cdras 'v1 l
)
1753 index2
(cdras 'v2 l
)
1757 (return (fractest rest arg1 arg2 index1 nil index2 nil
'2kti
))))
1759 ;; Laplace transform of Bessel K and Bessel Y functions
1760 (cond ((setq l
(m2-onekoney u
*var
*))
1761 (setq index1
(cdras 'v1 l
)
1762 index2
(cdras 'v2 l
)
1766 (return (fractest rest arg1 arg2 index1 nil index2 nil
'ktiytj
))))
1768 ;; Laplace transform of Bessel I and Bessel J functions
1769 (cond ((setq l
(m2-oneionej u
*var
*))
1770 (setq index1
(cdras 'v1 l
)
1771 index2
(cdras 'v2 l
)
1772 index21
(cdras 'v21 l
)
1773 arg1
(mul '$%i
(cdras 'w1 l
))
1775 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1776 (return (lt2j rest arg1 arg2 index1 index2
))))
1778 ;; Laplace transform of Bessel I and Hankel 1 functions
1779 (cond ((setq l
(m2-bessel_i*hankel_1 u
*var
*))
1780 (setq index1
(cdras 'v1 l
)
1781 index2
(cdras 'v2 l
)
1782 arg1
(mul '$%i
(cdras 'w1 l
))
1784 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1785 (return (fractest1 rest arg1 arg2 index1 index2
1 'besshtjory
))))
1787 ;; Laplace transform of Bessel I and Hankel 2 functions
1788 (cond ((setq l
(m2-bessel_i*hankel_2 u
*var
*))
1789 (setq index1
(cdras 'v1 l
)
1790 index2
(cdras 'v2 l
)
1791 arg1
(mul '$%i
(cdras 'w1 l
))
1793 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1794 (return (fractest1 rest arg1 arg2 index1 index2
2 'besshtjory
))))
1796 ;; Laplace transform of Bessel Y and Bessel J functions
1797 (cond ((setq l
(m2-oneyonej u
*var
*))
1798 (setq index1
(cdras 'v1 l
)
1799 index2
(cdras 'v2 l
)
1803 (return (fractest1 rest arg2 arg1 index2 index1 nil
'bessytj
))))
1805 ;; Laplace transform of Bessel K and Bessel J functions
1806 (cond ((setq l
(m2-onekonej u
*var
*))
1807 (setq index1
(cdras 'v1 l
)
1808 index2
(cdras 'v2 l
)
1812 (return (fractest1 rest arg2 arg1 index2 index1 nil
'besskti
))))
1814 ;; Laplace transform of Hankel 1 and Bessel J functions
1815 (cond ((setq l
(m2-hankel_1*bessel_j u
*var
*))
1816 (setq index1
(cdras 'v1 l
)
1817 index2
(cdras 'v2 l
)
1821 (return (fractest1 rest arg2 arg1 index2 index1
1 'besshtjory
))))
1823 ;; Laplace transform of Hankel 2 and Bessel J functions
1824 (cond ((setq l
(m2-hankel_2*bessel_j u
*var
*))
1825 (setq index1
(cdras 'v1 l
)
1826 index2
(cdras 'v2 l
)
1830 (return (fractest1 rest arg2 arg1 index2 index1
2 'besshtjory
))))
1832 ;; Laplace transform of Bessel Y and Hankel 1 functions
1833 (cond ((setq l
(m2-bessel_y*hankel_1 u
*var
*))
1834 (setq index1
(cdras 'v1 l
)
1835 index2
(cdras 'v2 l
)
1839 (return (fractest1 rest arg2 arg1 index2 index1
1 'htjoryytj
))))
1841 ;; Laplace transform of Bessel Y and Hankel 2 functions
1842 (cond ((setq l
(m2-bessel_y*hankel_2 u
*var
*))
1843 (setq index1
(cdras 'v1 l
)
1844 index2
(cdras 'v2 l
)
1848 (return (fractest1 rest arg2 arg1 index2 index1
2 'htjoryytj
))))
1850 ;; Laplace transform of Bessel K and Hankel 1 functions
1851 (cond ((setq l
(m2-bessel_k*hankel_1 u
*var
*))
1852 (setq index1
(cdras 'v1 l
)
1853 index2
(cdras 'v2 l
)
1857 (return (fractest1 rest arg2 arg1 index2 index1
1 'htjorykti
))))
1859 ;; Laplace transform of Bessel K and Hankel 2 functions
1860 (cond ((setq l
(m2-bessel_k*hankel_2 u
*var
*))
1861 (setq index1
(cdras 'v1 l
)
1862 index2
(cdras 'v2 l
)
1866 (return (fractest1 rest arg2 arg1 index2 index1
2 'htjorykti
))))
1868 ;; Laplace transform of Bessel I and Bessel Y functions
1869 (cond ((setq l
(m2-oneioney u
*var
*))
1870 (setq index1
(cdras 'v1 l
)
1871 index2
(cdras 'v2 l
)
1872 arg1
(mul '$%i
(cdras 'w1 l
))
1874 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1875 (return (fractest1 rest arg1 arg2 index1 index2 nil
'bessytj
))))
1877 ;; Laplace transform of Bessel I and Bessel K functions
1878 (cond ((setq l
(m2-oneionek u
*var
*))
1879 (setq index1
(cdras 'v1 l
)
1880 index2
(cdras 'v2 l
)
1881 arg1
(mul '$%i
(cdras 'w1 l
))
1883 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
1884 (return (fractest1 rest arg1 arg2 index1 index2 nil
'besskti
))))
1886 ;; Laplace transform of Struve H function
1887 (cond ((setq l
(m2-struve_h u
*var
*))
1888 (setq index1
(cdras 'v l
)
1891 (return (lt1hstruve rest arg1 index1
))))
1893 ;; Laplace transform of Struve L function
1894 (cond ((setq l
(m2-struve_l u
*var
*))
1895 (setq index1
(cdras 'v l
)
1898 (return (lt1lstruve rest arg1 index1
))))
1900 ;; Laplace transform of little Lommel s function
1901 (cond ((setq l
(m2-ones u
*var
*))
1902 (setq index1
(cdras 'v1 l
)
1903 index2
(cdras 'v2 l
)
1906 (return (lt1s rest arg1 index1 index2
))))
1908 ;; Laplace transform of Lommel S function
1909 (cond ((setq l
(m2-oneslommel u
*var
*))
1910 (setq index1
(cdras 'v1 l
)
1911 index2
(cdras 'v2 l
)
1914 (return (fractest2 rest arg1 index1 index2
'slommel
))))
1916 ;; Laplace transform of Bessel Y function
1917 (cond ((setq l
(m2-oney u
*var
*))
1918 (setq index1
(cdras 'v l
)
1921 (return (lt1yref rest arg1 index1
))))
1923 ;; Laplace transform of Bessel K function
1924 (cond ((setq l
(m2-onek u
*var
*))
1925 (setq index1
(cdras 'v l
)
1928 (cond ((zerop1 index1
)
1929 ;; Special case for a zero index
1930 (return (lt-bessel_k0 rest arg1
)))
1932 (return (fractest2 rest arg1 index1 nil
'kti
))))))
1934 ;; Laplace transform of Parabolic Cylinder function
1935 (cond ((setq l
(m2-parabolic_cylinder_d u
*var
*))
1936 (setq index1
(cdras 'v l
)
1939 (return (fractest2 rest arg1 index1 nil
'd
))))
1941 ;; Laplace transform of Incomplete Gamma function
1942 (cond ((setq l
(m2-onegammaincomplete u
*var
*))
1943 (setq arg1
(cdras 'w1 l
)
1946 (return (fractest2 rest arg1 arg2 nil
'gamma_incomplete
))))
1948 ;; Laplace transform of Batemann function
1949 (cond ((setq l
(m2-onekbateman u
*var
*))
1950 (setq index1
(cdras 'v l
)
1953 (return (fractest2 rest arg1 index1 nil
'kbateman
))))
1955 ;; Laplace transform of Bessel J function
1956 (cond ((setq l
(m2-onej u
*var
*))
1957 (setq index1
(cdras 'v l
)
1960 (return (lt1j rest arg1 index1
))))
1962 ;; Laplace transform of lower incomplete Gamma function
1963 (cond ((setq l
(m2-onegamma-incomplete-lower u
*var
*))
1964 (setq arg1
(cdras 'w1 l
)
1967 (return (lt1gamma-incomplete-lower rest arg1 arg2
))))
1969 ;; Laplace transform of Hankel 1 function
1970 (cond ((setq l
(m2-hankel_1 u
*var
*))
1971 (setq index1
(cdras 'v l
)
1974 (return (fractest2 rest arg1 index1
1 'htjory
))))
1976 ;; Laplace transform of Hankel 2 function
1977 (cond ((setq l
(m2-hankel_2 u
*var
*))
1978 (setq index1
(cdras 'v l
)
1981 (return (fractest2 rest arg1 index1
2 'htjory
))))
1983 ;; Laplace transform of Whittaker M function
1984 (cond ((setq l
(m2-onem u
*var
*))
1985 (setq index1
(cdras 'v1 l
)
1986 index11
(cdras 'v2 l
)
1989 (return (lt1m rest arg1 index1 index11
))))
1991 ;; Laplace transform of Whittaker M function
1992 (cond ((setq l
(m2-whittaker_m u
*var
*))
1993 (setq index1
(cdras 'v1 l
)
1994 index2
(cdras 'v2 l
)
1997 (return (lt1m rest arg1 index1 index2
))))
1999 ;; Laplace transform of the Generalized Laguerre function, %l[v1,v2](w)
2000 (cond ((setq l
(m2-onel u
*var
*))
2001 (setq index1
(cdras 'v1 l
)
2002 index11
(cdras 'v2 l
)
2005 (return (integertest rest arg1 index1 index11
'l
))))
2007 ;; Laplace transform for the Generalized Laguerre function
2008 ;; We call the routine for %l[v1,v2](w).
2009 (cond ((setq l
(m2-one-gen-laguerre u
*var
*))
2010 (setq index1
(cdras 'v1 l
)
2011 index2
(cdras 'v2 l
)
2014 (return (integertest rest arg1 index1 index2
'l
))))
2016 ;; Laplace transform for the Laguerre function
2017 ;; We call the routine for %l[v1,0](w).
2018 (cond ((setq l
(m2-one-laguerre u
*var
*))
2019 (setq index1
(cdras 'v1 l
)
2022 (return (integertest rest arg1 index1
0 'l
))))
2024 ;; Laplace transform of Gegenbauer function
2025 (cond ((setq l
(m2-onec u
*var
*))
2026 (setq index1
(cdras 'v1 l
)
2027 index11
(cdras 'v2 l
)
2030 (return (integertest rest arg1 index1 index11
'c
))))
2032 ;; Laplace transform of Chebyshev function of the first kind
2033 (cond ((setq l
(m2-onet u
*var
*))
2034 (setq index1
(cdras 'v1 l
)
2037 (return (integertest rest arg1 index1 nil
't
))))
2039 ;; Laplace transform of Chebyshev function of the second kind
2040 (cond ((setq l
(m2-oneu u
*var
*))
2041 (setq index1
(cdras 'v1 l
)
2044 (return (integertest rest arg1 index1 nil
'u
))))
2046 ;; Laplace transform for the Hermite function, hermite(index1,arg1)
2047 (cond ((setq l
(m2-one-hermite u
*var
*))
2048 (setq index1
(cdras 'v1 l
)
2052 (cond ((and (maxima-integerp index1
)
2055 ;; When index1 is an even or odd integer, we transform
2056 ;; directly to a hypergeometric function. For this case we
2057 ;; get a Laplace transform when the arg is the
2058 ;; square root of the variable.
2059 (sendexec rest
(hermite-to-hypergeometric index1 arg1
)))
2061 (integertest rest arg1 index1 nil
'he
))))))
2063 ;; Laplace transform of %p[v1,v2](w), Associated Legendre P function
2064 (cond ((setq l
(m2-hyp-onep u
*var
*))
2065 (setq index1
(cdras 'v1 l
)
2066 index11
(cdras 'v2 l
)
2069 (return (lt1p rest arg1 index1 index11
))))
2071 ;; Laplace transform of Associated Legendre P function
2072 (cond ((setq l
(m2-assoc_legendre_p u
*var
*))
2073 (setq index1
(cdras 'v1 l
)
2074 index2
(cdras 'v2 l
)
2077 (return (lt1p rest arg1 index1 index2
))))
2079 ;; Laplace transform of %p[v1,v2,v3](w), Jacobi function
2080 (cond ((setq l
(m2-onepjac u
*var
*))
2081 (setq index1
(cdras 'v1 l
)
2082 index2
(cdras 'v2 l
)
2083 index21
(cdras 'v3 l
)
2086 (return (pjactest rest arg1 index1 index2 index21
))))
2088 ;; Laplace transform of Jacobi P function
2089 (cond ((setq l
(m2-jacobi_p u
*var
*))
2090 (setq index1
(cdras 'v1 l
)
2091 index2
(cdras 'v2 l
)
2092 index21
(cdras 'v3 l
)
2095 (return (pjactest rest arg1 index1 index2 index21
))))
2097 ;; Laplace transform of Associated Legendre function of the second kind
2098 (cond ((setq l
(m2-oneq u
*var
*))
2099 (setq index1
(cdras 'v1 l
)
2100 index11
(cdras 'v2 l
)
2103 (return (lt1q rest arg1 index1 index11
))))
2105 ;; Laplace transform of Associated Legendre function of the second kind
2106 (cond ((setq l
(m2-assoc_legendre_q u
*var
*))
2107 (setq index1
(cdras 'v1 l
)
2108 index2
(cdras 'v2 l
)
2111 (return (lt1q rest arg1 index1 index2
))))
2113 ;; Laplace transform of %p[v1](w), Legendre P function
2114 (cond ((setq l
(m2-onep0 u
*var
*))
2115 (setq index1
(cdras 'v1 l
)
2119 (return (lt1p rest arg1 index1 index11
))))
2121 ;; Laplace transform of Legendre P function
2122 (cond ((setq l
(m2-legendre_p u
*var
*))
2123 (setq index1
(cdras 'v1 l
)
2126 (return (lt1p rest arg1 index1
0))))
2128 ;; Laplace transform of Whittaker W function
2129 (cond ((setq l
(m2-onew u
*var
*))
2130 (setq index1
(cdras 'v1 l
)
2131 index11
(cdras 'v2 l
)
2134 (return (whittest rest arg1 index1 index11
))))
2136 ;; Laplace transform of Whittaker W function
2137 (cond ((setq l
(m2-whittaker_w u
*var
*))
2138 (setq index1
(cdras 'v1 l
)
2139 index2
(cdras 'v2 l
)
2142 (return (whittest rest arg1 index1 index2
))))
2144 ;; Laplace transform of square of Bessel J function
2145 (cond ((setq l
(m2-onej^
2 u
*var
*))
2146 (setq index1
(cdras 'v l
)
2149 (return (lt1j^
2 rest arg1 index1
))))
2151 ;; Laplace transform of square of Hankel 1 function
2152 (cond ((setq l
(m2-hankel_1^
2 u
*var
*))
2153 (setq index1
(cdras 'v l
)
2156 (return (fractest rest arg1 arg1 index1
1 index1
1 '2htjory
))))
2158 ;; Laplace transform of square of Hankel 2 function
2159 (cond ((setq l
(m2-hankel_2^
2 u
*var
*))
2160 (setq index1
(cdras 'v l
)
2163 (return (fractest rest arg1 arg1 index1
2 index1
2 '2htjory
))))
2165 ;; Laplace transform of square of Bessel Y function
2166 (cond ((setq l
(m2-oney^
2 u
*var
*))
2167 (setq index1
(cdras 'v l
)
2170 (return (fractest rest arg1 arg1 index1 nil index1 nil
'2ytj
))))
2172 ;; Laplace transform of square of Bessel K function
2173 (cond ((setq l
(m2-onek^
2 u
*var
*))
2174 (setq index1
(cdras 'v l
)
2177 (return (fractest rest arg1 arg1 index1 nil index1 nil
'2kti
))))
2179 ;; Laplace transform of two Bessel I functions
2180 (cond ((setq l
(m2-twoi u
*var
*))
2181 (setq index1
(cdras 'v1 l
)
2182 index2
(cdras 'v2 l
)
2183 arg1
(mul '$%i
(cdras 'w1 l
))
2184 arg2
(mul '$%i
(cdras 'w2 l
))
2185 rest
(mul (power '$%i
(neg index1
))
2186 (power '$%i
(neg index1
))
2188 (return (lt2j rest arg1 arg2 index1 index2
))))
2190 ;; Laplace transform of Bessel I. We use I[v](w)=%i^n*J[n](%i*w).
2191 (cond ((setq l
(m2-onei u
*var
*))
2192 (setq index1
(cdras 'v l
)
2193 arg1
(mul '$%i
(cdras 'w l
))
2194 rest
(mul (power '$%i
(neg index1
)) (cdras 'u l
)))
2195 (return (lt1j rest arg1 index1
))))
2197 ;; Laplace transform of square of Bessel I function
2198 (cond ((setq l
(m2-onei^
2 u
*var
*))
2199 (setq index1
(cdras 'v l
)
2200 arg1
(mul '$%i
(cdras 'w l
))
2201 rest
(mul (power '$%i
(neg index1
))
2202 (power '$%i
(neg index1
))
2204 (return (lt1j^
2 rest arg1 index1
))))
2206 ;; Laplace transform of Erf function
2207 (cond ((setq l
(m2-onerf u
*var
*))
2208 (setq arg1
(cdras 'w l
)
2210 (return (lt1erf rest arg1
))))
2212 ;; Laplace transform of the logarithmic function.
2213 ;; We add an algorithm for the Laplace transform and call the routine
2214 ;; lt-log. The old code is still present, but isn't called.
2215 (cond ((setq l
(m2-onelog u
*var
*))
2216 (setq arg1
(cdras 'w l
)
2218 (return (lt-log rest arg1
))))
2220 ;; Laplace transform of Erfc function
2221 (cond ((setq l
(m2-onerfc u
*var
*))
2222 (setq arg1
(cdras 'w l
)
2224 (return (fractest2 rest arg1 nil nil
'erfc
))))
2226 ;; Laplace transform of expintegral_ei.
2227 ;; Maxima uses the build in transformation to the gamma_incomplete
2228 ;; function and simplifies the log functions of the transformation. We do
2229 ;; not use the dispatch mechanism of fractest2, but call sendexec directly
2230 ;; with the transformed function.
2231 (cond ((setq l
(m2-oneexpintegral_ei u
*var
*))
2232 (setq arg1
(cdras 'w l
)
2234 (let (($expintrep
'%gamma_incomplete
)
2236 (return (sratsimp (sendexec rest
(ftake '%expintegral_ei arg1
)))))))
2238 ;; Laplace transform of expintegral_e1
2239 (cond ((setq l
(m2-oneexpintegral_e1 u
*var
*))
2240 (setq arg1
(cdras 'w l
)
2242 (let (($expintrep
'%gamma_incomplete
)
2244 (return (sratsimp (sendexec rest
(ftake '%expintegral_e1 arg1
)))))))
2246 ;; Laplace transform of expintegral_e
2247 (cond ((setq l
(m2-oneexpintegral_e u
*var
*))
2248 (setq arg1
(cdras 'v l
)
2251 (let (($expintrep
'%gamma_incomplete
)
2253 (return (sratsimp (sendexec rest
(ftake '%expintegral_e arg1 arg2
)))))))
2255 ;; Laplace transform of expintegral_si
2256 (cond ((setq l
(m2-oneexpintegral_si u
*var
*))
2257 (setq arg1
(cdras 'w l
)
2259 ;; We transform to the hypergeometric representation.
2261 (sendexec rest
(expintegral_si-to-hypergeometric arg1
)))))
2263 ;; Laplace transform of expintegral_shi
2264 (cond ((setq l
(m2-oneexpintegral_shi u
*var
*))
2265 (setq arg1
(cdras 'w l
)
2267 ;; We transform to the hypergeometric representation.
2269 (sendexec rest
(expintegral_shi-to-hypergeometric arg1
)))))
2271 ;; Laplace transform of expintegral_ci
2272 (cond ((setq l
(m2-oneexpintegral_ci u
*var
*))
2273 (setq arg1
(cdras 'w l
)
2275 ;; We transform to the hypergeometric representation.
2276 ;; Because we have Logarithmic terms in the transformation,
2277 ;; we switch on the flag $logexpand and do a ratsimp.
2278 (let (($logexpand
'$super
))
2281 (sendexec rest
(expintegral_ci-to-hypergeometric arg1
)))))))
2283 ;; Laplace transform of expintegral_chi
2284 (cond ((setq l
(m2-oneexpintegral_chi u
*var
*))
2285 (setq arg1
(cdras 'w l
)
2287 ;; We transform to the hypergeometric representation.
2288 ;; Because we have Logarithmic terms in the transformation,
2289 ;; we switch on the flag $logexpand and do a ratsimp.
2290 (let (($logexpand
'$super
))
2293 (sendexec rest
(expintegral_chi-to-hypergeometric arg1
)))))))
2295 ;; Laplace transform of Complete elliptic integral of the first kind
2296 (cond ((setq l
(m2-onekelliptic u
*var
*))
2297 (setq arg1
(cdras 'w l
)
2299 (return (lt1kelliptic rest arg1
))))
2301 ;; Laplace transform of Complete elliptic integral of the first kind
2302 (cond ((setq l
(m2-elliptic_kc u
*var
*))
2303 (setq arg1
(cdras 'w l
)
2305 (return (lt1kelliptic rest arg1
))))
2307 ;; Laplace transform of Complete elliptic integral of the second kind
2308 (cond ((setq l
(m2-onee u
*var
*))
2309 (setq arg1
(cdras 'w l
)
2311 (return (lt1e rest arg1
))))
2313 ;; Laplace transform of Complete elliptic integral of the second kind
2314 (cond ((setq l
(m2-elliptic_ec u
*var
*))
2315 (setq arg1
(cdras 'w l
)
2317 (return (lt1e rest arg1
))))
2319 ;; Laplace transform of %f[v1,v2](w1,w2,w3), Hypergeometric function
2320 ;; We support the Laplace transform of the build in symbol %f. We do
2321 ;; not use the mechanism of defining an "Expert on Laplace transform",
2322 ;; the expert function does a call to lt-ltp. We do this call directly.
2323 (cond ((setq l
(m2-onef u
*var
*))
2324 (setq rest
(cdras 'u l
)
2326 index1
(list (cdras 'w1 l
) (cdras 'w2 l
)))
2327 (return (lt-ltp 'f rest arg1 index1
))))
2329 ;; Laplace transform of Hypergeometric function
2330 (cond ((setq l
(m2-hypergeometric u
*var
*))
2331 (setq rest
(cdras 'u l
)
2333 index1
(list (cdras 'w1 l
) (cdras 'w2 l
)))
2334 (return (lt-ltp 'f rest arg1 index1
))))
2336 ;; Laplace transform of c * t^v * (a+t)^w
2337 ;; It is possible to combine arbpow2 and arbpow.
2338 (cond ((setq l
(m2-arbpow2 u
*var
*))
2339 (setq rest
(cdras 'c l
)
2343 index2
(cdras 'w l
))
2344 (return (lt-arbpow2 rest arg1 arg2 index1 index2
))))
2346 ;; Laplace transform of c * t^v
2347 (cond ((setq l
(m2-arbpow1 u
*var
*))
2348 (setq arg1
(cdras 'u l
)
2350 index1
(cdras 'v l
))
2351 (return (mul arg2
(lt-arbpow arg1 index1
)))))
2353 ;; We have specialized the pattern for arbpow1. Now a lot of integrals
2354 ;; will fail correctly and we have to return a noun form.
2355 (return (setq *hyp-return-noun-flag
* 'other-j-cases-next
))))
2357 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2359 ;;; Algorithm 2.1: Laplace transform of c*t^u
2361 ;;; Table of Integral Transforms
2363 ;;; p. 137, formula 1:
2366 ;;; -> gamma(u+1)*p^(-u-1)
2369 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2371 (defun lt-arbpow (expr pow
)
2372 (cond ((or (eq expr
*var
*) (equal pow
0))
2375 (setq *hyp-return-noun-flag
* 'lt-arbow-failed
))))
2377 ;; Check if conditions for f1p137 hold
2378 (defun f1p137test (pow)
2379 (cond ((eq ($asksign
($realpart
(add pow
1))) '$pos
)
2382 (setq *hyp-return-noun-flag
* 'fail-in-arbpow
))))
2385 (mul (take '(%gamma
) (add pow
1))
2386 (power *par
* (sub (mul -
1 pow
) 1))))
2388 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2390 ;;; Algorithm 2.2: Laplace transform of c*t^v*(1+t)^w
2392 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2394 (defun lt-arbpow2 (c a b pow1 pow2
)
2395 (if (eq ($asksign
(add pow1
1)) '$pos
)
2398 ;; The Laplace transform is an Incomplete Gamma function.
2400 (power a
(add pow2
1))
2402 (power (mul *par
* a
(inv b
)) (mul -
1 (add pow2
1)))
2403 (power '$%e
(mul *par
* a
(inv b
)))
2404 (take '(%gamma_incomplete
) (add pow2
1) (mul *par
* a
(inv b
)))))
2405 ((not (maxima-integerp (add pow1 pow2
2)))
2406 ;; The general result is a Hypergeometric U function U(a,b,z) which can
2407 ;; be represented by two Hypergeometic 1F1 functions for the special
2408 ;; case that the index b is not an integer value.
2410 (power a
(add pow1 pow2
1))
2411 (inv (power b
(add pow1
1)))
2412 (take '(%gamma
) (add pow1 pow2
1))
2413 (power (mul *par
* a
(inv b
)) (mul -
1 (add pow1 pow2
1)))
2414 (hgfsimp-exec (list (mul -
1 pow2
))
2415 (list (mul -
1 (add pow1 pow2
)))
2416 (mul *par
* a
(inv b
))))
2418 (power a
(add pow1 pow2
1))
2419 (inv (power b
(add pow1
1)))
2420 (take '(%gamma
) (add pow1
1))
2421 (take '(%gamma
) (mul -
1 (add pow1 pow2
1)))
2422 (inv (take '(%gamma
) (mul -
1 pow2
)))
2423 (hgfsimp-exec (list (add pow1
1))
2424 (list (add pow1 pow2
2))
2425 (mul *par
* a
(inv b
))))))
2427 ;; The most general case is a result with the Hypergeometric U function.
2429 (power a
(add pow1 pow2
1))
2430 (inv (power b
(add pow1
1)))
2431 (take '(%gamma
) (add pow1
1))
2432 (take '($hypergeometric_u
)
2435 (mul *par
* a
(inv b
))))))
2436 (setq *hyp-return-noun-flag
* 'lt-arbpow2-failed
)))
2438 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2440 ;;; Algorithm 2.3: Laplace transform of the Logarithmic function
2442 ;;; c*t^(v-1)*log(a*t)
2443 ;;; -> c*gamma(v)*s^(-v)*(psi[0](v)-log(s/a))
2445 ;;; This is the formula for an expression with log(t) scaled like 1/a*F(s/a).
2447 ;;; For the following cases we have to add further algorithm:
2448 ;;; log(1+a*x), log(x+a), log(x)^2.
2449 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2451 (defun lt-log (rest arg
)
2452 (let* ((l (m2-c*t^v rest
*var
*))
2454 (v (add (cdras 'v l
) 1))) ; because v -> v-1
2456 ((and l
(eq ($asksign v
) '$pos
))
2457 (let* ((l1 (m2-a*t arg
*var
*))
2462 (inv (power *par
* v
))
2463 (sub (take '(mqapply) (list '($psi array
) 0) v
)
2464 (take '(%log
) (div *par
* a
)))))
2466 (setq *hyp-return-noun-flag
* 'lt-log-failed
)))))
2468 (setq *hyp-return-noun-flag
* 'lt-log-failed
)))))
2470 ;; Pattern for lt-log.
2471 ;; Extract the argument of a function: a*t+c for c=0.
2472 (defun m2-a*t
(expr var
)
2475 ((mtimes) (u alike1
,var
) (a free
,var
))
2476 ((coeffpp) (c equal
0)))))
2478 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2479 ;;; Algorithm 2.4: Laplace transform of the Whittaker function
2481 ;;; Test for Whittaker W function. Simplify this if possible, or
2482 ;;; convert to Whittaker M function.
2484 ;;; We have r * %w[i1,i2](a)
2486 ;;; Formula 16, p. 217
2488 ;;; t^(v-1)*%w[k,u](a*t)
2489 ;;; -> gamma(u+v+1/2)*gamma(v-u+1/2)*a^(u+1/2)/
2490 ;;; (gamma(v-k+1)*(p+a/2)^(u+v+1/2)
2491 ;;; *2f1(u+v+1/2,u-k+1/2;v-k+1;(p-a/2)/(p+a/2))
2493 ;;; For Re(v +/- mu) > -1/2
2494 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2496 (defun whittest (r a i1 i2
)
2498 (cond ((f16p217test r a i1 i2
))
2501 (and (maxima-integerp (setq n
(sub (sub '((rat simp
) 1 2) i2
) i1
)))
2502 (member ($sign n
) '($zero $neg $nz
))))
2504 (and (maxima-integerp (setq m
(sub (add '((rat simp
) 1 2) i2
) i1
)))
2505 (member ($sign m
) '($zero $neg $nz
)))))
2506 ;; 1/2-u-k and 1/2+u-k are not zero or a negative integer
2507 ;; Transform to Whittaker M and try again.
2508 (distrexecinit ($expand
(mul (init r
) (wtm a i1 i2
)))))
2510 ;; Both conditions fails, return a noun form.
2511 (setq *hyp-return-noun-flag
* 'whittest-failed
)))))
2513 (defun f16p217test (r a i1 i2
)
2514 ;; We have r*%w[i1,i2](a)
2515 (let ((l (m2-c*t^v r
*var
*)))
2516 ;; Make sure r is of the form c*t^v
2518 (let* ((v (add (cdras 'v l
) 1))
2520 ;; Check that v + i2 + 1/2 > 0 and v - i2 + 1/2 > 0.
2521 (when (and (eq ($asksign
(add (add v i2
) '((rat simp
) 1 2))) '$pos
)
2522 (eq ($asksign
(add (sub v i2
) '((rat simp
) 1 2))) '$pos
))
2523 ;; Ok, we satisfy the conditions. Now extract the arg.
2524 ;; The transformation is only valid for an argument a*t. We have
2525 ;; to special the pattern to make sure that we satisfy the condition.
2526 (let ((l (m2-a*t a
*var
*)))
2528 (let ((a (cdras 'a l
)))
2529 ;; We're ready now to compute the transform.
2531 (power a
(add i2
'((rat simp
) 1 2)))
2532 (take '(%gamma
) (add (add v i2
) '((rat simp
) 1 2)))
2533 (take '(%gamma
) (add (sub v i2
) '((rat simp
) 1 2)))
2534 (inv (mul (take '(%gamma
) (add (sub v i1
) 1))
2535 (power (add *par
* (div a
2))
2536 (add (add i2 v
) '((rat simp
) 1 2)))))
2537 (hgfsimp-exec (list (add (add i2 v
'((rat simp
) 1 2)))
2538 (add (sub i2 i1
) '((rat simp
) 1 2)))
2539 (list (add (sub v i1
) 1))
2540 (div (sub *par
* (div a
2))
2541 (add *par
* (div a
2)))))))))))))
2543 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2544 ;;; Algorithm 2.5: Laplace transform of bessel_k(0,a*t)
2546 ;;; The general algorithm handles the Bessel K function for an order |v|<1.
2547 ;;; but does not include the special case v=0. Return the Laplace transform:
2549 ;;; bessel_k(0,a*t) --> acosh(s/a)/sqrt(s^2-a^2)
2551 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2553 (defun lt-bessel_k0 (rest arg
)
2554 (let* ((l (m2-c*t^v rest
*var
*))
2557 (l (m2-a*t arg
*var
*))
2559 (cond ((and l
(zerop1 v
))
2561 (take '(%acosh
) (div *par
* a
))
2562 (inv (power (sub (mul *par
* *par
*) (mul a a
))
2563 '((rat simp
) 1 2)))))
2565 (setq *hyp-return-noun-flag
* 'lt-bessel_k-failed
)))))
2567 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2569 ;;; DISPATCH FUNCTIONS TO CHANGE THE REPRESENTATION OF SPECIAL FUNCTIONS
2571 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2573 ;; Laplace transform of a product of Bessel functions. A1, A2 are
2574 ;; the args of the two functions. I1, I2 are the indices of each
2575 ;; function. I11, I21 are secondary indices of each function, if any.
2576 ;; FLG is a symbol indicating how we should handle the special
2577 ;; functions (and also indicates what the special functions are.)
2579 ;; I11 and I21 are for the Hankel functions.
2581 (defun fractest (r a1 a2 i1 i11 i2 i21 flg
)
2582 (cond ((or (and (listp i1
) (equal (caar i1
) 'rat
)
2583 (listp i2
) (equal (caar i2
) 'rat
))
2585 ;; We have to Bessel or Hankel functions. Both indizes have to be
2586 ;; rational numbers or we have two Hankel functions.
2588 (cond ((eq flg
'2ytj
)
2592 (mul (htjory i1 i11 a1
)
2593 (htjory i2 i21 a2
)))
2601 (setq *hyp-return-noun-flag
* 'product-of-y-with-nofract-indices
))))
2603 ;; Laplace transform of a product of Bessel functions. A1, A2 are
2604 ;; the args of the two functions. I1, I2 are the indices of each
2605 ;; function. I is a secondary index to one function, if any.
2606 ;; FLG is a symbol indicating how we should handle the special
2607 ;; functions (and also indicates what the special functions are.)
2609 ;; I is for the kind of Hankel function.
2611 (defun fractest1 (r a1 a2 i1 i2 i flg
)
2612 (cond ((or (and (listp i2
)
2613 (equal (caar i2
) 'rat
))
2614 (eq flg
'besshtjory
))
2615 ;; We have two Bessel or Hankel functions. The second index has to
2616 ;; be a rational number or one of the functions is a Hankel function
2617 ;; and the second function is Bessel J or Bessel I
2619 (cond ((eq flg
'bessytj
)
2620 (mul (take '(%bessel_j
) i1 a1
)
2622 ((eq flg
'besshtjory
)
2623 (mul (take '(%bessel_j
) i1 a1
)
2625 ((eq flg
'htjoryytj
)
2626 (mul (htjory i1 i a1
)
2629 (mul (take '(%bessel_j
) i1 a1
)
2631 ((eq flg
'htjorykti
)
2632 (mul (htjory i1 i a1
)
2635 (setq *hyp-return-noun-flag
* 'product-of-i-y-of-nofract-index
))))
2637 ;; Laplace transform of a single special function. A is the arg of
2638 ;; the special function. I1, I11 are the indices of the function. FLG
2639 ;; is a symbol indicating how we should handle the special functions
2640 ;; (and also indicates what the special functions are.)
2642 ;; I11 is the kind of Hankel function
2644 (defun fractest2 (r a1 i1 i11 flg
)
2645 (cond ((or (and (listp i1
)
2646 (equal (caar i1
) 'rat
))
2649 (eq flg
'gamma_incomplete
)
2654 ;; The index is a rational number or flg has the value of one of the
2655 ;; above special functions.
2657 (cond ((eq flg
'ytj
)
2665 ((eq flg
'gamma_incomplete
)
2666 (gamma_incomplete-to-gamma-incomplete-lower a1 i1
))
2672 (slommeltjandy i1 i11 a1
)))))
2674 (setq *hyp-return-noun-flag
* 'y-of-nofract-index
))))
2676 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2678 (defun integertest (r arg i1 i2 flg
)
2679 (cond ((maxima-integerp i1
)
2680 (dispatchpoltrans r arg i1 i2 flg
))
2682 (setq *hyp-return-noun-flag
* 'index-should-be-an-integer-in-polys
))))
2684 (defun dispatchpoltrans (r x i1 i2 flg
)
2686 (cond ((eq flg
'l
)(ltw x i1 i2
))
2687 ((eq flg
'he
)(hetd x i1
))
2688 ((eq flg
'c
)(ctpjac x i1 i2
))
2689 ((eq flg
't
)(ttpjac x i1
))
2690 ((eq flg
'u
)(utpjac x i1
)))))
2692 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2693 ;; (-1)^n*n!*laguerre(n,a,x) = U(-n,a+1,x)
2695 ;; W[k,u](z) = exp(-z/2)*z^(u+1/2)*U(1/2+u-k,1+2*u,z)
2699 ;; laguerre(n,a,x) = (-1)^n*U(-n,a+1,x)/n!
2701 ;; U(-n,a+1,x) = exp(z/2)*z^(-a/2-1/2)*W[1/2+a/2+n,a/2](z)
2705 ;; laguerre(n,a,x) = (-1)^n/n!*exp(z/2)*z^(-a/2-1/2)*M[1/2+a/2+n,a/2](z)
2707 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2710 (let ((diva2 (div a
2)))
2711 (mul (take '(%gamma
) (add n a
1))
2712 (inv (take '(%gamma
) (add a
1)))
2713 (inv (take '(%gamma
) (add n
1)))
2714 (power x
(sub '((rat simp
) -
1 2) diva2
))
2715 (power '$%e
(div x
2))
2716 (list '(mqapply simp
)
2717 (list '($%m simp array
)
2718 (add '((rat simp
) 1 2) diva2 n
) diva2
) x
))))
2720 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2721 ;; Hermite He function as a parabolic cylinder function
2723 ;; Tables of Integral Transforms
2727 ;; D[n](z) = (-1)^n*exp(z^2/4)*diff(exp(-z^2/2),z,n);
2731 ;; He[n](x) = (-1)^n*exp(x^2/2)*diff(exp(-x^2/2),x,n)
2732 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2735 (mul (power '$%e
(mul x x
'((rat simp
) 1 4)))
2736 ;; At this time the Parabolic Cylinder D function is not implemented
2737 ;; as a simplifying function. We call nevertheless the simplifer.
2738 (take '($parabolic_cylinder_d
) n x
)))
2740 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2742 ;; Transform Gegenbauer function to Jacobi P function
2743 ;; ultraspherical(n,v,x) = gamma(2*v+n)*gamma(v+1/2)/gamma(2*v)/gamma(v+n+1/2)
2744 ;; *jacobi_p(n,v-1/2,v-1/2,x)
2745 (defun ctpjac (x n v
)
2746 (let ((inv2 '((rat simp
) 1 2)))
2747 (mul (take '(%gamma
) (add v v n
))
2748 (inv (take '(%gamma
) (add v v
)))
2749 (take '(%gamma
) (add inv2 v
))
2750 (inv (take '(%gamma
) (add v inv2 n
)))
2751 (pjac x n
(sub v inv2
) (sub v inv2
)))))
2753 ;; Transform Chebyshev T function to Jacobi P function
2754 ;; chebyshev_t(n,x) = gamma(n+1)*sqrt(%pi)/gamma(n+1/2)*jacobi_p(n,-1/2,-1/2,x)
2756 (let ((inv2 '((rat simp
) 1 2)))
2757 (mul (take '(%gamma
) n
1)
2758 (power '$%pi inv2
) ; gamma(1/2)
2759 (inv (take '(%gamma
) (add inv2 n
)))
2760 (pjac x n
(mul -
1 inv2
) (mul -
1 inv2
)))))
2762 ;; Transform Chebyshev U function to Jacobi P function
2763 ;; chebyshev_u(n,x) = gamma(n+2)*sqrt(%pi)/2/gamma(3/2+n)*jacobi_p(n,1/2,1/2,x)
2765 (let ((inv2 '((rat simp
) 1 2)))
2766 (mul (take '(%gamma
) (add n
2))
2768 (power '$%pi inv2
) ; gamma(1/2)
2769 (inv (take '(%gamma
) (add inv2 n
1)))
2770 (pjac x n inv2 inv2
))))
2772 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2774 (defun pjactest (rest arg index1 index2 index3
)
2775 (cond ((maxima-integerp index1
)
2779 (list index1 index2 index3
)))
2781 (setq *hyp-return-noun-flag
* 'ind-should-be-an-integer-in-polys
))))
2783 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2785 ;;; Laplace transform of a single Bessel Y function.
2787 ;;; REST is the multiplier, ARG1 is the arg, and INDEX1 is the order of
2788 ;;; the Bessel Y function.
2789 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2791 (defun lt1yref (rest arg1 index1
)
2792 ;; If the index is an integer, use LT1Y. Otherwise, convert Bessel
2793 ;; Y to Bessel J and compute the transform of that. We do this
2794 ;; because converting Y to J for an integer index doesn't work so
2795 ;; well without taking limits.
2796 (cond ((maxima-integerp index1
)
2797 ;; Do not call lt1y but lty directly.
2798 ;; lt1y calls lt-ltp with the flag 'oney. lt-ltp checks this flag
2799 ;; and calls lty. So we can do it at this place and the algorithm is
2801 (lty rest arg1 index1
))
2802 (t (fractest2 rest arg1 index1 nil
'ytj
))))
2804 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2806 ;;; TRANSFORMATIONS TO CHANGE THE REPRESENTATION OF SPECIAL FUNCTIONS
2808 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
2810 ;; erfc in terms of D, parabolic cylinder function
2812 ;; Tables of Integral Transforms
2815 ;; erfc(x) = (%pi*x)^(-1/2)*exp(-x^2/2)*W[-1/4,1/4](x^2)
2818 ;; D[v](z) = 2^(v/2+1/2)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
2822 ;; erfc(x) = %pi^(-1/2)*2^(1/4)*exp(-x^2/2)*D[-1](x*sqrt(2))
2825 (let ((inv2 '((rat simp
) 1 2)))
2826 (mul (power 2 inv2
) ; Should this be 2^(1/4)?
2827 (inv (power '$%pi inv2
))
2828 (power '$%e
(mul -
1 inv2 x x
))
2829 ;; At this time the Parabolic Cylinder D function is not implemented
2830 ;; as a simplifying function. We call nevertheless the simplifer
2831 ;; to simplify the arguments. When we implement the function
2832 ;; The symbol has to be changed to the noun form.
2833 (take '($parabolic_cylinder_d
) -
1 (mul (power 2 inv2
) x
)))))
2835 ;; Lommel S function in terms of Bessel J and Bessel Y.
2838 ;; S[u,v](z) = s[u,v](z) + {2^(u-1)*gamma((u-v+1)/2)*gamma((u+v+1)/2)}
2839 ;; * {sin[(u-v)*%pi/2]*bessel_j(v,z)
2840 ;; - cos[(u-v)*%pi/2]*bessel_y(v,z)
2842 (defun slommeltjandy (m n z
)
2843 (let ((arg (mul '((rat simp
) 1 2) '$%pi
(sub m n
))))
2844 (add (littleslommel m n z
)
2845 (mul (power 2 (sub m
1))
2846 (take '(%gamma
) (div (sub (add m
1) n
) 2))
2847 (take '(%gamma
) (div (add m n
1) 2))
2848 (sub (mul (take '(%sin
) arg
)
2849 (take '(%bessel_j
) n z
))
2850 (mul (take '(%cos
) arg
)
2851 (take '(%bessel_y
) n z
)))))))
2853 ;; Whittaker W function in terms of Whittaker M function
2857 ;; W[k,u](z) = gamma(-2*u)/gamma(1/2-u-k)*M[k,u](z)
2858 ;; + gamma(2*u)/gamma(1/2+u-k)*M[k,-u](z)
2860 (defun wtm (a i1 i2
)
2861 (add (mul (take '(%gamma
) (mul -
2 i2
))
2863 (inv (take '(%gamma
) (sub (sub '((rat simp
) 1 2) i2
) i1
))))
2864 (mul (take '(%gamma
) (add i2 i2
))
2865 (mwhit a i1
(mul -
1 i2
))
2866 (inv (take '(%gamma
) (sub (add '((rat simp
) 1 2) i2
) i1
))))))
2868 ;; Incomplete gamma function in terms of Whittaker W function
2870 ;; Tables of Integral Transforms, p. 387
2872 ;; gamma_incomplete(a,x) = x^((a-1)/2)*exp(-x/2)*W[(a-1)/2,a/2](x)
2874 (defun gammaincompletetw (a x
)
2875 (mul (power x
(div (sub a
1) 2))
2876 (power '$%e
(div x -
2))
2877 (wwhit x
(div (sub a
1) 2)(div a
2))))
2879 ;;; Incomplete Gamma function in terms of lower incomplete Gamma function
2881 ;;; Only for a=0 we use the general representation as a Whittaker W function:
2883 ;;; gamma_incomplete(a,x) = x^((a-1)/2)*exp(-x/2)*W[(a-1)/2,a/2](x)
2885 ;;; In all other cases we transform to a lower incomplete Gamma function:
2887 ;;; gamma_incomplete(a,x) = gamma(a)- gamma_incomplete_lower(a,x)
2889 ;;; The lower incomplete Gamma function will be further transformed to a Hypergeometric 1F1
2890 ;;; representation. With this change we get more simple and correct results for
2891 ;;; the Laplace transform of the Incomplete Gamma function.
2893 (defun gamma_incomplete-to-gamma-incomplete-lower (a x
)
2894 (if (or (eq ($sign a
) '$zero
)
2895 (and (integerp a
) (< a
0)))
2896 ;; The representation as a Whittaker W function for a=0 or a negative
2897 ;; integer (The gamma function is not defined for this case.)
2898 (mul (power x
(div (sub a
1) 2))
2899 (power '$%e
(div x -
2))
2900 (wwhit x
(div (sub a
1) 2) (div a
2)))
2901 ;; In all other cases the representation as a lower incomplete Gamma function
2902 (sub (take '(%gamma
) a
)
2903 (list '(%gamma_incomplete_lower simp
) a x
))))
2905 ;; Bessel Y in terms of Bessel J
2909 ;; bessel_y(v,z) = bessel_j(v,z)*cot(v*%pi)-bessel_j(-v,z)/sin(v*%pi)
2912 (sub (mul (take '(%bessel_j
) i a
)
2913 (take '(%cot
) (mul i
'$%pi
)))
2914 (mul (take '(%bessel_j
) (mul -
1 i
) a
)
2915 (inv (take '(%sin
) (mul i
'$%pi
))))))
2917 ;; Parabolic cylinder function in terms of Whittaker W function.
2919 ;; See Table of Integral Transforms, p.386:
2921 ;; D[v](z) = 2^(v/2+1/4)*z^(-1/2)*W[v/2+1/4,1/4](z^2/2)
2924 (mul (power 2 (add (div i
2) '((rat simp
) 1 4)))
2925 (power a
'((rat simp
) -
1 2))
2926 (wwhit (mul a a
'((rat simp
) 1 2))
2927 (add (div i
2) '((rat simp
) 1 4))
2928 '((rat simp
) 1 4))))
2930 ;; Bateman's function in terms of Whittaker W function
2932 ;; See Table of Integral Transforms, p.386:
2934 ;; k[2*v](z) = 1/gamma(v+1)*W[v,1/2](2*z)
2936 (defun kbatemantw (v a
)
2937 (div (wwhit (add a a
) (div v
2) '((rat simp
) 1 2))
2938 (take '(%gamma
) (add (div v
2) 1))))
2940 ;; Bessel K in terms of Bessel I
2944 ;; bessel_k(v,z) = %pi/2*(bessel_i(-v,z)-bessel_i(v,z))/sin(v*%pi)
2949 (inv (take '(%sin
) (mul i
'$%pi
)))
2950 (sub (take '(%bessel_i
) (mul -
1 i
) a
)
2951 (take '(%bessel_i
) i a
))))
2953 ;; Express Hankel function in terms of Bessel J or Y function.
2957 ;; H[v,1](z) = %i*csc(v*%pi)*(exp(-v*%pi*%i)*bessel_j(v,z) - bessel_j(-v,z))
2960 ;; H[v,2](z) = %i*csc(v*%pi)*(bessel_j(-v,z) - exp(-v*%pi*%i)*bessel_j(v,z))
2962 ;; Both formula are not valid for v an integer.
2963 ;; For this case use the definitions of the Hankel functions:
2964 ;; H[v,1](z) = bessel_j(v,z) + %i* bessel_y(v,z)
2965 ;; H[v,2](z) = bessel_j(v,z) - %i* bessel_y(v,z)
2967 ;; All this can be implemented more simple.
2968 ;; We do not need the transformation to bessel_j for rational numbers,
2969 ;; because the correct transformation for bessel_y is already implemented.
2970 ;; It is enough to use the definitions for the Hankel functions.
2972 (defun htjory (v sort z
)
2973 ;; V is the order, SORT is the kind of Hankel function (1 or 2), Z
2975 (cond ((and (consp v
)
2977 (equal (caar v
) 'rat
))
2978 ;; If the order is a rational number of some sort,
2980 ;; (bessel_j(-v,z) - bessel_j(v,z)*exp(-v*%pi*%i))/(%i*sin(v*%pi*%i))
2981 (div (numjory v sort z
'j
)
2982 (mul '$%i
(take '(%sin
) (mul v
'$%pi
)))))
2984 ;; Transform hankel_1(v,z) to bessel_j(v,z)+%i*bessel_y(v,z)
2985 (add (take '(%bessel_j
) v z
)
2986 (mul '$%i
(take '(%bessel_y
) v z
))))
2988 ;; Transform hankel_2(v,z) to bessel_j(v,z)-%i*bessel_y(v,z)
2989 (sub (take '(%bessel_j
) v z
)
2990 (mul '$%i
(take '(%bessel_y
) v z
))))
2992 ;; We should never reach this point of code.
2993 ;; Problem: The user input for the symbol %h[v,sort](t) is not checked.
2994 ;; Therefore the user can generate this error as long as we do not cut
2995 ;; out the support for the Laplace transform of the symbol %h.
2996 (merror "htjory: Called with wrong sort of Hankel functions."))))
2998 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3000 ;;; Three helper functions only used by htjory
3002 ;; Bessel J or Y, depending on if FLG is 'J or not.
3003 (defun desjy (v z flg
)
3005 (take '(%bessel_j
) v z
))
3007 (take '(%bessel_y
) v z
))))
3009 (defun numjory (v sort z flg
)
3010 (cond ((equal sort
1)
3011 ;; bessel(-v, z) - exp(-v*%pi*%i)*bessel(v, z)
3013 ;; Where bessel is bessel_j if FLG is 'j. Otherwise, bessel
3016 ;; bessel_y(-v, z) - exp(-v*%pi*%i)*bessel_y(v, z)
3017 (sub (desjy (mul -
1 v
) z flg
)
3018 (mul (power '$%e
(mul -
1 v
'$%pi
'$%i
))
3021 ;; exp(-v*%pi*%i)*bessel(v,z) - bessel(-v,z), where bessel is
3022 ;; bessel_j or bessel_y, depending on if FLG is 'j or not.
3023 (sub (mul (power '$%e
(mul v
'$%pi
'$%i
))
3025 (desmjy (mul -
1 v
) z flg
)))))
3027 (defun desmjy (v z flg
)
3030 (take '(%bessel_j
) v z
))
3033 (mul -
1 (take '(%bessel_y
) v z
)))))
3035 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3037 ;;; TRANSFORM TO HYPERGEOMETRIC FUNCTION WITHOUT USING THE ROUTINE REF
3039 ;;; This functions are called in the routine lt-sf-log to get the
3040 ;;; representation in terms of a hypergeometric function.
3041 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3043 ;; The algorithm of the implemented Hermite function %he does not work for
3044 ;; the known Laplace transforms. For an even or odd integer order, we
3045 ;; can represent the Hermite function by the Hypergeometric function 1F1.
3046 ;; With this representations we get the expected Laplace transforms.
3047 (defun hermite-to-hypergeometric (order arg
)
3049 ((and (maxima-integerp order
)
3051 ;; Transform to 1F1 for order an even integer
3052 (mul (power 2 order
)
3053 (power '$%pi
(div 1 2))
3054 (inv (take '(%gamma
) (div (sub 1 order
) 2)))
3055 (list '(mqapply simp
)
3056 (list '($%f array simp
) 1 1)
3057 (list '(mlist) (div order -
2))
3058 (list '(mlist) '((rat simp
) 1 2))
3061 ((and (maxima-integerp order
)
3063 ;; Transform to 1F1 for order an odd integer
3066 (power '$%pi
'((rat simp
) 1 2))
3067 (inv (take '(%gamma
) (div order -
2)))
3068 (list '(mqapply simp
)
3069 (list '($%f simp array
) 1 1)
3070 (list '(mlist simp
) (div (sub 1 order
) 2))
3071 (list '(mlist simp
) '((rat simp
) 3 2))
3074 ;; The general case, transform to 2F0
3075 ;; For this case we have no Laplace transform.
3076 (mul (power (mul 2 arg
) order
)
3077 (list '(mqapply simp
)
3078 (list '($%f array simp
) 2 0)
3079 (list '(mlist simp
) (div order
2) (div (sub 1 order
) 2))
3080 (list '(mlist simp
))
3081 (div -
1 (mul arg arg
)))))))
3083 ;;; Hypergeometric representation of the Exponential Integral Si
3084 ;;; Si(z) = z*1F2([1/2],[3/2,3/2],-z^2/4)
3085 (defun expintegral_si-to-hypergeometric (arg)
3087 (list '(mqapply simp
)
3088 (list '($%f array simp
) 1 2)
3089 (list '(mlist simp
) '((rat simp
) 1 2))
3090 (list '(mlist simp
) '((rat simp
) 3 2) '((rat simp
) 3 2))
3091 (div (mul -
1 arg arg
) 4))))
3093 ;;; Hypergeometric representation of the Exponential Integral Shi
3094 ;;; Shi(z) = z*1F2([1/2],[3/2,3/2],z^2/4)
3095 (defun expintegral_shi-to-hypergeometric (arg)
3097 (list '(mqapply simp
)
3098 (list '($%f simp array
) 1 2)
3099 (list '(mlist simp
) '((rat simp
) 1 2))
3100 (list '(mlist simp
) '((rat simp
) 3 2) '((rat simp
) 3 2))
3101 (div (mul arg arg
) 4))))
3103 ;;; Hypergeometric representation of the Exponential Integral Ci
3104 ;;; Ci(z) = -z^2/4*2F3([1,1],[2,2,3/2],-z^2/4)+log(z)+%gamma
3105 (defun expintegral_ci-to-hypergeometric (arg)
3106 (add (mul (div (mul -
1 arg arg
) 4)
3107 (list '(mqapply simp
)
3108 (list '($%f simp array
) 2 3)
3109 (list '(mlist simp
) 1 1)
3110 (list '(mlist simp
) 2 2 '((rat simp
) 3 2))
3111 (div (mul -
1 arg arg
) 4)))
3115 ;;; Hypergeometric representation of the Exponential Integral Chi
3116 ;;; Chi(z) = z^2/4*2F3([1,1],[2,2,3/2],z^2/4)+log(z)+%gamma
3117 (defun expintegral_chi-to-hypergeometric (arg)
3118 (add (mul (div (mul arg arg
) 4)
3119 (list '(mqapply simp
)
3120 (list '($%f array simp
) 2 3)
3121 (list '(mlist simp
) 1 1)
3122 (list '(mlist simp
) 2 2 '((rat simp
) 3 2))
3123 (div (mul arg arg
) 4)))
3127 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3129 ;;; EXPERTS ON LAPLACE TRANSFORMS
3131 ;;; LT<foo> functions are various experts on Laplace transforms of the
3132 ;;; function <foo>. The expression being transformed is
3133 ;;; r*<foo>(args). The first arg of each expert is r, The remaining
3134 ;;; args are the arg(s) and/or parameters.
3135 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3137 ;; Laplace transform of one Bessel J
3138 (defun lt1j (rest arg index
)
3139 (lt-ltp 'onej rest arg index
))
3141 ;; Laplace transform of two Bessel J functions.
3142 ;; The argument of each must be the same, but the orders may be different.
3143 (defun lt2j (rest arg1 arg2 index1 index2
)
3144 (cond ((not (equal arg1 arg2
))
3145 (setq *hyp-return-noun-flag
* 'product-of-bessel-with-different-args
))
3147 ;; Call lt-ltp two transform and integrate two Bessel J functions.
3148 (lt-ltp 'twoj rest arg1
(list 'list index1 index2
)))))
3150 ;; Laplace transform of a square of a Bessel J function
3151 (defun lt1j^
2 (rest arg index
)
3152 (cond ((alike1 index
'((rat) -
1 2))
3153 ;; Special case: Laplace transform of bessel_j(v,arg)^2, v = -1/2.
3154 ;; For this case the algorithm for the product of two Bessel functions
3155 ;; does not work. Call the integrator with the hypergeometric
3156 ;; representation: 2/%pi/arg*1/2*(1+0F1([], [1/2], -arg*arg)).
3157 (sendexec (mul (div 2 '$%pi
)
3160 (add '((rat simp
) 1 2)
3161 (mul '((rat simp
) 1 2)
3162 (list '(mqapply simp
)
3163 (list '($%f simp array
) 1 0)
3164 (list '(mlist simp
) )
3165 (list '(mlist simp
) '((rat simp
) 1 2))
3166 (mul -
1 (mul arg arg
)))))))
3168 (lt-ltp 'twoj rest arg
(list 'list index index
)))))
3170 ;; Laplace transform of Incomplete Gamma function
3171 (defun lt1gamma-incomplete-lower (rest arg1 arg2
)
3172 (lt-ltp 'gamma-incomplete-lower rest arg2 arg1
))
3174 ;; Laplace transform of Whittaker M function
3175 (defun lt1m (r a i1 i2
)
3176 (lt-ltp 'onem r a
(list i1 i2
)))
3178 ;; Laplace transform of Jacobi function
3179 (defun lt1p (r a i1 i2
)
3180 (lt-ltp 'hyp-onep r a
(list i1 i2
)))
3182 ;; Laplace transform of Associated Legendre function of the second kind
3183 (defun lt1q (r a i1 i2
)
3184 (lt-ltp 'oneq r a
(list i1 i2
)))
3186 ;; Laplace transform of Erf function
3187 (defun lt1erf (rest arg
)
3188 (lt-ltp 'onerf rest arg nil
))
3190 ;; Laplace transform of Log function
3191 (defun lt1log (rest arg
)
3192 (lt-ltp 'onelog rest arg nil
))
3194 ;; Laplace transform of Complete elliptic integral of the first kind
3195 (defun lt1kelliptic (rest arg
)
3196 (lt-ltp 'onekelliptic rest arg nil
))
3198 ;; Laplace transform of Complete elliptic integral of the second kind
3199 (defun lt1e (rest arg
)
3200 (lt-ltp 'onee rest arg nil
))
3202 ;; Laplace transform of Struve H function
3203 (defun lt1hstruve (rest arg1 index1
)
3204 (lt-ltp 'hs rest arg1 index1
))
3206 ;; Laplace transform of Struve L function
3207 (defun lt1lstruve (rest arg1 index1
)
3208 (lt-ltp 'hl rest arg1 index1
))
3210 ;; Laplace transform of Lommel s function
3211 (defun lt1s (rest arg1 index1 index2
)
3212 (lt-ltp 's rest arg1
(list index1 index2
)))
3214 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3216 ;;; TRANSFORM TO HYPERGEOMETRIC FUNCTION AND DO THE INTEGRATION
3218 ;;; FLG = special function we're transforming
3219 ;;; REST = other stuff
3220 ;;; ARG = arg of special function
3221 ;;; INDEX = index of special function.
3223 ;;; So we're transforming REST*FLG(INDEX, ARG).
3224 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3226 (defun lt-ltp (flg rest arg index
)
3228 (when (and (listp index
)
3229 (eq (car index
) 'list
))
3230 (setq index
(list (cadr index
) (caddr index
))))
3234 (car (setq l1
(ref flg index arg
)))))
3236 ;; Convert the special function to a hypgergeometric
3237 ;; function. L1 is the special function converted to the
3238 ;; hypergeometric function. d*x^m*%e^a*x looks for that
3239 ;; factor in the expanded form.
3242 ;; We currently don't know how to handle this yet.
3243 ;; We add the return of a noun form.
3244 (setq *hyp-return-noun-flag
* 'other-ca-later
)))))
3246 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3247 ;; Dispatch function to convert the given function to a hypergeometric
3250 ;; The first arg is a symbol naming the function; the last is the
3251 ;; argument to the function. The second arg is the index (or list of
3252 ;; indices) to the function. Not used if the function doesn't have
3255 ;; The result is a list of 2 elements: The first element is a
3256 ;; multiplier; the second, the hypergeometric function itself.
3257 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3259 (defun ref (flg index arg
)
3261 (onej (j1tf index arg
))
3262 (twoj (j2tf (car index
) (cadr index
) arg
))
3263 (hs (hstf index arg
))
3264 (hl (lstf index arg
))
3265 (s (stf (car index
) (cadr index
) arg
))
3267 (onelog (logtf arg
))
3268 (onekelliptic (kelliptictf arg
)) ; elliptic_kc
3269 (onee (etf arg
)) ; elliptic_ec
3270 (onem (mtf (car index
) (cadr index
) arg
))
3271 (hyp-onep (ptf (car index
) (cadr index
) arg
))
3272 (oneq (qtf (car index
) (cadr index
) arg
))
3273 (gamma-incomplete-lower (gamma-incomplete-lower-tf index arg
))
3274 (onepjac (pjactf (car index
) (cadr index
) (caddr index
) arg
))
3278 ;; Transform %f to internal representation FPQ
3279 (list 1 (ref-fpq (rest (car index
)) (rest (cadr index
)) arg
)))))
3281 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3283 ;;; TRANSFORM FUNCTION IN TERMS OF HYPERGEOMETRIC FUNCTION
3285 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3287 ;; Create a hypergeometric form that we recognize. The function is
3288 ;; %f[n,m](p; q; arg). We represent this as a list of the form
3289 ;; (fpq (<length p> <length q>) <p> <q> <arg>)
3291 (defun ref-fpq (p q arg
)
3292 (list 'fpq
(list (length p
) (length q
))
3295 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3297 ;; Whittaker M function in terms of hypergeometric function
3301 ;; M[k,u](z) = exp(-z/2)*z^(1/2+u)*M(1/2+u-k,1+2*u,z)
3303 (defun mtf (i1 i2 arg
)
3304 (list (mul (power arg
(add i2
'((rat simp
) 1 2)))
3305 (power '$%e
(div arg -
2)))
3306 (ref-fpq (list (add '((rat simp
) 1 2) i2
(mul -
1 i1
)))
3307 (list (add i2 i2
1))
3310 ;; Jacobi P in terms of hypergeometric function
3314 ;; F(-n,a+1+b+n; a+1; x) = n!/poch(a+1,n)*jacobi_p(n,a,b,1-2*x)
3316 ;; jacobi_p(n,a,b,x) = poch(a+1,n)/n!*F(-n,a+1+b+n; a+1; (1-x)/2)
3317 ;; = gamma(a+n+1)/gamma(a+1)/n!*F(-n,a+1+b+n; a+1; (1-x)/2)
3319 ;; We have a problem:
3320 ;; We transform the argument x -> (1-x)/2. But an argument with a constant
3321 ;; part is not integrable with the implemented algorithm for the hypergeometric
3322 ;; function. It might be possible to get a result for an argument y=1-2*x
3323 ;; with x=t^-2. But for this case the routine lt-ltp fails. The routine
3324 ;; recognize a constant term in the argument, but does not take into account
3325 ;; that the constant term might vanish, when we transform to a hypergeometric
3327 ;; Because of this the Laplace transform for the following functions does not
3328 ;; work too: Legendre P, Chebyshev T, Chebyshev U, and Gegenbauer.
3330 (defun pjactf (n a b x
)
3331 (list (mul (take '(%gamma
) (add n a
1))
3332 (inv (take '(%gamma
) (add a
1)))
3333 (inv (take '(%gamma
) (add n
1))))
3334 (ref-fpq (list (mul -
1 n
) (add n a b
1))
3336 (sub '((rat simp
) 1 2) (div x
2)))))
3338 ;; ArcSin in terms of hypergeometric function
3342 ;; F(1/2,1/2; 3/2; z^2) = asin(z)/z
3344 ;; asin(z) = z*F(1/2,1/2; 3/2; z^2)
3347 (let ((inv2 '((rat simp
) 1 2)))
3349 (ref-fpq (list inv2 inv2
)
3350 (list '((rat simp
) 3 2))
3353 ;; ArcTan in terms of hypergeometric function
3357 ;; F(1/2,1; 3/2; -z^2) = atan(z)/z
3359 ;; atan(z) = z*F(1/2,1; 3/2; -z^2)
3363 (ref-fpq (list '((rat simp
) 1 2) 1)
3364 (list '((rat simp
) 3 2))
3367 ;; Associated Legendre function P in terms of hypergeometric function
3371 ;; assoc_legendre_p(v,u,z) = ((z+1)/(z-2))^(u/2)/gamma(1-u)*F(-v,v+1;1-u,(1-z)/2)
3373 ;; FIXME: What about the branch cut? 8.1.2 is for z not on the real
3374 ;; line with -1 < z < 1.
3377 (list (mul (inv (take '(%gamma
) (sub 1 m
)))
3378 (power (div (add z
1)
3381 (ref-fpq (list (mul -
1 n
) (add n
1))
3383 (sub '((rat simp
) 1 2) (div z
2)))))
3385 ;; Associated Legendre function Q in terms of hypergeometric function
3389 ;; assoc_legendre_q(v,u,z)
3390 ;; = exp(%i*u*%pi)*2^(-v-1)*sqrt(%pi) *
3391 ;; gamma(v+u+1)/gamma(v+3/2)*z^(-v-u-1)*(z^2-1)^(u/2) *
3392 ;; F(1+v/2+u/2, 1/2+v/2+u/2; v+3/2; 1/z^2)
3394 ;; FIXME: What about the branch cut?
3397 (list (mul (power '$%e
(mul m
'$%pi
'$%i
))
3398 (power '$%pi
'((rat simp
) 1 2))
3399 (take '(%gamma
) (add m n
1))
3400 (power 2 (sub -
1 n
))
3401 (inv (take '(%gamma
) (add n
'((rat simp
) 3 2))))
3402 (power z
(mul -
1 (add m n
1)))
3403 (power (sub (mul z z
) 1) (div m
2)))
3404 (ref-fpq (list (div (add m n
1) 2) (div (add m n
2) 2))
3405 (list (add n
'((rat simp
) 3 2)))
3408 ;; Incomplete lower Gamma in terms of hypergeometric function
3412 ;; M(a,a+1,-x) = a*x^(-a)*gamma_incomplete_lower(a,x)
3414 ;; gamma_incomplete_lower(a,x) = x^a/a*M(a,a+1,-x)
3416 (defun gamma-incomplete-lower-tf (a x
)
3417 (list (mul (inv a
) (power x a
))
3422 ;; Complete elliptic K in terms of hypergeometric function
3426 ;; K(k) = %pi/2*2F1(1/2,1/2; 1; k^2)
3428 (defun kelliptictf (k)
3429 (let ((inv2 '((rat simp
) 1 2)))
3430 (list (mul inv2
'$%pi
)
3431 (ref-fpq (list inv2 inv2
)
3435 ;; Complete elliptic E in terms of hypergeometric function
3439 ;; E(k) = %pi/2*2F1(-1/2,1/2;1;k^2)
3442 (let ((inv2 '((rat simp
) 1 2)))
3443 (list (mul inv2
'$%pi
)
3444 (ref-fpq (list (mul -
1 inv2
) inv2
)
3448 ;; erf in terms of hypgeometric function.
3452 ;; erf(z) = 2*z/sqrt(%pi)*M(1/2,3/2,-z^2)
3453 ;; = 2*z/sqrt(%pi)*exp(-z^2)*M(1,3/2,z^2)
3456 (list (mul 2 arg
(power '$%pi
'((rat simp
) -
1 2)))
3457 (ref-fpq (list '((rat simp
) 1 2))
3458 (list '((rat simp
) 3 2))
3461 ;; log in terms of hypergeometric function
3463 ;; We know from A&S 15.1.3 that
3465 ;; F(1,1;2;z) = -log(1-z)/z.
3467 ;; So log(z) = (z-1)*F(1,1;2;1-z)
3475 ;; Bessel J function expressed as a hypergeometric function.
3479 ;; bessel_j(v,z) = (z/2)^v*sum (-z^2/4)^k/k!/gamma(v+k+1)
3482 ;; = (z/2)^v/gamma(v+1)*sum 1/poch(v+1,k)*(-z^2/4)^k/k!
3484 ;; = (z/2)^v/gamma(v+1) * 0F1(; v+1; -z^2/4)
3487 (list (mul (inv (power 2 v
))
3489 (inv (take '(%gamma
) (add v
1))))
3492 (mul '((rat simp
) -
1 4) (power z
2)))))
3494 ;; Product of 2 Bessel J functions in terms of hypergeometric function
3496 ;; See Y. L. Luke, formula 39, page 216:
3498 ;; bessel_j(u,z)*bessel_j(v,z)
3499 ;; = (z/2)^(u+v)/gamma(u+1)/gamma(v+1) *
3500 ;; 2F3((u+v+1)/2, (u+v+2)/2; u+1, v+1, u+v+1; -z^2)
3502 (defun j2tf (n m arg
)
3503 (list (mul (inv (take '(%gamma
) (add n
1)))
3504 (inv (take '(%gamma
) (add m
1)))
3505 (inv (power 2 (add n m
)))
3506 (power arg
(add n m
)))
3507 (ref-fpq (list (add '((rat simp
) 1 2) (div n
2) (div m
2))
3508 (add 1 (div n
2) (div m
2)))
3509 (list (add 1 n
) (add 1 m
) (add 1 n m
))
3510 (mul -
1 (power arg
2)))))
3512 ;; Struve H function in terms of hypergeometric function.
3514 ;; A&S 12.1.2 gives the following series for the Struve H function:
3517 ;; H[v](z) = (z/2)^(v+1)*sum (-1)^k*(z/2)^(2*k)/gamma(k+3/2)/gamma(k+v+3/2)
3520 ;; We can write this in the form
3522 ;; H[v](z) = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2)
3525 ;; * sum n!/poch(3/2,n)/poch(v+3/2,n)*(-z^2/4)^n/n!
3528 ;; = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2) * 1F2(1;3/2,v+3/2;(-z^2/4))
3530 ;; See also A&S 12.1.21.
3533 (let ((d32 (div 3 2)))
3534 (list (mul (power (div z
2)(add v
1))
3535 (inv (take '(%gamma
) d32
))
3536 (inv (take '(%gamma
) (add v d32
))))
3538 (list d32
(add v d32
))
3539 (mul '((rat simp
) -
1 4) z z
)))))
3541 ;; Struve L function in terms of hypergeometric function
3545 ;; L[v](z) = -%i*exp(-v*%i*%pi/2)*H[v](%i*z)
3547 ;; This function computes exactly this way. (But why is %i written as
3548 ;; exp(%i*%pi/2) instead of just %i)
3550 ;; A&S 12.2.1 gives the series expansion as
3553 ;; L[v](z) = (z/2)^(v+1)*sum (z/2)^(2*k)/gamma(k+3/2)/gamma(k+v+3/2)
3556 ;; It's quite easy to derive
3558 ;; L[v](z) = 2/sqrt(%pi)*(z/2)^(v+1)/gamma(v+3/2) * 1F2(1;3/2,v+3/2;(z^2/4))
3561 (let ((d32 (div 3 2)))
3562 (list (mul (power (div z
2) (add v
1))
3563 (inv (take '(%gamma
) d32
))
3564 (inv (take '(%gamma
) (add v d32
))))
3566 (list d32
(add v d32
))
3567 (mul '((rat simp
) 1 4) z z
)))))
3569 ;; Lommel s function in terms of hypergeometric function
3571 ;; See Y. L. Luke, p 217, formula 1
3573 ;; 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)
3576 (list (mul (power z
(add m
1))
3577 (inv (sub (add m
1) n
))
3580 (list (div (sub (add m
3) n
) 2)
3581 (div (add m n
3) 2))
3582 (mul '((rat simp
) -
1 4) z z
))))
3584 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3586 ;;; Algorithm 3: Laplace transform of a hypergeometric function
3588 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3590 (defun %$etest
(l l1
)
3592 (setq q
(cdras 'q l
))
3593 (cond ((equal q
1)(setq a
0)(go loop
)))
3594 (setq a
(cdras 'a l
))
3596 (return (substl (sub *par
* a
)
3598 (execf19 l
(cadr l1
))))))
3600 (defun execf19 (l1 l2
)
3602 (setq ans
(execargmatch (car (cddddr l2
))))
3603 (cond ((eq (car ans
) 'dionimo
)
3604 (return (dionarghyp l1 l2
(cadr ans
)))))
3605 ;; We specialized the pattern for the argument. Now the test fails
3606 ;; correctly and we have to add the return of a noun form.
3608 (setq *hyp-return-noun-flag
* 'next-for-other-args
))))
3610 ;; Executive for recognizing the sort of argument to the
3611 ;; hypergeometric function. We look to see if the arg is of the form
3612 ;; a*x^m + c. Return a list of 'dionimo (what does that mean?) and
3615 (defun execargmatch (arg)
3617 (cond ((setq l1
(m2-a*x^m
+c
($factor arg
) *var
*))
3618 (return (list 'dionimo l1
))))
3619 (cond ((setq l1
(m2-a*x^m
+c
($expand arg
) *var
*))
3620 (return (list 'dionimo l1
))))
3621 ;; The return value has to be a list.
3622 (return (list 'other-case-args-to-follow
))))
3624 ;; We have hypergeometric function whose arg looks like a*x^m+c. L1
3625 ;; matches the d*x^m... part, L2 is the hypergeometric function and
3626 ;; arg is the match for a*x^m+c.
3628 (defun dionarghyp (l1 l2 arg
)
3630 (setq a
(cdras 'a arg
)
3634 ((and (integerp m
) ; The power of the argument has to be
3635 (plusp m
) ; a positive integer.
3636 (equal 0 c
)) ; No const term to the argument.
3637 (return (f19cond a m l1 l2
)))
3640 (setq *hyp-return-noun-flag
* 'prop4-and-other-cases-to-follow
))))))
3642 (defun f19cond (a m l1 l2
)
3650 ;; At this point, we have the function d*x^s*%f[p,q](l1, l2, (a*t)^m).
3651 ;; Check to see if Formula 19, p 220 applies.
3652 (cond ((and (not (eq ($asksign
(sub (add p m -
1) q
)) '$pos
))
3653 (eq ($asksign
(add s
1)) '$pos
))
3655 (f19p220-simp (add s
1)
3660 (return (setq *hyp-return-noun-flag
*
3661 'failed-on-f19cond-multiply-the-other-cases-with-d
))))
3663 ;; Table of Laplace transforms, p 220, formula 19:
3665 ;; If m + k <= n + 1, and Re(s) > 0, the Laplace transform of
3667 ;; t^(s-1)*%f[m,n]([a1,...,am],[p1,...,pn],(c*t)^k)
3670 ;; 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)
3672 ;; with Re(p) > 0 if m + k <= n, Re(p+k*c*exp(2*%pi*%i*r/k)) > 0 for r
3673 ;; = 0, 1,...,k-1, if m + k = n + 1.
3675 ;; The args below are s, [a's], [p's], c^k, k.
3677 (defun f19p220-simp (s l1 l2 cf k
)
3678 (mul (take '(%gamma
) s
)
3679 (inv (power *par
* s
))
3680 (hgfsimp-exec (append l1
(addarglist s k
))
3684 (power (inv *par
*) k
)))))
3686 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3688 ;; Return a list of s/k, (s+1)/k, ..., (s+|k|-1)/k
3689 (defun addarglist (s k
)
3690 (let ((abs-k (abs k
))
3693 (push (div (add s n
) k
) res
))
3696 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3698 ;;; Pattern for the Laplace transform of the hypergeometric function
3700 ;; Match d*x^m*%e^(a*x). If we match, Q is the e^(a*x) part, A is a,
3701 ;; M is M, and D is d.
3702 (defun m2-d*x^m
*%e^a
*x
(expr var par
)
3705 ((coefftt) (d free2
,var
,par
))
3706 ((mexpt) (x alike1
,var
) (m free2
,var
,par
))
3710 ((coefftt) (a free2
,var
,par
))
3711 (x alike1
,var
))))))
3714 (defun m2-f+c
(expr var
)
3717 ((coeffpt) (f has
,var
))
3718 ((coeffpp) (c free
,var
)))))
3721 ;; The pattern was too general. We match also a*t^2+b*t. But that's not correct.
3722 (defun m2-a*x^m
+c
(expr var
)
3725 ((coefft) ; more special (not coeffpt)
3727 ((mexpt) (x alike1
,var
) (m free-not-zero-p
,var
)))
3728 ((coeffpp) (c free
,var
)))))
3730 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3732 ;;; Algorithm 4: SPECIAL HANDLING OF Bessel Y for an integer order
3734 ;;; This is called for one Bessel Y function, when the order is an integer.
3735 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3737 (defun lty (rest arg index
)
3739 (cond ((setq l
(m2-d*x^m
*%e^a
*x rest
*var
* *par
*))
3740 (return (execfy l arg index
))))
3741 (return (setq *hyp-return-noun-flag
* 'fail-in-lty
))))
3743 (defun execfy (l arg index
)
3745 (setq ans
(execargmatch arg
))
3746 (cond ((eq (car ans
) 'dionimo
)
3747 (return (dionarghyp-y l index
(cadr ans
)))))
3748 (return (setq *hyp-return-noun-flag
* 'fail-in-execfy
))))
3750 (defun dionarghyp-y (l index arg
)
3752 (setq a
(cdras 'a arg
)
3755 (cond ((and (zerp c
) (equal m
1))
3756 (let ((ans (f2p105v2cond a l index
)))
3757 (unless (symbolp ans
)
3759 (cond ((and (zerp c
) (equal m
'((rat simp
) 1 2)))
3760 (let ((ans (f50cond a l index
)))
3761 (unless (symbolp ans
)
3763 (return (setq *hyp-return-noun-flag
* 'fail-in-dionarghyp-y
))))
3765 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3767 ;;; Algorithm 4.1: Laplace transform of t^n*bessel_y(v,a*t)
3768 ;;; v is an integer and n>=v
3770 ;;; Table of Integral Transforms
3772 ;;; Volume 2, p 105, formula 2 is a formula for the Y-transform of
3774 ;;; f(x) = x^(u-3/2)*exp(-a*x)
3776 ;;; where the Y-transform is defined by
3778 ;;; integrate(f(x)*bessel_y(v,x*y)*sqrt(x*y), x, 0, inf)
3782 ;;; -2/%pi*gamma(u+v)*sqrt(y)*(y^2+a^2)^(-u/2)
3783 ;;; *assoc_legendre_q(u-1,-v,a/sqrt(y^2+a^2))
3785 ;;; with a > 0, Re u > |Re v|.
3787 ;;; In particular, with a slight change of notation, we have
3789 ;;; integrate(x^(u-1)*exp(-p*x)*bessel_y(v,a*x)*sqrt(a), x, 0, inf)
3791 ;;; which is the Laplace transform of x^(u-1/2)*bessel_y(v,x).
3793 ;;; Thus, the Laplace transform is
3795 ;;; -2/%pi*gamma(u+v)*sqrt(a)*(a^2+p^2)^(-u/2)
3796 ;;; *assoc_legendre_q(u-1,-v,p/sqrt(a^2+p^2))
3797 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3799 (defun f2p105v2cond (a l index
)
3801 (setq d
(cdras 'd l
) ; contains constant part of integrand
3804 (cond ((eq ($asksign
($realpart
(sub m index
))) '$pos
)
3805 (return (mul d
(f2p105v2cond-simp m index a
)))))
3806 (return (setq *hyp-return-noun-flag
* 'fail-in-f2p105v2cond
))))
3808 (defun f2p105v2cond-simp (m v a
)
3811 (take '(%gamma
) (add m v
))
3812 (power (add (mul a a
) (mul *par
* *par
*))
3813 (mul -
1 '((rat simp
) 1 2) m
))
3814 ;; Call Associated Legendre Q function, which simplifies accordingly.
3815 ;; We have to do a Maxima function call, because $assoc_legendre_q is
3816 ;; not in Maxima core and has to be autoloaded.
3817 (mfuncall '$assoc_legendre_q
3821 (power (add (mul a a
) (mul *par
* *par
*))
3822 '((rat simp
) -
1 2))))))
3824 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3826 ;;; Algorithm 4.2: Laplace transform of t^n*bessel_y(v, a*sqrt(t))
3828 ;;; Table of Integral Transforms
3830 ;;; p. 188, formula 50:
3832 ;;; t^(u-1/2)*bessel_y(2*v,2*sqrt(a)*sqrt(t))
3833 ;;; -> a^(-1/2)*p^(-u)*exp(-a/2/p)
3834 ;;; * [tan((u-v)*%pi)*gamma(u+v+1/2)/gamma(2*v+1)*M[u,v](a/p)
3835 ;;; -sec((u-v)*%pi)*W[u,v](a/p)]
3836 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
3838 (defun f50cond (a l v
)
3840 (setq d
(cdras 'd l
)
3841 m
(add (cdras 'm l
) '((rat simp
) 1 2))
3844 ((and (eq ($asksign
($realpart
(add m v
'((rat simp
) 1 2)))) '$pos
)
3845 (eq ($asksign
($realpart
(sub (add m
'((rat simp
) 1 2)) v
))) '$pos
)
3846 (not (maxima-integerp (mul (sub (add m m
) (add v v
1))
3847 '((rat simp
) 1 2)))))
3848 (setq a
(mul a a
'((rat simp
) 1 4)))
3849 (return (f50p188-simp d m v a
))))
3850 (return (setq *hyp-return-noun-flag
* 'fail-in-f50cond
))))
3852 (defun f50p188-simp (d u v a
)
3854 (power a
'((rat simp
) -
1 2))
3855 (power *par
* (mul -
1 u
))
3856 (power '$%e
(div a
(mul -
2 *par
*)))
3857 (sub (mul (take '(%tan
) (mul '$%pi
(sub u v
)))
3858 (take '(%gamma
) (add u v
'((rat simp
) 1 2)))
3859 (inv (take '(%gamma
) (add v v
1)))
3860 (mwhit (div a
*par
*) u v
))
3861 (mul (take '(%sec
) (mul '$%pi
(sub u v
)))
3862 (wwhit (div a
*par
*) u v
)))))
3864 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;