forked from pastfalk/LEOPARD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmpfunbq.f90
executable file
·3295 lines (2526 loc) · 71 KB
/
mpfunbq.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
! Basic function module (module MPFUNB); includes real*16 support.
! Search for !> for version differences.
! Revision date: 5 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 (MPFUNB):
! This module contains routines for: add, subtract, multiply, divide;
! comparison; double/multi conversion; double/multi multiplication/division;
! integer and fractional parts; nearest integer; nth power; nth root;
! square root; and rounding and normalization. Routines in this package
! are not intended to be directly called by the user; instead, use the
! high-level language interface modules.
module mpfunb
use mpfuna
contains
subroutine mpabrt (ier)
! This routine terminates execution. Users may wish to replace the
! default STOP with a call to a system routine that provides a traceback.
implicit none
integer ier
! End of declaration
write (mpldb, 1) ier
1 format ('*** MPABRT: Execution terminated, error code =',i4)
stop
end subroutine mpabrt
subroutine mpadd (a, b, c, mpnw)
! This routine adds MPR numbers A and B to yield C.
implicit none
double precision a(0:mpnw+5), b(0:mpnw+5), c(0:mpnw+5), d(0:mpnw+6), db
integer i, ia, ib, ic, ish, ixa, ixb, ixd, mpnw, m1, m2, m3, m4, m5, na, nb, &
nd, nsh
! End of declaration
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. b(0) < abs (b(2)) + 4 .or. &
c(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPADD: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia = sign (1.d0, a(2))
ib = sign (1.d0, b(2))
na = min (int (abs (a(2))), mpnw)
nb = min (int (abs (b(2))), mpnw)
! Check for zero inputs.
if (na == 0) then
! A is zero -- the result is B.
c(1) = mpnw
c(2) = sign (nb, ib)
do i = 2, nb + 2
c(i+1) = b(i+1)
enddo
goto 100
elseif (nb == 0) then
! B is zero -- the result is A.
c(1) = mpnw
c(2) = sign (na, ia)
do i = 2, na + 2
c(i+1) = a(i+1)
enddo
goto 100
endif
if (ia == ib) then
db = 1.d0
else
db = -1.d0
endif
ixa = a(3)
ixb = b(3)
ish = ixa - ixb
if (ish >= 0) then
! A has greater exponent than B, so B must be shifted to the right.
m1 = min (na, ish)
m2 = min (na, nb + ish)
m3 = na
m4 = min (max (na, ish), mpnw + 1)
m5 = min (max (na, nb + ish), mpnw + 1)
do i = 1, m1
d(i+3) = a(i+3)
enddo
do i = m1 + 1, m2
d(i+3) = a(i+3) + db * b(i+2-ish+1)
enddo
do i = m2 + 1, m3
d(i+3) = a(i+3)
enddo
do i = m3 + 1, m4
d(i+3) = 0.d0
enddo
do i = m4 + 1, m5
d(i+3) = db * b(i+2-ish+1)
enddo
nd = m5
ixd = ixa
d(nd+4) = 0.d0
d(nd+5) = 0.d0
else
! B has greater exponent than A, so A must be shifted to the right.
nsh = - ish
m1 = min (nb, nsh)
m2 = min (nb, na + nsh)
m3 = nb
m4 = min (max (nb, nsh), mpnw + 1)
m5 = min (max (nb, na + nsh), mpnw + 1)
do i = 1, m1
d(i+3) = db * b(i+3)
enddo
do i = m1 + 1, m2
d(i+3) = a(i+2-nsh+1) + db * b(i+3)
enddo
do i = m2 + 1, m3
d(i+3) = db * b(i+3)
enddo
do i = m3 + 1, m4
d(i+3) = 0.d0
enddo
do i = m4 + 1, m5
d(i+3) = a(i+2-nsh+1)
enddo
nd = m5
ixd = ixb
d(nd+4) = 0.d0
d(nd+5) = 0.d0
endif
! Call mpnorm to fix up result and store in c.
d(0) = mpnw + 7
d(1) = mpnw
d(2) = sign (nd, ia)
d(3) = ixd
call mpnorm (d, c, mpnw)
100 continue
return
end subroutine mpadd
subroutine mpcabs (a, b, mpnw)
! This routine returns the absolute value of the MPC argument A (the
! result is of type MPR).
implicit none
integer la, lb, mpnw, mpnw1
double precision a(0:2*mpnw+11), b(0:2*mpnw+5), s0(0:mpnw+6), s1(0:mpnw+6), &
s2(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) then
write (mpldb, 1)
1 format ('*** MPCABS: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
mpnw1 = mpnw + 1
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
call mpmul (a, a, s0, mpnw1)
call mpmul (a(la), a(la), s1, mpnw1)
call mpadd (s0, s1, s2, mpnw1)
call mpsqrt (s2, s0, mpnw1)
call mproun (s0, mpnw)
call mpeq (s0, b, mpnw)
return
end subroutine mpcabs
subroutine mpcadd (a, b, c, mpnw)
! This routine adds the MPC numbers A and B.
implicit none
integer la, lb, lc, mpnw
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(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 ('*** MPCADD: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
call mpadd (a, b, c, mpnw)
call mpadd (a(la), b(lb), c(lc), mpnw)
return
end subroutine mpcadd
subroutine mpcdiv (a, b, c, mpnw)
! This routine divides the MPC numbers A and B.
implicit none
integer la, lb, lc, mpnw, mpnw1
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(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)
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 ('*** MPCDIV: uninitialized or inadequately sized arrays')
call mpabrt (99)
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, b, s0, mpnw1)
call mpmul (a(la), b(lb), s1, mpnw1)
call mpadd (s0, s1, s2, mpnw1)
call mpmul (a, b(lb), s0, mpnw1)
s0(2) = - s0(2)
call mpmul (a(la), b, s1, mpnw1)
call mpadd (s0, s1, s3, mpnw1)
call mpmul (b, b, s0, mpnw1)
call mpmul (b(lb), b(lb), s1, mpnw1)
call mpadd (s0, s1, s4, mpnw1)
call mpdiv (s2, s4, s0, mpnw1)
call mpdiv (s3, s4, s1, mpnw1)
call mproun (s0, mpnw)
call mproun (s1, mpnw)
call mpeq (s0, c, mpnw)
call mpeq (s1, c(lc), mpnw)
return
end subroutine mpcdiv
subroutine mpceq (a, b, mpnw)
! Sets the MPC number B equal to A.
implicit none
integer i, ia, la, lb, mpnw, na
double precision a(0:2*mpnw+11), b(0:2*mpnw+11)
! 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 ('*** MPCEQ: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia = sign (1.d0, a(2))
na = min (int (abs (a(2))), mpnw)
if (na == 0) then
b(1) = mpnw
b(2) = 0.d0
b(3) = 0.d0
goto 110
endif
b(1) = mpnw
b(2) = sign (na, ia)
do i = 2, na + 2
b(i+1) = a(i+1)
enddo
b(na+4) = 0.d0
b(na+5) = 0.d0
110 continue
ia = sign (1.d0, a(la+2))
na = min (int (abs (a(la+2))), mpnw)
if (na == 0) then
b(lb+1) = mpnw
b(lb+2) = 0.d0
b(lb+3) = 0.d0
goto 120
endif
b(lb+1) = mpnw
b(lb+2) = sign (na, ia)
do i = 2, na + 2
b(i+lb+1) = a(i+la+1)
enddo
b(na+lb+4) = 0.d0
b(na+lb+5) = 0.d0
120 continue
return
end subroutine mpceq
subroutine mpcmul (a, b, c, mpnw)
! This routine multiplies the MPC numbers A and B.
implicit none
integer la, lb, lc, mpnw, mpnw1
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(0:2*mpnw+11), &
s0(0:mpnw+6), s1(0:mpnw+6), s2(0:mpnw+6), s3(0:mpnw+6)
! 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 ('*** MPCMUL: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
mpnw1 = mpnw + 1
s0(0) = mpnw + 7
s1(0) = mpnw + 7
s2(0) = mpnw + 7
s3(0) = mpnw + 7
call mpmul (a, b, s0, mpnw1)
call mpmul (a(la), b(lb), s1, mpnw1)
call mpsub (s0, s1, s2, mpnw1)
call mpmul (a, b(lb), s0, mpnw1)
call mpmul (a(la), b, s1, mpnw1)
call mpadd (s0, s1, s3, mpnw1)
call mproun (s2, mpnw)
call mproun (s3, mpnw)
call mpeq (s2, c, mpnw)
call mpeq (s3, c(lc), mpnw)
return
end subroutine mpcmul
subroutine mpcnpwr (a, n, b, mpnw)
! This computes the N-th power of the MPC number A and returns the MPC result
! in B. When N is zero, 1 is returned. When N is negative, the reciprocal
! of A ^ |N| is returned.
implicit none
integer i, j, k, kk, kn, k0, k1, k2, la, lb, lc, mn, mpnw, mpnw1, n, na, nn, nws
double precision cl2, t1, mprxx
parameter (cl2 = 1.4426950408889633d0, mprxx = 1d-14)
double precision a(0:mpnw+11), b(0:mpnw+11), s0(0:2*mpnw+13), s1(0:2*mpnw+13), &
s2(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 ('*** MPCNPWR: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
na = min (int (abs (a(2))), mpnw)
if (na == 0) then
if (n >= 0) then
b(1) = mpnw
b(2) = 0.d0
b(3) = 0.d0
goto 120
else
write (mpldb, 2)
2 format ('*** MPCNPWR: Argument is zero and N is negative or zero.')
call mpabrt (57)
endif
endif
mpnw1 = mpnw + 1
lc = mpnw + 7
s0(0) = mpnw + 7
s0(lc) = mpnw + 7
s1(0) = mpnw + 7
s1(lc) = mpnw + 7
s2(0) = mpnw + 7
s2(lc) = mpnw + 7
nn = abs (n)
if (nn == 0) then
call mpdmc (1.d0, 0, b, mpnw)
call mpdmc (0.d0, 0, b(lb), mpnw)
goto 120
elseif (nn == 1) then
call mpceq (a, s2, mpnw1)
goto 110
elseif (nn == 2) then
call mpcmul (a, a, s2, mpnw1)
goto 110
endif
! Determine the least integer MN such that 2 ^ MN .GT. NN.
t1 = nn
mn = cl2 * log (t1) + 1.d0 + mprxx
call mpdmc (1.d0, 0, s2, mpnw1)
call mpceq (a, s0, mpnw1)
kn = nn
! Compute B ^ N using the binary rule for exponentiation.
do j = 1, mn
kk = kn / 2
if (kn /= 2 * kk) then
call mpcmul (s2, s0, s1, mpnw1)
call mpceq (s1, s2, mpnw1)
endif
kn = kk
if (j < mn) then
call mpcmul (s0, s0, s1, mpnw1)
call mpceq (s1, s0, mpnw1)
endif
enddo
! Compute reciprocal if N is negative.
110 continue
if (n < 0) then
call mpdmc (1.d0, 0, s1, mpnw1)
call mpdmc (0.d0, 0, s1(lc), mpnw1)
call mpcdiv (s1, s2, s0, mpnw1)
call mpceq (s0, s2, mpnw1)
endif
! Restore original precision level.
call mproun (s2, mpnw)
call mproun (s2(lc), mpnw)
call mpceq (s2, b, mpnw)
120 continue
return
end subroutine mpcnpwr
subroutine mpconjg (a, b, mpnw)
! This routine returns the conjugate of the MPC argument A.
implicit none
integer la, lb, mpnw
double precision a(0:2*mpnw+11), b(0:2*mpnw+11)
! 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 ('*** MPCONJ: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
call mpceq (a, b, mpnw)
b(lb+2) = - b(lb+2)
return
end subroutine mpconjg
subroutine mpcsqrt (a, b, mpnw)
! This routine returns the square root of the MPC argument A.
! The formula is:
! 1/Sqrt[2] * (Sqrt[r + a1] + I * a2 / Sqrt[r + a1]) if a1 >= 0, or
! 1/Sqrt[2] * (|a2| / Sqrt[r - a1] + I * Sgn[a2] * Sqrt[r - a1]) if a1 < 0,
! where r = Sqrt[a1^2 + a2^2], and a1 and a2 are the real and imaginary
! parts of A.
implicit none
integer la, lb, mpnw, mpnw1
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 ('*** MPCSQRT: uninitialized or inadequately sized arrays')
call mpabrt (99)
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 mpsqrt (s2, s0, mpnw1)
if (a(2) >= 0.d0) then
call mpadd (s0, a, s1, mpnw1)
call mpsqrt (s1, s0, mpnw1)
call mpdiv (a(la), s0, s1, mpnw1)
else
call mpsub (s0, a, s2, mpnw1)
call mpsqrt (s2, s1, mpnw1)
call mpdiv (a(la), s1, s0, mpnw1)
s0(2) = abs (s0(2))
s1(2) = sign (s1(2), a(la+2))
endif
call mpeq (mpsqrt22con, s2, mpnw1)
call mpmul (s0, s2, s3, mpnw1)
call mpmul (s1, s2, s4, mpnw1)
call mproun (s3, mpnw)
call mproun (s4, mpnw)
call mpeq (s3, b, mpnw)
call mpeq (s4, b(lb), mpnw)
return
end subroutine mpcsqrt
subroutine mpcsub (a, b, c, mpnw)
! This routine subtracts the MPC numbers A and B.
implicit none
integer la, lb, lc, mpnw
double precision a(0:2*mpnw+11), b(0:2*mpnw+11), c(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 ('*** MPCSUB: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
call mpsub (a, b, c, mpnw)
call mpsub (a(la), b(lb), c(lc), mpnw)
return
end subroutine mpcsub
subroutine mpcpr (a, b, ic, mpnw)
! This routine compares the MPR numbers A and B and returns in IC the value
! -1, 0, or 1 depending on whether A < B, A = B, or A > B.
! Note that the first and second words do NOT need to be the same for the
! result to be "equal".
implicit none
double precision a(0:mpnw+5), b(0:mpnw+5), s0(0:mpnw+5)
integer i, ia, ib, ic, ma, mb, mpnw, na, nb
! End of declaration
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. b(0) < abs (b(2)) + 4) then
write (mpldb, 1)
1 format ('*** MPCPR: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
s0(0) = mpnw + 6
call mpsub (a, b, s0, mpnw)
if (s0(2) < 0.d0) then
ic = -1
elseif (s0(2) == 0.d0) then
ic = 0
else
ic = 1
endif
return
end subroutine mpcpr
subroutine mpdiv (a, b, c, mpnw)
! This divides the MPR number A by the MP number B to yield C.
implicit none
integer i, i1, i2, i3, ia, ib, ij, is, j, j3, md, mpnw, mpnwx, na, nb, nc
parameter (mpnwx = 200)
double precision a1, a2, b1, b2, c1, c2, dc, rb, ss, t0, t1, t2
double precision a(0:mpnw+5), b(0:mpnw+5), c(0:mpnw+5), d(0:mpnw+5)
! End of declaration
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. b(0) < abs (b(2)) + 4 .or. &
c(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPDIV: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia = sign (1.d0, a(2))
ib = sign (1.d0, b(2))
na = min (int (abs (a(2))), mpnw)
nb = min (int (abs (b(2))), mpnw)
! Check if dividend is zero.
if (na .eq. 0) then
c(1) = mpnw
c(2) = 0.d0
c(3) = 0.d0
c(4) = 0.d0
goto 190
endif
! Check if divisor is zero.
if (nb == 0) then
write (mpldb, 2)
2 format ('*** MPDIV: Divisor is zero.')
call mpabrt (31)
endif
if (na > mpnwx .and. nb > mpnwx) then
! Precision is very high, so call mpdivx.
call mpdivx (a, b, c, mpnw)
goto 200
endif
! Compute double precision approximation to trial divisor and its reciprocal.
t0 = mpbdx * b(4)
if (nb >= 2) t0 = t0 + b(5)
if (nb >= 3) t0 = t0 + mprdx * b(6)
rb = 1.d0 / t0
md = min (na + nb, mpnw)
d(0) = mpnw + 6
d(1) = mpnw
d(2) = 0.d0
do i = 2, na + 1
d(i+1) = a(i+2)
enddo
do i = na + 2, md + 4
d(i+1) = 0.d0
enddo
! Perform ordinary long division algorithm.
do j = 2, mpnw + 3
! Compute trial dividend and trial quotient term.
t1 = mpbdx**2 * d(j) + mpbdx * d(j+1) + d(j+2)
if (j <= mpnw + 2) t1 = t1 + mprdx * d(j+3)
t0 = anint (rb * t1)
! Split trial quotient term into high and low order halves, 24-bits each.
a1 = mpb24x * aint (mpr24x * t0)
a2 = t0 - a1
j3 = j - 3
i2 = min (nb, mpnw + 2 - j3) + 2
ij = i2 + j3
! Multiply trial quotient term by each term of divisor, saving high and low-
! order parts into appropriate terms of the running dividend.
do i = 3, i2
i3 = i + j3
b1 = mpb24x * aint (mpr24x * b(i+1))
b2 = b(i+1) - b1
dc = a1 * b2 + a2 * b1
c1 = mpbdx * aint (mprdx * dc)
c2 = dc - c1
d(i3) = d(i3) - mprdx * (a1 * b1 + c1)
d(i3+1) = d(i3+1) - a2 * b2 - c2
enddo
! Release carries periodically to avoid overflowing the exact integer
! capacity of double precision floating point words in D.
if (mod (j - 1, mpnpr) .eq. 0) then
do i = j + 1, ij
t1 = d(i)
t2 = int (mprdx * t1)
d(i) = t1 - mpbdx * t2
d(i-1) = d(i-1) + t2
enddo
endif
d(j+1) = d(j+1) + mpbdx * d(j)
d(j) = t0
! Continue computing quotient terms past the end of the input dividend, until
! trial running dividend is zero.
if (j >= na + 2) then
if (ij <= mpnw + 1) d(ij+4) = 0.d0
endif
enddo
! Final bookkeeping.
j = mpnw + 3
170 continue
d(j+1) = 0.d0
if (d(2) == 0.d0) then
is = 1
else
is = 2
endif
nc = min (j - 1, mpnw)
d(nc+4) = 0.d0
d(nc+5) = 0.d0
do i = j + 1, 3, -1
d(i+1) = d(i-is+1)
enddo
d(0) = mpnw + 6
d(1) = mpnw
d(2) = sign (nc, ia * ib)
d(3) = a(3) - b(3) + is - 2
c(1) = mpnw
! Call mpnorm to fix up any remaining bugs and perform rounding.
call mpnorm (d, c, mpnw)
190 continue
200 continue
return
end subroutine mpdiv
subroutine mpdivd (a, b, c, mpnw)
! This routine divides the MPR number A by the DP number B to yield C.
! NOTE however that if A = 0.1D0, for example, then C will NOT be the true
! multiprecision equivalent of the quotient, since 0.1d0 is not an exact
! binary value.
! Examples of exact binary values (good): 123456789.d0, 0.25d0, -5.3125d0.
! Examples of inexact binary values (bad): 0.1d0, 1234567.8d0, -3333.3d0.
implicit none
integer i, ij, is, i1, i2, i3, ia, ib, j, j3, k, md, mpnw, n, na, nb, nc, n1
double precision a1, a2, b, bb, b1, b2, br, c1, c2, dc, dd, d1, d2, rb, &
t0, t1
double precision a(0:mpnw+5), c(0:mpnw+5), d(0:mpnw+5)
! End of declaration
if (mpnw < 4 .or. a(0) < abs (a(2)) + 4 .or. c(0) < mpnw + 6) then
write (mpldb, 1)
1 format ('*** MPDIVD: uninitialized or inadequately sized arrays')
call mpabrt (99)
endif
ia = sign (1.d0, a(2))
na = min (int (abs (a(2))), mpnw)
ib = sign (1.d0, b)
! Check if dividend is zero.
if (na == 0) then
c(1) = mpnw
c(2) = 0.d0
c(3) = 0.d0
goto 190
endif
! Check if divisor is zero.
if (b == 0.d0) then
write (mpldb, 3)
3 format ('*** MPDIVD: Divisor is zero.')
call mpabrt (32)
endif
n1 = 0
bb = abs (b)
! Reduce BB to within 1 and MPBDX.
if (bb >= mpbdx) then
do k = 1, 100
bb = mprdx * bb
if (bb < mpbdx) then
n1 = n1 + k
goto 120
endif
enddo
elseif (bb < 1.d0) then
do k = 1, 100
bb = mpbdx * bb
if (bb >= 1.d0) then
n1 = n1 - k
goto 120
endif
enddo
endif
120 continue
! If B cannot be represented exactly in a single mantissa word, use MPDIV.
if (bb /= aint (bb)) then
bb = sign (bb, b)
d(0) = mpnw + 6
d(1) = mpnw
call mpdmc (bb, n1 * mpnbt, d, mpnw)
call mpdiv (a, d, c, mpnw)
goto 190
endif
! Compute double precision approximation to trial divisor and its reciprocal.
t0 = mpbdx * bb
rb = 1.d0 / t0
b1 = mpb24x * aint (mpr24x * bb)
b2 = bb - b1
md = min (na + 1, mpnw)
d(0) = mpnw + 6
d(1) = mpnw
d(2) = 0.d0
do i = 2, na + 1
d(i+1) = a(i+2)
enddo
do i = na + 2, md + 4
d(i+1) = 0.d0
enddo
! Perform ordinary short division algorithm.
do j = 2, mpnw + 3
! Compute trial dividend and trial quotient term.
t1 = mpbdx**2 * d(j) + mpbdx * d(j+1) + d(j+2)
if (j <= mpnw + 2) t1 = t1 + mprdx * d(j+3)
t0 = anint (rb * t1)
! Split trial quotient term into high and low order halves, 24-bits each.
a1 = mpb24x * aint (mpr24x * t0)
a2 = t0 - a1
! Multiply trial quotient term by each term of divisor, saving high and low-
! order parts into appropriate terms of the running dividend.
dc = a1 * b2 + a2 * b1
c1 = mpbdx * aint (mprdx * dc)
c2 = dc - c1
d(j) = d(j) - mprdx * (a1 * b1 + c1)
d(j+1) = d(j+1) - a2 * b2 - c2
d(j+1) = d(j+1) + mpbdx * d(j)
d(j) = t0
! Continue computing quotient terms past the end of the input dividend, until
! trial running dividend is zero.
if (j >= na + 2) then
if (j <= mpnw + 1) d(j+4) = 0.d0
endif
enddo
! Final bookkeeping.
j = mpnw + 3
170 continue
d(j+1) = 0.d0
if (d(2) == 0.d0) then
is = 1
else
is = 2
endif
nc = min (j - 1, mpnw)
d(nc+4) = 0.d0
d(nc+5) = 0.d0
do i = j + 1, 3, -1
d(i+1) = d(i-is+1)
enddo