-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathspecial.f90
780 lines (580 loc) · 15.9 KB
/
special.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
MODULE special
IMPLICIT NONE
CONTAINS
!These functions come from:
!http://people.sc.fsu.edu/~jburkardt/f_src/special_functions/special_functions.html
subroutine hygfx ( a, b, c, x, hf )
!*****************************************************************************80
!
!! HYGFX evaluates the hypergeometric function F(A,B,C,X).
!
! Licensing:
!
! The original FORTRAN77 version of this routine is copyrighted by
! Shanjie Zhang and Jianming Jin. However, they give permission to
! incorporate this routine into a user program that the copyright
! is acknowledged.
!
! Modified:
!
! 08 September 2007
!
! Author:
!
! Original FORTRAN77 version by Shanjie Zhang, Jianming Jin.
! FORTRAN90 version by John Burkardt.
!
! Reference:
!
! Shanjie Zhang, Jianming Jin,
! Computation of Special Functions,
! Wiley, 1996,
! ISBN: 0-471-11963-6,
! LC: QA351.C45
!
! Parameters:
!
! Input, real ( kind = 8 ) A, B, C, X, the arguments of the function.
! C must not be equal to a nonpositive integer.
! X < 1.
!
! Output, real HF, the value of the function.
!
implicit none
real ( kind = 8 ) a
real ( kind = 8 ) a0
real ( kind = 8 ) aa
real ( kind = 8 ) b
real ( kind = 8 ) bb
real ( kind = 8 ) c
real ( kind = 8 ) c0
real ( kind = 8 ) c1
real ( kind = 8 ), parameter :: el = 0.5772156649015329D+00
real ( kind = 8 ) eps
real ( kind = 8 ) f0
real ( kind = 8 ) f1
real ( kind = 8 ) g0
real ( kind = 8 ) g1
real ( kind = 8 ) g2
real ( kind = 8 ) g3
real ( kind = 8 ) ga
real ( kind = 8 ) gabc
real ( kind = 8 ) gam
real ( kind = 8 ) gb
real ( kind = 8 ) gbm
real ( kind = 8 ) gc
real ( kind = 8 ) gca
real ( kind = 8 ) gcab
real ( kind = 8 ) gcb
real ( kind = 8 ) gm
real ( kind = 8 ) hf
real ( kind = 8 ) hw
integer ( kind = 4 ) j
integer ( kind = 4 ) k
logical l0
logical l1
logical l2
logical l3
logical l4
logical l5
integer ( kind = 4 ) m
integer ( kind = 4 ) nm
real ( kind = 8 ) pa
real ( kind = 8 ) pb
real ( kind = 8 ), parameter :: pi = 3.141592653589793D+00
real ( kind = 8 ) r
real ( kind = 8 ) r0
real ( kind = 8 ) r1
real ( kind = 8 ) rm
real ( kind = 8 ) rp
real ( kind = 8 ) sm
real ( kind = 8 ) sp
real ( kind = 8 ) sp0
real ( kind = 8 ) x
real ( kind = 8 ) x1
l0 = ( c == aint ( c ) ) .and. ( c < 0.0D+00 )
l1 = ( 1.0D+00 - x < 1.0D-15 ) .and. ( c - a - b <= 0.0D+00 )
l2 = ( a == aint ( a ) ) .and. ( a < 0.0D+00 )
l3 = ( b == aint ( b ) ) .and. ( b < 0.0D+00 )
l4 = ( c - a == aint ( c - a ) ) .and. ( c - a <= 0.0D+00 )
l5 = ( c - b == aint ( c - b ) ) .and. ( c - b <= 0.0D+00 )
if ( l0 .or. l1 ) then
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'HYGFX - Fatal error!'
write ( *, '(a)' ) ' The hypergeometric series is divergent.'
return
end if
if ( 0.95D+00 < x ) then
eps = 1.0D-08
else
eps = 1.0D-15
end if
if ( x == 0.0D+00 .or. a == 0.0D+00 .or. b == 0.0D+00 ) then
hf = 1.0D+00
return
else if ( 1.0D+00 - x == eps .and. 0.0D+00 < c - a - b ) then
call gamma_f ( c, gc )
call gamma_f ( c - a - b, gcab )
call gamma_f ( c - a, gca )
call gamma_f ( c - b, gcb )
hf = gc * gcab /( gca *gcb )
return
else if ( 1.0D+00 + x <= eps .and. abs ( c - a + b - 1.0D+00 ) <= eps ) then
g0 = sqrt ( pi ) * 2.0D+00**( - a )
call gamma_f ( c, g1 )
call gamma_f ( 1.0D+00 + a / 2.0D+00 - b, g2 )
call gamma_f ( 0.5D+00 + 0.5D+00 * a, g3 )
hf = g0 * g1 / ( g2 * g3 )
return
else if ( l2 .or. l3 ) then
if ( l2 ) then
nm = int ( abs ( a ) )
end if
if ( l3 ) then
nm = int ( abs ( b ) )
end if
hf = 1.0D+00
r = 1.0D+00
do k = 1, nm
r = r * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) &
/ ( k * ( c + k - 1.0D+00 ) ) * x
hf = hf + r
end do
return
else if ( l4 .or. l5 ) then
if ( l4 ) then
nm = int ( abs ( c - a ) )
end if
if ( l5 ) then
nm = int ( abs ( c - b ) )
end if
hf = 1.0D+00
r = 1.0D+00
do k = 1, nm
r = r * ( c - a + k - 1.0D+00 ) * ( c - b + k - 1.0D+00 ) &
/ ( k * ( c + k - 1.0D+00 ) ) * x
hf = hf + r
end do
hf = ( 1.0D+00 - x )**( c - a - b ) * hf
return
end if
aa = a
bb = b
x1 = x
!
! WARNING: ALTERATION OF INPUT ARGUMENTS A AND B, WHICH MIGHT BE CONSTANTS.
!
if ( x < 0.0D+00 ) then
x = x / ( x - 1.0D+00 )
if ( a < c .and. b < a .and. 0.0D+00 < b ) then
a = bb
b = aa
end if
b = c - b
end if
if ( 0.75D+00 <= x ) then
gm = 0.0D+00
if ( abs ( c - a - b - aint ( c - a - b ) ) < 1.0D-15 ) then
m = int ( c - a - b )
call gamma_f ( a, ga )
call gamma_f ( b, gb )
call gamma_f ( c, gc )
call gamma_f ( a + m, gam )
call gamma_f ( b + m, gbm )
call psi ( a, pa )
call psi ( b, pb )
if ( m /= 0 ) then
gm = 1.0D+00
end if
do j = 1, abs ( m ) - 1
gm = gm * j
end do
rm = 1.0D+00
do j = 1, abs ( m )
rm = rm * j
end do
f0 = 1.0D+00
r0 = 1.0D+00
r1 = 1.0D+00
sp0 = 0.0D+00
sp = 0.0D+00
if ( 0 <= m ) then
c0 = gm * gc / ( gam * gbm )
c1 = - gc * ( x - 1.0D+00 )**m / ( ga * gb * rm )
do k = 1, m - 1
r0 = r0 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) &
/ ( k * ( k - m ) ) * ( 1.0D+00 - x )
f0 = f0 + r0
end do
do k = 1, m
sp0 = sp0 + 1.0D+00 / ( a + k - 1.0D+00 ) &
+ 1.0D+00 / ( b + k - 1.0D+00 ) - 1.0D+00 / real ( k, kind = 8 )
end do
f1 = pa + pb + sp0 + 2.0D+00 * el + log ( 1.0D+00 - x )
hw = f1
do k = 1, 250
sp = sp + ( 1.0D+00 - a ) / ( k * ( a + k - 1.0D+00 ) ) &
+ ( 1.0D+00 - b ) / ( k * ( b + k - 1.0D+00 ) )
sm = 0.0D+00
do j = 1, m
sm = sm + ( 1.0D+00 - a ) &
/ ( ( j + k ) * ( a + j + k - 1.0D+00 ) ) &
+ 1.0D+00 / ( b + j + k - 1.0D+00 )
end do
rp = pa + pb + 2.0D+00 * el + sp + sm + log ( 1.0D+00 - x )
r1 = r1 * ( a + m + k - 1.0D+00 ) * ( b + m + k - 1.0D+00 ) &
/ ( k * ( m + k ) ) * ( 1.0D+00 - x )
f1 = f1 + r1 * rp
if ( abs ( f1 - hw ) < abs ( f1 ) * eps ) then
exit
end if
hw = f1
end do
hf = f0 * c0 + f1 * c1
else if ( m < 0 ) then
m = - m
c0 = gm * gc / ( ga * gb * ( 1.0D+00 - x )**m )
c1 = - ( - 1 )**m * gc / ( gam * gbm * rm )
do k = 1, m - 1
r0 = r0 * ( a - m + k - 1.0D+00 ) * ( b - m + k - 1.0D+00 ) &
/ ( k * ( k - m ) ) * ( 1.0D+00 - x )
f0 = f0 + r0
end do
do k = 1, m
sp0 = sp0 + 1.0D+00 / real ( k, kind = 8 )
end do
f1 = pa + pb - sp0 + 2.0D+00 * el + log ( 1.0D+00 - x )
do k = 1, 250
sp = sp + ( 1.0D+00 - a ) &
/ ( k * ( a + k - 1.0D+00 ) ) &
+ ( 1.0D+00 - b ) / ( k * ( b + k - 1.0D+00 ) )
sm = 0.0D+00
do j = 1, m
sm = sm + 1.0D+00 / real ( j + k, kind = 8 )
end do
rp = pa + pb + 2.0D+00 * el + sp - sm + log ( 1.0D+00 - x )
r1 = r1 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) &
/ ( k * ( m + k ) ) * ( 1.0D+00 - x )
f1 = f1 + r1 * rp
if ( abs ( f1 - hw ) < abs ( f1 ) * eps ) then
exit
end if
hw = f1
end do
hf = f0 * c0 + f1 * c1
end if
else
call gamma_f ( a, ga )
call gamma_f ( b, gb )
call gamma_f ( c, gc )
call gamma_f ( c - a, gca )
call gamma_f ( c - b, gcb )
call gamma_f ( c - a - b, gcab )
call gamma_f ( a + b - c, gabc )
c0 = gc * gcab / ( gca * gcb )
c1 = gc * gabc / ( ga * gb ) * ( 1.0D+00 - x )**( c - a - b )
hf = 0.0D+00
r0 = c0
r1 = c1
do k = 1, 250
r0 = r0 * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) &
/ ( k * ( a + b - c + k ) ) * ( 1.0D+00 - x )
r1 = r1 * ( c - a + k - 1.0D+00 ) * ( c - b + k - 1.0D+00 ) &
/ ( k * ( c - a - b + k ) ) * ( 1.0D+00 - x )
hf = hf + r0 + r1
if ( abs ( hf - hw ) < abs ( hf ) * eps ) then
exit
end if
hw = hf
end do
hf = hf + c0 + c1
end if
else
a0 = 1.0D+00
if ( a < c .and. c < 2.0D+00 * a .and. b < c .and. c < 2.0D+00 * b ) then
a0 = ( 1.0D+00 - x )**( c - a - b )
a = c - a
b = c - b
end if
hf = 1.0D+00
r = 1.0D+00
do k = 1, 250
r = r * ( a + k - 1.0D+00 ) * ( b + k - 1.0D+00 ) &
/ ( k * ( c + k - 1.0D+00 ) ) * x
hf = hf + r
if ( abs ( hf - hw ) <= abs ( hf ) * eps ) then
exit
end if
hw = hf
end do
hf = a0 * hf
end if
if ( x1 < 0.0D+00 ) then
x = x1
c0 = 1.0D+00 / ( 1.0D+00 - x )**aa
hf = c0 * hf
end if
a = aa
b = bb
if ( 120 < k ) then
write ( *, '(a)' ) ' '
write ( *, '(a)' ) 'HYGFX - Warning!'
write ( *, '(a)' ) ' A large number of iterations were needed.'
write ( *, '(a)' ) ' The accuracy of the results should be checked.'
end if
return
end subroutine hygfx
subroutine gamma_f ( x, ga )
!*****************************************************************************80
!
!! GAMMA evaluates the Gamma function.
!
! Licensing:
!
! The original FORTRAN77 version of this routine is copyrighted by
! Shanjie Zhang and Jianming Jin. However, they give permission to
! incorporate this routine into a user program that the copyright
! is acknowledged.
!
! Modified:
!
! 08 September 2007
!
! Author:
!
! Original FORTRAN77 version by Shanjie Zhang, Jianming Jin.
! FORTRAN90 version by John Burkardt.
!
! Reference:
!
! Shanjie Zhang, Jianming Jin,
! Computation of Special Functions,
! Wiley, 1996,
! ISBN: 0-471-11963-6,
! LC: QA351.C45
!
! Parameters:
!
! Input, real ( kind = 8 ) X, the argument.
! X must not be 0, or any negative integer.
!
! Output, real ( kind = 8 ) GA, the value of the Gamma function.
!
implicit none
real ( kind = 8 ), dimension ( 26 ) :: g = (/ &
1.0D+00, &
0.5772156649015329D+00, &
-0.6558780715202538D+00, &
-0.420026350340952D-01, &
0.1665386113822915D+00, &
-0.421977345555443D-01, &
-0.96219715278770D-02, &
0.72189432466630D-02, &
-0.11651675918591D-02, &
-0.2152416741149D-03, &
0.1280502823882D-03, &
-0.201348547807D-04, &
-0.12504934821D-05, &
0.11330272320D-05, &
-0.2056338417D-06, &
0.61160950D-08, &
0.50020075D-08, &
-0.11812746D-08, &
0.1043427D-09, &
0.77823D-11, &
-0.36968D-11, &
0.51D-12, &
-0.206D-13, &
-0.54D-14, &
0.14D-14, &
0.1D-15 /)
real ( kind = 8 ) ga
real ( kind = 8 ) gr
integer ( kind = 4 ) k
integer ( kind = 4 ) m
integer ( kind = 4 ) m1
real ( kind = 8 ), parameter :: pi = 3.141592653589793D+00
real ( kind = 8 ) r
real ( kind = 8 ) x
real ( kind = 8 ) z
if ( x == aint ( x ) ) then
if ( 0.0D+00 < x ) then
ga = 1.0D+00
m1 = int ( x ) - 1
do k = 2, m1
ga = ga * k
end do
else
ga = 1.0D+300
end if
else
if ( 1.0D+00 < abs ( x ) ) then
z = abs ( x )
m = int ( z )
r = 1.0D+00
do k = 1, m
r = r * ( z - real ( k, kind = 8 ) )
end do
z = z - real ( m, kind = 8 )
else
z = x
end if
gr = g(26)
do k = 25, 1, -1
gr = gr * z + g(k)
end do
ga = 1.0D+00 / ( gr * z )
if ( 1.0D+00 < abs ( x ) ) then
ga = ga * r
if ( x < 0.0D+00 ) then
ga = - pi / ( x* ga * sin ( pi * x ) )
end if
end if
end if
return
end subroutine gamma_f
subroutine psi ( x, ps )
!*****************************************************************************80
!
!! PSI computes the PSI function.
!
! Licensing:
!
! The original FORTRAN77 version of this routine is copyrighted by
! Shanjie Zhang and Jianming Jin. However, they give permission to
! incorporate this routine into a user program that the copyright
! is acknowledged.
!
! Modified:
!
! 08 September 2007
!
! Author:
!
! Original FORTRAN77 by Shanjie Zhang, Jianming Jin.
! FORTRAN90 version by John Burkardt.
!
! Reference:
!
! Shanjie Zhang, Jianming Jin,
! Computation of Special Functions,
! Wiley, 1996,
! ISBN: 0-471-11963-6,
! LC: QA351.C45
!
! Parameters:
!
! Input, real ( kind = 8 ) X, the argument.
!
! Output, real ( kind = 8 ) PS, the value of the PSI function.
!
implicit none
real ( kind = 8 ), parameter :: a1 = -0.83333333333333333D-01
real ( kind = 8 ), parameter :: a2 = 0.83333333333333333D-02
real ( kind = 8 ), parameter :: a3 = -0.39682539682539683D-02
real ( kind = 8 ), parameter :: a4 = 0.41666666666666667D-02
real ( kind = 8 ), parameter :: a5 = -0.75757575757575758D-02
real ( kind = 8 ), parameter :: a6 = 0.21092796092796093D-01
real ( kind = 8 ), parameter :: a7 = -0.83333333333333333D-01
real ( kind = 8 ), parameter :: a8 = 0.4432598039215686D+00
real ( kind = 8 ), parameter :: el = 0.5772156649015329D+00
integer ( kind = 4 ) k
integer ( kind = 4 ) n
real ( kind = 8 ), parameter :: pi = 3.141592653589793D+00
real ( kind = 8 ) ps
real ( kind = 8 ) s
real ( kind = 8 ) x
real ( kind = 8 ) x2
real ( kind = 8 ) xa
xa = abs ( x )
s = 0.0D+00
if ( x == aint ( x ) .and. x <= 0.0D+00 ) then
ps = 1.0D+300
return
else if ( xa == aint ( xa ) ) then
n = int ( xa )
do k = 1, n - 1
s = s + 1.0D+00 / real ( k, kind = 8 )
end do
ps = - el + s
else if ( xa + 0.5D+00 == aint ( xa + 0.5D+00 ) ) then
n = int ( xa - 0.5D+00 )
do k = 1, n
s = s + 1.0D+00 / real ( 2 * k - 1, kind = 8 )
end do
ps = - el + 2.0D+00 * s - 1.386294361119891D+00
else
if ( xa < 10.0D+00 ) then
n = 10 - int ( xa )
do k = 0, n - 1
s = s + 1.0D+00 / ( xa + real ( k, kind = 8 ) )
end do
xa = xa + real ( n, kind = 8 )
end if
x2 = 1.0D+00 / ( xa * xa )
ps = log ( xa ) - 0.5D+00 / xa + x2 * ((((((( &
a8 &
* x2 + a7 ) &
* x2 + a6 ) &
* x2 + a5 ) &
* x2 + a4 ) &
* x2 + a3 ) &
* x2 + a2 ) &
* x2 + a1 )
ps = ps - s
end if
if ( x < 0.0D+00 ) then
ps = ps - pi * cos ( pi * x ) / sin ( pi * x ) - 1.0D+00 / x
end if
return
end subroutine psi
!Fits data (vx,vy) with a polynomial of degree d
!Source: http://rosettacode.org/wiki/Polynomial_regression
function polyfit(vx, vy, d)
implicit none
integer, intent(in) :: d
double precision, dimension(d+1) :: polyfit
double precision, dimension(:), intent(in) :: vx, vy
double precision, dimension(:,:), allocatable :: X
double precision, dimension(:,:), allocatable :: XT
double precision, dimension(:,:), allocatable :: XTX
integer :: i, j
integer :: n, lda, lwork
integer :: info
integer, dimension(:), allocatable :: ipiv
double precision, dimension(:), allocatable :: work
n = d+1
lda = n
lwork = n
allocate(ipiv(n))
allocate(work(lwork))
allocate(XT(n, size(vx)))
allocate(X(size(vx), n))
allocate(XTX(n, n))
! prepare the matrix
do i = 0, d
do j = 1, size(vx)
X(j, i+1) = vx(j)**i
end do
end do
XT = transpose(X)
XTX = matmul(XT, X)
! calls to LAPACK subs DGETRF and DGETRI
call DGETRF(n, n, XTX, lda, ipiv, info)
if ( info /= 0 ) then
print *, "problem"
return
end if
call DGETRI(n, XTX, lda, ipiv, work, lwork, info)
if ( info /= 0 ) then
print *, "problem"
return
end if
polyfit = matmul( matmul(XTX, XT), vy)
deallocate(ipiv)
deallocate(work)
deallocate(X)
deallocate(XT)
deallocate(XTX)
end function
end MODULE special