-
Notifications
You must be signed in to change notification settings - Fork 7
/
mpfund.f90
executable file
·2495 lines (1957 loc) · 59.8 KB
/
mpfund.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
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!*****************************************************************************
! MPFUN-Fort: A thread-safe arbitrary precision computation package
! Transcendental function module (module MPFUND)
! Revision date: 6 Feb 2016
! AUTHOR:
! David H. Bailey
! Lawrence Berkeley National Lab (retired) and University of California, Davis
! Email: dhbailey@lbl.gov
! COPYRIGHT AND DISCLAIMER:
! All software in this package (c) 2015 David H. Bailey.
! By downloading or using this software you agree to the copyright, disclaimer
! and license agreement in the accompanying file DISCLAIMER.txt.
! PURPOSE OF PACKAGE:
! This package permits one to perform floating-point computations (real and
! complex) to arbitrarily high numeric precision, by making only relatively
! minor changes to existing Fortran-90 programs. All basic arithmetic
! operations and transcendental functions are supported, together with several
! special functions.
! In addition to fast execution times, one key feature of this package is a
! 100% THREAD-SAFE design, which means that user-level applications can be
! easily converted for parallel execution, say using a threaded parallel
! environment such as OpenMP. There are NO global shared variables (except
! static compile-time data), and NO initialization is necessary unless
! extremely high precision (> 19,500 digits) is required.
! DOCUMENTATION:
! A detailed description of this package, and instructions for compiling
! and testing this program on various specific systems are included in the
! README file accompanying this package, and, in more detail, in the
! following technical paper:
! David H. Bailey, "MPFUN2015: A thread-safe arbitrary precision package,"
! available at http://www.davidhbailey.com/dhbpapers/mpfun2015.pdf.
! DESCRIPTION OF THIS MODULE (MPFUND):
! This module contains subroutines for basic transcendental functions,
! including routines for cos, sin, inverse cos/sin, hyperbolic cos/sin,
! and inverse hyperbolic cos, sin, as well as routines to compute pi and log(2).
module mpfund
use mpfuna
use mpfunb
use mpfunc
contains
subroutine mpagmr (a, b, c, mpnw)
! This performs the arithmetic-geometric mean (AGM) iterations on A and B.
! The AGM algorithm is as follows: Set a_0 = a and b_0 = b, then iterate
! a_{k+1} = (a_k + b_k)/2
! b_{k+1} = sqrt (a_k * b_k)
! until convergence (i.e., until a_k = b_k to available precision).
! The result is returned in C.
implicit none
integer i, itrmx, j, mpnw, mpnw1
parameter (itrmx = 50)
double precision a(0:mpnw+5), b(0:mpnw+5), c(0:mpnw+5), s0(0:mpnw+6), &
s1(0:mpnw+6), s2(0:mpnw+6), s3(0:mpnw+6)
if (mpnw < 4 .or. a(0) < mpnw + 4 .or. b(0) < abs (a(2)) + 4 .or. &
c(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPAGMR: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
s3(0) = mpnw + 7
mpnw1 = mpnw + 1
call mpeq (a, s1, mpnw1)
call mpeq (b, s2, mpnw1)
do j = 1, itrmx
call mpadd (s1, s2, s0, mpnw1)
call mpmuld (s0, 0.5d0, s3, mpnw1)
call mpmul (s1, s2, s0, mpnw1)
call mpsqrt (s0, s2, mpnw1)
call mpeq (s3, s1, mpnw1)
call mpsub (s1, s2, s0, mpnw1)
! Check for convergence.
if (s0(2) .eq. 0.d0 .or. s0(3) < 1.d0 - mpnw1) goto 100
enddo
write (mpldb, 2)
2 format ('*** MPAGMR: Iteration limit exceeded.')
call mpabrt (5)
100 continue
call mproun (s1, mpnw)
call mpeq (s1, c, mpnw)
return
end subroutine mpagmr
subroutine mpang (x, y, a, mpnw)
! This computes the MPR angle A subtended by the MPR pair (X, Y) considered as
! a point in the x-y plane. This is more useful than an arctan or arcsin
! routine, since it places the result correctly in the full circle, i.e.
! -Pi < A <= Pi. Pi must be precomputed to at least MPNW words precision
! and the stored in the array in module MPMODA.
! The Taylor series for Sin converges much more slowly than that of Arcsin.
! Thus this routine does not employ Taylor series, but instead computes
! Arccos or Arcsin by solving Cos (a) = x or Sin (a) = y using one of the
! following Newton iterations, both of which converge to a:
! z_{k+1} = z_k - [x - Cos (z_k)] / Sin (z_k)
! z_{k+1} = z_k + [y - Sin (z_k)] / Cos (z_k)
! The first is selected if Abs (x) <= Abs (y); otherwise the second is used.
! These iterations are performed with a maximum precision level MPNW that
! is dynamically changed, approximately doubling with each iteration.
! If the precision level MPNW exceeds MPNWX words, this subroutine calls
! MPANGX instead. By default, MPNWX = 100 (approx. 1450 digits).
implicit none
integer i, iq, ix, iy, k, kk, mpnw, mpnwx, mpnw1, mq, nit, nx, ny, n1, n2
double precision cl2, cpi, mprxx, t1, t2, t3
parameter (cl2 = 1.4426950408889633d0, cpi = 3.141592653589793d0, &
mprxx = 1d-14, mpnwx = 100, nit = 3)
double precision a(0:mpnw+5), x(0:mpnw+5), y(0:mpnw+5), s0(0:mpnw+6), &
s1(0:mpnw+6), s2(0:mpnw+6), s3(0:mpnw+6), s4(0:mpnw+6), s5(0:mpnw+6)
! End of declaration
if (mpnw < 4 .or. x(0) < mpnw + 4 .or. x(0) < abs (x(2)) + 4 .or. &
y(0) < mpnw + 4 .or. y(0) < abs (y(2)) + 4 .or. a(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPANG: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ix = sign (1.d0, x(2))
nx = min (int (abs (x(2))), mpnw)
iy = sign (1.d0, y(2))
ny = min (int (abs (y(2))), mpnw)
mpnw1 = mpnw + 1
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
s3(0) = mpnw + 7
s4(0) = mpnw + 7
s5(0) = mpnw + 7
! Check if both X and Y are zero.
if (nx .eq. 0 .and. ny .eq. 0) then
write (mpldb, 2)
2 format ('*** MPANG: Both arguments are zero.')
call mpabrt (7)
endif
! Check if Pi has been precomputed.
if (mpnw1 > mppicon(1)) then
write (mpldb, 3) mpnw1
3 format ('*** MPANG: Pi must be precomputed to precision',i9,' words.'/ &
'See documentation for details.')
call mpabrt (8)
endif
! If the precision level mpnw exceeds mpnwx words, call mpangx.
if (mpnw > mpnwx) then
call mpangx (x, y, a, mpnw)
goto 120
endif
! Check if one of X or Y is zero.
if (nx .eq. 0) then
call mpeq (mppicon, s0, mpnw1)
if (iy .gt. 0) then
call mpmuld (s0, 0.5d0, a, mpnw)
else
call mpmuld (s0, -0.5d0, a, mpnw)
endif
goto 120
elseif (ny .eq. 0) then
if (ix .gt. 0) then
a(1) = mpnw
a(2) = 0.d0
a(3) = 0.d0
else
call mpeq (s0, a, mpnw)
endif
goto 120
endif
! Determine the least integer MQ such that 2 ^ MQ .GE. MPNW.
t1 = mpnw1
mq = cl2 * log (t1) + 1.d0 - mprxx
! Normalize x and y so that x^2 + y^2 = 1.
call mpmul (x, x, s0, mpnw1)
call mpmul (y, y, s1, mpnw1)
call mpadd (s0, s1, s2, mpnw1)
call mpsqrt (s2, s3, mpnw1)
call mpdiv (x, s3, s1, mpnw1)
call mpdiv (y, s3, s2, mpnw1)
! Compute initial approximation of the angle.
call mpmdc (s1, t1, n1, mpnw1)
call mpmdc (s2, t2, n2, mpnw1)
n1 = max (n1, -mpnbt)
n2 = max (n2, -mpnbt)
t1 = t1 * 2.d0 ** n1
t2 = t2 * 2.d0 ** n2
t3 = atan2 (t2, t1)
call mpdmc (t3, 0, s5, mpnw1)
! The smaller of x or y will be used from now on to measure convergence.
! This selects the Newton iteration (of the two listed above) that has the
! largest denominator.
if (abs (t1) .le. abs (t2)) then
kk = 1
call mpeq (s1, s0, mpnw1)
else
kk = 2
call mpeq (s2, s0, mpnw1)
endif
mpnw1 = 4
iq = 0
! Perform the Newton-Raphson iteration described above with a dynamically
! changing precision level MPNW (one greater than powers of two).
do k = 1, mq
mpnw1 = min (2 * mpnw1 - 2, mpnw) + 1
100 continue
call mpcssnr (s5, s1, s2, mpnw1)
if (kk .eq. 1) then
call mpsub (s0, s1, s3, mpnw1)
call mpdiv (s3, s2, s4, mpnw1)
call mpsub (s5, s4, s1, mpnw1)
else
call mpsub (s0, s2, s3, mpnw1)
call mpdiv (s3, s1, s4, mpnw1)
call mpadd (s5, s4, s1, mpnw1)
endif
call mpeq (s1, s5, mpnw1)
if (k == mq - nit .and. iq == 0) then
iq = 1
goto 100
endif
enddo
! Restore original precision level.
call mproun (s5, mpnw)
call mpeq (s5, a, mpnw)
120 continue
return
end subroutine mpang
subroutine mpangx (x, y, a, mpnw)
! This computes the MPR angle A subtended by the MPR pair (X, Y) considered as
! a point in the x-y plane. This is more useful than an arctan or arcsin
! routine, since it places the result correctly in the full circle, i.e.
! -Pi < A <= Pi. Pi and Log(2) must be precomputed to at least MPNW words
! precision and the stored in the array in module MPMODA.
! This routine simply calls mpclogx. For modest precision, use mpang.
implicit none
integer i, mpnw, mp7
double precision a(0:mpnw+5), x(0:mpnw+5), y(0:mpnw+5), s0(0:2*mpnw+13), &
s1(0:2*mpnw+13)
! End of declaration
if (mpnw < 4 .or. x(0) < mpnw + 4 .or. x(0) < abs (x(2)) + 4 .or. &
y(0) < mpnw + 4 .or. y(0) < abs (y(2)) + 4 .or. a(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPANGX: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
mp7 = mpnw + 7
s0(0) = mp7
s0(mp7) = mp7
s1(0) = mp7
s1(mp7) = mp7
call mpeq (x, s0, mpnw)
call mpeq (y, s0(mp7), mpnw)
call mpclogx (s0, s1, mpnw)
call mpeq (s1(mp7), a, mpnw)
return
end subroutine mpangx
subroutine mpcagm (a, b, c, mpnw)
! This performs the arithmetic-geometric mean (AGM) iterations on A and B
! for MPC arguments A and B.
! The AGM algorithm is as follows: Set a_0 = a and b_0 = b, then iterate
! a_{k+1} = (a_k + b_k)/2
! b_{k+1} = sqrt (a_k * b_k)
! until convergence (i.e., until a_k = b_k to available precision).
! The result is returned in C.
implicit none
integer i, itrmx, j, la, lb, lc, mp7, mpnw, mpnw1
parameter (itrmx = 50)
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(0:2*mpnw+11), &
s0(0:2*mpnw+13), s1(0:2*mpnw+13), s2(0:2*mpnw+13), s3(0:2*mpnw+13)
la = a(0)
lb = b(0)
lc = c(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < abs (b(2)) + 4 .or. b(lb) < abs (b(lb+2)) + 4 .or. &
c(0) < mpnw + 6 .or. c(lc) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCAGM: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
mp7 = mpnw + 7
s0(0) = mp7
s0(mp7) = mp7
s1(0) = mp7
s1(mp7) = mp7
s2(0) = mp7
s2(mp7) = mp7
s3(0) = mp7
s3(mp7) = mp7
mpnw1 = mpnw + 1
call mpceq (a, s1, mpnw1)
call mpceq (b, s2, mpnw1)
do j = 1, itrmx
call mpcadd (s1, s2, s0, mpnw1)
call mpmuld (s0, 0.5d0, s3, mpnw1)
call mpmuld (s0(mp7), 0.5d0, s3(mp7), mpnw1)
call mpcmul (s1, s2, s0, mpnw1)
call mpcsqrt (s0, s2, mpnw1)
call mpceq (s3, s1, mpnw1)
call mpcsub (s1, s2, s0, mpnw1)
! Check for convergence.
if ((s0(2) .eq. 0.d0 .or. s0(3) < 1.d0 - mpnw1) .and. &
(s0(mp7+2) .eq. 0.d0 .or. s0(mp7+3) < 1.d0 - mpnw1)) goto 100
enddo
write (mpldb, 2)
2 format ('*** MPCAGM: Iteration limit exceeded.')
call mpabrt (5)
100 continue
call mproun (s1, mpnw)
call mproun (s1(mp7), mpnw)
call mpceq (s1, c, mpnw)
return
end subroutine mpcagm
subroutine mpcexp (a, b, mpnw)
! This computes Exp[A], for MPC A.
! The formula is: E^a1 * (Cos[a2] + I * Sin[a2]), where a1 and a2 are
! the real and imaginary parts of A.
! If the precision level MPNW exceeds MPNWX words, this subroutine calls
! MPCEXP instead. By default, MPNWX = 300 (about 4300 digits).
implicit none
integer la, lb, mpnw, mpnwx, mpnw1
parameter (mpnwx = 300)
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), s0(0:mpnw+6), &
s1(0:mpnw+6), s2(0:mpnw+6), s3(0:mpnw+6), s4(0:mpnw+6)
! End of declaration
la = a(0)
lb = b(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < mpnw + 6 .or. b(lb) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCEXP: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
! If the precision level mpnw exceeds mpnwx, call mpcexpx.
if (mpnw > mpnwx) then
call mpcexpx (a, b, mpnw)
goto 100
endif
mpnw1 = mpnw + 1
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
s3(0) = mpnw + 7
s4(0) = mpnw + 7
call mpexp (a, s0, mpnw1)
call mpcssnr (a(la), s1, s2, mpnw1)
call mpmul (s0, s1, s3, mpnw1)
call mpmul (s0, s2, s4, mpnw1)
call mproun (s3, mpnw)
call mproun (s4, mpnw)
call mpeq (s3, b, mpnw)
call mpeq (s4, b(lb), mpnw)
100 continue
return
end subroutine mpcexp
subroutine mpcexpx (a, b, mpnw)
! This computes the exponential of the MPC number A and returns the MPC
! result in B.
! This routine employs the following Newton iteration, which converges to b:
! x_{k+1} = x_k + x_k * [a - Log (x_k)]
! These iterations are performed with a maximum precision level MPNW that
! is dynamically changed, approximately doubling with each iteration.
! For modest levels of precision, use mpcexp.
implicit none
integer i, ia, iq, k, la, lb, mpnw, mpnw1, mp7, mq, na, nb, nit, n0, n1, n2
double precision cl2, t0, t1, t2, mprxx
parameter (cl2 = 1.4426950408889633d0, nit = 3, mprxx = 1d-14)
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), s0(0:2*mpnw+13), &
s1(0:2*mpnw+13), s2(0:2*mpnw+13), s3(0:2*mpnw+13), s4(0:2*mpnw+13), &
r1(0:mpnw+6), r2(0:mpnw+6)
! End of declaration
la = a(0)
lb = b(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < mpnw + 6 .or. b(lb) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCEXPX: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia = sign (1.d0, a(2))
na = min (int (abs (a(2))), mpnw)
call mpmdc (a, t1, n1, mpnw)
! Check for overflows and underflows.
if (n1 > 30) then
if (t1 > 0.d0) then
write (mpldb, 2)
2 format ('*** MPCEXPX: Real part of argument is too large.')
call mpabrt (34)
else
b(1) = mpnw
b(2) = 0.d0
b(3) = 0.d0
b(nb+1) = mpnw
b(nb+2) = 0.d0
b(nb+3) = 0.d0
goto 130
endif
endif
t1 = t1 * 2.d0 ** n1
if (abs (t1) > 1488522236.d0) then
if (t1 > 0) then
write (mpldb, 2)
call mpabrt (34)
else
b(1) = mpnw
b(2) = 0.d0
b(3) = 0.d0
b(nb+1) = mpnw
b(nb+2) = 0.d0
b(nb+3) = 0.d0
goto 130
endif
endif
! Check if imaginary part is too large to compute meaningful cos/sin values.
call mpmdc (a(la), t1, n1, mpnw)
if (n1 >= mpnbt * (mpnw - 1)) then
write (mpldb, 3)
3 format ('*** MPCEXPX: imaginary part is too large to compute cos or sin.')
call mpabrt (28)
endif
mpnw1 = mpnw + 1
mp7 = mpnw + 7
s0(0) = mp7
s0(mp7) = mp7
s1(0) = mp7
s1(mp7) = mp7
s2(0) = mp7
s2(mp7) = mp7
s3(0) = mp7
s3(mp7) = mp7
s4(0) = mp7
s4(mp7) = mp7
r1(0) = mp7
r2(0) = mp7
! Check if Pi and Log(2) have been precomputed.
if (mpnw1 > mplog2con(1)) then
write (mpldb, 4) mpnw1
4 format ('*** MPCEXPX: Pi and Log(2) must be precomputed to precision',i9,' words.'/ &
'See documentation for details.')
call mpabrt (53)
endif
! Reduce imaginary part to within -pi and pi.
call mpeq (a, s4, mpnw1)
call mpmuld (mppicon, 2.d0, s0, mpnw1)
call mpdiv (a(la), s0, s1, mpnw1)
call mpnint (s1, s2, mpnw1)
call mpmul (s2, s0, s1, mpnw1)
call mpsub (a(la), s1, s4(mp7), mpnw1)
! Check if imaginary part is -pi; if so correct to +pi.
call mpadd (s4(mp7), mppicon, s2, mpnw1)
if (s2(2) <= 0.d0 .or. s2(3) < - mpnw) then
call mpeq (mppicon, s4(mp7), mpnw1)
endif
! Determine the least integer MQ such that 2 ^ MQ .GE. MPNW.
t2 = mpnw1
mq = cl2 * log (t2) + 1.d0 - mprxx
! Compute initial approximation of Exp (A) (DP accuracy is OK).
mpnw1 = 4
! The following code (between here and iq = 0) is the equivalent of:
! call mpcexp (s4, s3, mpnw1)
call mpdiv (s4, mplog2con, s0, mpnw1)
call mpinfr (s0, s1, s2, mpnw1)
call mpmdc (s1, t1, n1, mpnw1)
n1 = int (t1 * 2.d0**n1)
call mpmdc (s2, t0, n0, mpnw1)
n0 = min (max (n0, -100), 0)
t0 = 2.d0 ** (t0 * 2.d0**n0)
call mpdmc (t0, n1, s0, mpnw1)
call mpmdc (s4(mp7), t0, n0, mpnw1)
t0 = t0 * 2.d0 ** n0
t1 = cos (t0)
t2 = sin (t0)
if (abs (t1) < 1d-14) t1 = 0.d0
if (abs (t2) < 1d-14) t2 = 0.d0
call mpdmc (t1, 0, s1, mpnw1)
call mpdmc (t2, 0, s1(mp7), mpnw1)
call mpmul (s0, s1, s3, mpnw1)
call mpmul (s0, s1(mp7), s3(mp7), mpnw1)
iq = 0
! Perform the Newton-Raphson iteration described above with a dynamically
! changing precision level MPNW (one greater than powers of two).
do k = 0, mq
if (k > 1) mpnw1 = min (2 * mpnw1 - 2, mpnw) + 1
100 continue
call mpclogx (s3, s0, mpnw1)
! Check if we need to add or subtract 2*pi to the output of imaginary part,
! in order to remain consistent with previous iterations.
call mpsub (s4(mp7), s0(mp7), r1, 4)
if (r1(2) /= 0.d0 .and. r1(3) >= -1.d0) then
if (r1(2) > 0.d0) then
call mpadd (s0(mp7), mppicon, r1, mpnw1)
call mpadd (r1, mppicon, s0(mp7), mpnw1)
else
call mpsub (s0(mp7), mppicon, r1, mpnw1)
call mpsub (r1, mppicon, s0(mp7), mpnw1)
endif
endif
call mpcsub (s4, s0, s1, mpnw1)
call mpcmul (s3, s1, s2, mpnw1)
call mpcadd (s3, s2, s1, mpnw1)
call mpceq (s1, s3, mpnw1)
if (k .eq. mq - nit .and. iq .eq. 0) then
iq = 1
goto 100
endif
enddo
! Restore original precision level.
call mproun (s1, mpnw)
call mproun (s1(mp7), mpnw)
call mpceq (s1, b, mpnw)
130 continue
return
end subroutine mpcexpx
subroutine mpclog (a, b, mpnw)
! This computes Log[A], for MPC A.
! The formula is: 1/2 * Log[r] + I * Theta, where r = a1^2 + a2^2,
! Theta is the angle corresponding to (a1, a2), and a1 and a2 are the
! real and imaginary parts of A.
! If the precision level MPNW exceeds MPNWX words, this subroutine calls
! MPCLOGX instead. By default, MPNWX = 30 (about 400 digits).
implicit none
integer la, lb, mpnw, mpnwx, mpnw1
parameter (mpnwx = 30)
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), s0(0:mpnw+6), &
s1(0:mpnw+6), s2(0:mpnw+6), s3(0:mpnw+6), s4(0:mpnw+6)
! End of declaration
la = a(0)
lb = b(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < mpnw + 6 .or. b(lb) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCLOG: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
! If precision level mpnw exceeds mpnwx words, call mpclogx.
if (mpnw > mpnwx) then
call mpclogx (a, b, mpnw)
goto 100
endif
mpnw1 = mpnw + 1
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
s3(0) = mpnw + 7
s4(0) = mpnw + 7
call mpmul (a, a, s0, mpnw1)
call mpmul (a(la), a(la), s1, mpnw1)
call mpadd (s0, s1, s2, mpnw1)
call mplog (s2, s3, mpnw1)
call mpmuld (s3, 0.5d0, s0, mpnw1)
call mpang (a, a(la), s1, mpnw1)
call mproun (s0, mpnw)
call mproun (s1, mpnw)
call mpeq (s0, b, mpnw)
call mpeq (s1, b(lb), mpnw)
100 continue
return
end subroutine mpclog
subroutine mpclogx (a, b, mpnw)
! This computes the natural logarithm of the MP number A and returns the MP
! result in B. Pi and Log(2) must be precomputed to at least MPNW words
! precision and the stored in the arrays in module MPMODA.
! This uses the following algorithm, which is due to Salamin and Brent. If
! A is extremely close to 1, use a Taylor series. Otherwise select n such
! that z = a * 2^n is at least 2^m, where m is the number of bits of desired
! precision in the result. Then
! Log(x) = Pi / [2 AGM (1, 4/x)]
! For modest precision, or if A is close to 2, use mplog.
implicit none
integer i, ia1, ia2, is, itrmax, i1, i2, la, lb, k, mpnw, mpnw1, mp7, &
na1, na2, n1, n2
double precision st, rtol, tol, t1, t2, tn
parameter (itrmax = 1000000, rtol = 0.5d0**7)
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), f1(0:18), f4(0:18), &
s0(0:2*mpnw+13), s1(0:2*mpnw+13), s2(0:2*mpnw+13), s3(0:2*mpnw+13), &
s4(0:2*mpnw+13)
! End of declaration
la = a(0)
lb = b(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < mpnw + 6 .or. b(lb) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCLOGX: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia1 = sign (1.d0, a(2))
na1 = min (int (abs (a(2))), mpnw)
ia2 = sign (1.d0, a(la+2))
na2 = min (int (abs (a(la+2))), mpnw)
if (na1 .eq. 0 .and. na2 .eq. 0) then
write (mpldb, 2)
2 format ('*** MPCLOGX: Argument is zero.')
call mpabrt (50)
endif
! Check if input is exactly one.
if (a(2) == 1.d0 .and. a(3) == 0.d0 .and. a(4) == 1.d0 .and. &
a(la+2) == 0.d0) then
b(1) = mpnw
b(2) = 0.d0
b(3) = 0.d0
b(4) = 0.d0
b(lb+1) = mpnw
b(lb+2) = 0.d0
b(lb+3) = 0.d0
b(lb+4) = 0.d0
goto 120
endif
mpnw1 = mpnw + 1
mp7 = mpnw + 7
s0(0) = mp7
s0(mp7) = mp7
s1(0) = mp7
s1(mp7) = mp7
s2(0) = mp7
s2(mp7) = mp7
s3(0) = mp7
s3(mp7) = mp7
s4(0) = mp7
s4(mp7) = mp7
! Check if Pi and Log(2) have been precomputed.
if (mpnw1 > mplog2con(1)) then
write (mpldb, 3) mpnw1
3 format ('*** MPCLOGX: Pi and Log(2) must be precomputed to precision',i9,' words.'/ &
'See documentation for details.')
call mpabrt (53)
endif
f1(0) = 9.d0
f1(1) = mpnw1
f1(2) = 1.d0
f1(3) = 0.d0
f1(4) = 1.d0
f1(5) = 0.d0
f1(6) = 0.d0
f1(9) = 9.d0
f1(10) = mpnw1
f1(11) = 0.d0
f1(12) = 0.d0
f1(13) = 0.d0
f4(0) = 9.d0
f4(1) = mpnw1
f4(2) = 1.d0
f4(3) = 0.d0
f4(4) = 4.d0
f4(5) = 0.d0
f4(6) = 0.d0
f4(9) = 9.d0
f4(10) = mpnw1
f4(11) = 0.d0
f4(12) = 0.d0
f4(13) = 0.d0
! If the argument is sufficiently close to 1, employ a Taylor series.
call mpcsub (a, f1, s0, mpnw1)
if ((s0(2) == 0.d0 .or. s0(3) <= min (-2.d0, - rtol * mpnw1)) .and. &
(s0(mp7+2) == 0.d0 .or. s0(mp7+3) <= min (-2.d0, - rtol * mpnw1))) then
call mpceq (s0, s1, mpnw1)
call mpceq (s1, s2, mpnw1)
i1 = 1
is = 1
if (s0(2) == 0.d0) then
tol = s0(mp7+3) - mpnw1
elseif (s0(mp7+2) == 0.d0) then
tol = s0(3) - mpnw1
else
tol = max (s0(3), s0(mp7+3)) - mpnw1
endif
do i1 = 2, itrmax
is = - is
st = is * i1
call mpcmul (s1, s2, s3, mpnw1)
call mpceq (s3, s2, mpnw1)
call mpdivd (s3, st, s4, mpnw1)
call mpdivd (s3(mp7), st, s4(mp7), mpnw1)
call mpcadd (s0, s4, s3, mpnw1)
call mpceq (s3, s0, mpnw1)
if ((s4(2) == 0.d0 .or. s4(3) < tol) .and. &
(s4(mp7+2) == 0.d0 .or. s4(mp7+3) < tol)) goto 110
enddo
write (mpldb, 4) itrmax
4 format ('*** MPCLOGX: Iteration limit exceeded:',i10)
call mpabrt (54)
endif
! Multiply the input by a large power of two.
call mpmdc (a, t1, n1, mpnw1)
n2 = mpnbt * (mpnw1 / 2 + 2) - n1
tn = n2
call mpdmc (1.d0, n2, s1, mpnw1)
call mpmul (a, s1, s0, mpnw1)
call mpmul (a(la), s1, s0(mp7), mpnw1)
! Perform AGM iterations.
call mpceq (f1, s1, mpnw1)
call mpcdiv (f4, s0, s2, mpnw1)
call mpcagm (s1, s2, s3, mpnw1)
! Compute Pi / (2 * A), where A is the limit of the AGM iterations.
call mpmuld (s3, 2.d0, s0, mpnw1)
call mpmuld (s3(mp7), 2.d0, s0(mp7), mpnw1)
call mpeq (mppicon, s3, mpnw1)
call mpdmc (0.d0, 0, s3(mp7), mpnw1)
call mpcdiv (s3, s0, s1, mpnw1)
! Subtract TN * Log(2).
call mpeq (mplog2con, s3, mpnw1)
call mpmuld (s3, tn, s2, mpnw1)
call mpsub (s1, s2, s0, mpnw1)
! Check if imaginary part is -pi; if so correct to +pi.
call mpadd (s1(mp7), mppicon, s2, mpnw1)
if (s2(2) <= 0.d0 .or. s2(3) < - mpnw) then
call mpeq (mppicon, s0(mp7), mpnw1)
else
call mpeq (s1(mp7), s0(mp7), mpnw1)
endif
110 continue
! Restore original precision level.
call mproun (s0, mpnw)
call mproun (s0(mp7), mpnw)
call mpceq (s0, b, mpnw)
120 continue
return
end subroutine mpclogx
subroutine mpcpowcc (a, b, c, mpnw)
! This computes A^B, where A and B are MPC.
implicit none
integer la, lb, lc, l3, mpnw
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(0:2*mpnw+11), &
s1(0:2*mpnw+11), s2(0:2*mpnw+11)
! End of declaration
la = a(0)
lb = b(0)
lc = c(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < abs (b(2)) + 4 .or. b(lb) < abs (b(lb+2)) + 4 &
.or. c(0) < mpnw + 6 .or. c(lc) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCPOWCC: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
l3 = mpnw + 6
s1(0) = l3
s1(l3) = l3
s2(0) = l3
s2(l3) = l3
call mpclog (a, s1, mpnw)
call mpcmul (s1, b, s2, mpnw)
call mpcexp (s2, c, mpnw)
return
end subroutine mpcpowcc
subroutine mpcpowcr (a, b, c, mpnw)
! This computes A^B, where A is MPC and B is MPR.
implicit none
integer la, lb, lc, l3, mpnw
double precision a(0:2*mpnw+11), b(0:mpnw+5), c(0:2*mpnw+11), &
s1(0:2*mpnw+11), s2(0:2*mpnw+11)
! End of declaration
la = a(0)
lb = b(0)
lc = c(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. a(la) < abs (a(la+2)) + 4 &
.or. b(0) < abs (b(2)) + 4 &
.or. c(0) < mpnw + 6 .or. c(lc) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCPOWCR: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
l3 = mpnw + 6
s1(0) = l3
s1(l3) = l3
s2(0) = l3
s2(l3) = l3
call mpclog (a, s1, mpnw)
call mpmul (b, s1, s2, mpnw)
call mpmul (b, s1(l3), s2(l3), mpnw)
call mpcexp (s2, c, mpnw)
return
end subroutine mpcpowcr
subroutine mpcpowrc (a, b, c, mpnw)
! This computes A^B, where A is MPR and and B is MPC.
implicit none
integer la, lb, lc, l3, mpnw
double precision a(0:mpnw+5), b(0:2*mpnw+11), c(0:2*mpnw+11), &
s1(0:2*mpnw+11), s2(0:2*mpnw+11)
! End of declaration
la = a(0)
lb = b(0)
lc = c(0)
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 &
.or. b(0) < abs (b(2)) + 4 .or. b(lb) < abs (b(lb+2)) + 4 &
.or. c(0) < mpnw + 6 .or. c(lc) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPCPOWRC: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
l3 = mpnw + 6
s1(0) = l3
s2(0) = l3
s2(l3) = l3
call mplog (a, s1, mpnw)
call mpmul (s1, b, s2, mpnw)
call mpmul (s1, b(lb), s2(l3), mpnw)
call mpcexp (s2, c, mpnw)
return
end subroutine mpcpowrc