-
Notifications
You must be signed in to change notification settings - Fork 0
/
rksuite.f90
3149 lines (3087 loc) · 112 KB
/
rksuite.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
module rksuite_90_prec
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
integer, parameter :: wp = selected_real_kind(10,50)
end module rksuite_90_prec
module rksuite_90
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
use rksuite_90_prec, only:wp
implicit none
private
public :: wp, setup, range_integrate, step_integrate, interpolate, &
global_error, statistics, reset_t_end, collect_garbage, &
set_stop_on_fatal, get_saved_fatal
!starta!
public :: rk_comm_real_1d
type rk_comm_real_1d
private
real(kind=wp) :: t, t_old, t_start, t_end, dir !indep!
real(kind=wp) :: h, h_old, h_start, h_average !indep!
real(kind=wp) :: tol
integer :: f_count, full_f_count, step_count, bad_step_count
logical :: at_t_start, at_t_end
real(kind=wp), dimension(:), pointer :: thresh, weights, ymax !shp-dep!
real(kind=wp), dimension(:), pointer :: scratch, y, yp, y_new !dep!
real(kind=wp), dimension(:), pointer :: y_old, yp_old, v0, v1 !dep!
real(kind=wp), dimension(:), pointer :: err_estimates, v2, v3 !dep!
real(kind=wp), dimension(:), pointer :: vtemp !dep!
real(kind=wp), dimension(:,:), pointer :: stages !dep!
real(kind=wp) :: a(13,13), b(13), c(13), bhat(13), r(11,6), e(7)
integer :: ptr(13), no_of_stages, rk_method, intrp_degree
logical :: intrp_able, intrp_needs_stages
real(kind=wp) :: toosml, cost, safety, expon, stability_radius, tan_angle, &
rs, rs1, rs2, rs3, rs4
integer :: order, last_stage, max_stiff_iters, no_of_ge_steps
logical :: fsal
real(kind=wp) :: ge_max_contrib
real(kind=wp) :: t_ge_max_contrib !indep!
integer :: ge_f_count
real(kind=wp), dimension(:), pointer :: ge_assess !shp-dep!
real(kind=wp), dimension(:), pointer :: ge_y, ge_yp, ge_y_new !dep!
real(kind=wp), dimension(:), pointer :: ge_err_estimates !dep!
real(kind=wp), dimension(:,:), pointer :: ge_stages !dep!
logical :: erason, erasfl
real(kind=wp) :: mcheps, dwarf, round_off, sqrrmc, cubrmc, sqtiny
integer :: outch
logical :: print_message, use_range
character(len=80) :: rec(10)
real(kind=wp) :: tlast, range_t_end !indep!
real(kind=wp), dimension(:), pointer :: xstage, ytemp !dep!
real(kind=wp), dimension(:,:), pointer :: p !dep!
integer :: stiff_bad_step_count, hit_t_end_count
real(kind=wp) :: errold
logical :: chkeff, phase2
integer, dimension(7) :: save_states
logical :: stop_on_fatal, saved_fatal_err
end type rk_comm_real_1d
!enda!
interface setup
module procedure setup_r1
end interface
interface range_integrate
module procedure range_integrate_r1
end interface
interface step_integrate
module procedure step_integrate_r1
end interface
interface statistics
module procedure statistics_r1
end interface
interface global_error
module procedure global_error_r1
end interface
interface reset_t_end
module procedure reset_t_end_r1
end interface
interface interpolate
module procedure interpolate_r1
end interface
interface set_stop_on_fatal
module procedure set_stop_on_fatal_r1
end interface
interface get_saved_fatal
module procedure get_saved_fatal_r1
end interface
interface collect_garbage
module procedure collect_garbage_r1
end interface
contains
!startb!
subroutine machine_const(round_off,sqrrmc,cubrmc,sqtiny,outch)
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
real(kind=wp), intent(out) :: round_off, sqrrmc, cubrmc, sqtiny
integer, intent(out) :: outch
!
real(kind=wp) :: dummy
real(kind=wp), parameter :: third=1.0_wp/3.0_wp, ten=10.0_wp
!
outch = 6
!
round_off = ten*epsilon(dummy)
sqrrmc = sqrt(epsilon(dummy))
cubrmc = epsilon(dummy)**third
sqtiny = sqrt(tiny(dummy))
!
end subroutine machine_const
subroutine method_const(rk_method, a, b, c, bhat, r, e, ptr, no_of_stages, &
intrp_degree, intrp_able, intrp_needs_stages, &
cost, safety, expon, stability_radius, &
tan_angle, rs, rs1, rs2, rs3, rs4, order, last_stage, &
max_stiff_iters, no_of_ge_steps, fsal, cdiff)
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
integer, intent(in) :: rk_method
real(kind=wp), intent(out) :: a(13,13), b(13), c(13), bhat(13), r(11,6), e(7)
integer, intent(out) :: ptr(13), no_of_stages, intrp_degree
logical, intent(out) :: intrp_able, intrp_needs_stages
real(kind=wp), intent(out) :: cost, safety, expon, stability_radius, &
tan_angle, rs, rs1, rs2, rs3, rs4, cdiff
integer, intent(out) :: order, last_stage, max_stiff_iters, no_of_ge_steps
logical, intent(out) :: fsal
!
integer :: i
real(kind=wp), parameter :: fivepc=0.05_wp, one=1.0_wp, two=2.0_wp, &
fifty=50.0_wp
!
select case (rk_method)
case(1)
!
! METHD = 1.
! This pair is from "A 3(2) Pair of Runge-Kutta Formulas" by P. Bogacki
! and L.F. Shampine, Appl. Math. Lett., 2, pp. 321-325, 1989. The authors
! are grateful to P. Bogacki for his assistance in implementing the pair.
!
no_of_stages = 4; fsal = .true.; order = 2
tan_angle = 8.9_wp; stability_radius = 2.3_wp
safety = 0.8_wp; intrp_able = .true.; intrp_degree = 3
intrp_needs_stages = .false.; no_of_ge_steps = 3
!
ptr(1:4) = (/ 0, 1, 2, 3 /)
!
a(2,1) = 1.0_wp/2.0_wp
a(3,1) = 0.0_wp
a(3,2) = 3.0_wp/4.0_wp
a(4,1) = 2.0_wp/9.0_wp
a(4,2) = 1.0_wp/3.0_wp
a(4,3) = 4.0_wp/9.0_wp
!
! The coefficients BHAT refer to the formula used to advance the
! integration, here the one of order 3. The coefficients B refer
! to the other formula, here the one of order 2. For this pair, BHAT
! is not needed since FSAL = .TRUE.
!
b(1) = 7.0_wp/24.0_wp
b(2) = 1.0_wp/4.0_wp
b(3) = 1.0_wp/3.0_wp
b(4) = 1.0_wp/8.0_wp
!
c(1) = 0.0_wp
c(2) = 1.0_wp/2.0_wp
c(3) = 3.0_wp/4.0_wp
c(4) = 1.0_wp
!
case (2)
!
! METHD = 2
! This pair is from "An Efficient Runge-Kutta (4,5) Pair" by P. Bogacki
! and L.F. Shampine, Rept. 89-20, Math. Dept., Southern Methodist
! University, Dallas, Texas, USA, 1989. The authors are grateful to
! P. Bogacki for his assistance in implementing the pair. Shampine and
! Bogacki subsequently modified the formula to enhance the reliability of
! the pair. The original fourth order formula is used in an estimate of
! the local error. If the step fails, the computation is broken off. If
! the step is acceptable, the first evaluation of the next step is done,
! i.e., the pair is implemented as FSAL and the local error of the step
! is again estimated with a fourth order formula using the additional data.
! The step must succeed with both estimators to be accepted. When the
! second estimate is formed, it is used for the subsequent adjustment of
! the step size because it is of higher quality. The two fourth order
! formulas are well matched to leading order, and only exceptionally do
! the estimators disagree -- problems with discontinuous coefficients are
! handled more reliably by using two estimators as is global error
! estimation.
!
no_of_stages = 8; fsal = .true.; order = 4
tan_angle = 5.2_wp; stability_radius = 3.9_wp
safety = 0.8_wp; intrp_able = .true.
intrp_needs_stages = .true.; intrp_degree = 6
no_of_ge_steps = 2
!
ptr(1:8) = (/ 0, 1, 2, 3, 4, 5, 6, 7 /)
!
a(2,1) = 1.0_wp/6.0_wp
a(3,1) = 2.0_wp/27.0_wp
a(3,2) = 4.0_wp/27.0_wp
a(4,1) = 183.0_wp/1372.0_wp
a(4,2) = -162.0_wp/343.0_wp
a(4,3) = 1053.0_wp/1372.0_wp
a(5,1) = 68.0_wp/297.0_wp
a(5,2) = -4.0_wp/11.0_wp
a(5,3) = 42.0_wp/143.0_wp
a(5,4) = 1960.0_wp/3861.0_wp
a(6,1) = 597.0_wp/22528.0_wp
a(6,2) = 81.0_wp/352.0_wp
a(6,3) = 63099.0_wp/585728.0_wp
a(6,4) = 58653.0_wp/366080.0_wp
a(6,5) = 4617.0_wp/20480.0_wp
a(7,1) = 174197.0_wp/959244.0_wp
a(7,2) = -30942.0_wp/79937.0_wp
a(7,3) = 8152137.0_wp/19744439.0_wp
a(7,4) = 666106.0_wp/1039181.0_wp
a(7,5) = -29421.0_wp/29068.0_wp
a(7,6) = 482048.0_wp/414219.0_wp
a(8,1) = 587.0_wp/8064.0_wp
a(8,2) = 0.0_wp
a(8,3) = 4440339.0_wp/15491840.0_wp
a(8,4) = 24353.0_wp/124800.0_wp
a(8,5) = 387.0_wp/44800.0_wp
a(8,6) = 2152.0_wp/5985.0_wp
a(8,7) = 7267.0_wp/94080.0_wp
!
! The coefficients B refer to the formula of order 4.
!
b(1) = 2479.0_wp/34992.0_wp
b(2) = 0.0_wp
b(3) = 123.0_wp/416.0_wp
b(4) = 612941.0_wp/3411720.0_wp
b(5) = 43.0_wp/1440.0_wp
b(6) = 2272.0_wp/6561.0_wp
b(7) = 79937.0_wp/1113912.0_wp
b(8) = 3293.0_wp/556956.0_wp
!
! The coefficients E refer to an estimate of the local error based on
! the first formula of order 4. It is the difference of the fifth order
! result, here located in A(8,:), and the fourth order result. By
! construction both E(2) and E(7) are zero.
!
e(1) = -3.0_wp/1280.0_wp
e(2) = 0.0_wp
e(3) = 6561.0_wp/632320.0_wp
e(4) = -343.0_wp/20800.0_wp
e(5) = 243.0_wp/12800.0_wp
e(6) = -1.0_wp/95.0_wp
e(7) = 0.0_wp
!
c(1) = 0.0_wp
c(2) = 1.0_wp/6.0_wp
c(3) = 2.0_wp/9.0_wp
c(4) = 3.0_wp/7.0_wp
c(5) = 2.0_wp/3.0_wp
c(6) = 3.0_wp/4.0_wp
c(7) = 1.0_wp
c(8) = 1.0_wp
!
! To do interpolation with this pair, some extra stages have to be computed.
! The following additional A and C coefficients are for this purpose.
! In addition there is an array R that plays a role for interpolation
! analogous to that of BHAT for the basic step.
!
c(9) = 1.0_wp/2.0_wp
c(10) = 5.0_wp/6.0_wp
c(11) = 1.0_wp/9.0_wp
!
a(9,1) = 455.0_wp/6144.0_wp
a(10,1) = -837888343715.0_wp/13176988637184.0_wp
a(11,1) = 98719073263.0_wp/1551965184000.0_wp
a(9,2) = 0.0_wp
a(10,2) = 30409415.0_wp/52955362.0_wp
a(11,2) = 1307.0_wp/123552.0_wp
a(9,3) = 10256301.0_wp/35409920.0_wp
a(10,3) = -48321525963.0_wp/759168069632.0_wp
a(11,3) = 4632066559387.0_wp/70181753241600.0_wp
a(9,4) = 2307361.0_wp/17971200.0_wp
a(10,4) = 8530738453321.0_wp/197654829557760.0_wp
a(11,4) = 7828594302389.0_wp/382182512025600.0_wp
a(9,5) = -387.0_wp/102400.0_wp
a(10,5) = 1361640523001.0_wp/1626788720640.0_wp
a(11,5) = 40763687.0_wp/11070259200.0_wp
a(9,6) = 73.0_wp/5130.0_wp
a(10,6) = -13143060689.0_wp/38604458898.0_wp
a(11,6) = 34872732407.0_wp/224610586200.0_wp
a(9,7) = -7267.0_wp/215040.0_wp
a(10,7) = 18700221969.0_wp/379584034816.0_wp
a(11,7) = -2561897.0_wp/30105600.0_wp
a(9,8) = 1.0_wp/32.0_wp
a(10,8) = -5831595.0_wp/847285792.0_wp
a(11,8) = 1.0_wp/10.0_wp
a(10,9) = -5183640.0_wp/26477681.0_wp
a(11,9) = -1.0_wp/10.0_wp
a(11,10) = -1403317093.0_wp/11371610250.0_wp
!
r(1:11,1) = 0.0_wp; r(2,1:6) = 0.0_wp
r(1,6) = -12134338393.0_wp/1050809760.0_wp
r(1,5) = -1620741229.0_wp/50038560.0_wp
r(1,4) = -2048058893.0_wp/59875200.0_wp
r(1,3) = -87098480009.0_wp/5254048800.0_wp
r(1,2) = -11513270273.0_wp/3502699200.0_wp
!
r(3,6) = -33197340367.0_wp/1218433216.0_wp
r(3,5) = -539868024987.0_wp/6092166080.0_wp
r(3,4) = -39991188681.0_wp/374902528.0_wp
r(3,3) = -69509738227.0_wp/1218433216.0_wp
r(3,2) = -29327744613.0_wp/2436866432.0_wp
!
r(4,6) = -284800997201.0_wp/19905339168.0_wp
r(4,5) = -7896875450471.0_wp/165877826400.0_wp
r(4,4) = -333945812879.0_wp/5671036800.0_wp
r(4,3) = -16209923456237.0_wp/497633479200.0_wp
r(4,2) = -2382590741699.0_wp/331755652800.0_wp
!
r(5,6) = -540919.0_wp/741312.0_wp
r(5,5) = -103626067.0_wp/43243200.0_wp
r(5,4) = -633779.0_wp/211200.0_wp
r(5,3) = -32406787.0_wp/18532800.0_wp
r(5,2) = -36591193.0_wp/86486400.0_wp
!
r(6,6) = 7157998304.0_wp/374350977.0_wp
r(6,5) = 30405842464.0_wp/623918295.0_wp
r(6,4) = 183022264.0_wp/5332635.0_wp
r(6,3) = -3357024032.0_wp/1871754885.0_wp
r(6,2) = -611586736.0_wp/89131185.0_wp
!
r(7,6) = -138073.0_wp/9408.0_wp
r(7,5) = -719433.0_wp/15680.0_wp
r(7,4) = -1620541.0_wp/31360.0_wp
r(7,3) = -385151.0_wp/15680.0_wp
r(7,2) = -65403.0_wp/15680.0_wp
!
r(8,6) = 1245.0_wp/64.0_wp
r(8,5) = 3991.0_wp/64.0_wp
r(8,4) = 4715.0_wp/64.0_wp
r(8,3) = 2501.0_wp/64.0_wp
r(8,2) = 149.0_wp/16.0_wp
r(8,1) = 1.0_wp
!
r(9,6) = 55.0_wp/3.0_wp
r(9,5) = 71.0_wp
r(9,4) = 103.0_wp
r(9,3) = 199.0_wp/3.0_wp
r(9,2) = 16.0d0
!
r(10,6) = -1774004627.0_wp/75810735.0_wp
r(10,5) = -1774004627.0_wp/25270245.0_wp
r(10,4) = -26477681.0_wp/359975.0_wp
r(10,3) = -11411880511.0_wp/379053675.0_wp
r(10,2) = -423642896.0_wp/126351225.0_wp
!
r(11,6) = 35.0_wp
r(11,5) = 105.0_wp
r(11,4) = 117.0_wp
r(11,3) = 59.0_wp
r(11,2) = 12.0_wp
!
case (3)
!
! METHD = 3
! This pair is from "High Order Embedded Runge-Kutta Formulae" by P.J.
! Prince and J.R. Dormand, J. Comp. Appl. Math.,7, pp. 67-75, 1981. The
! authors are grateful to P. Prince and J. Dormand for their assistance in
! implementing the pair.
!
no_of_stages = 13; fsal = .false.; order = 7
tan_angle = 11.0_wp; stability_radius = 5.2_wp
safety = 0.8_wp; intrp_able = .false.
intrp_needs_stages = .false.; intrp_degree = 0
no_of_ge_steps = 2
!
ptr(1:13) = (/ 0, 1, 2, 1, 3, 2, 4, 5, 6, 7, 8, 9, 1 /)
!
a(2,1) = 5.55555555555555555555555555556e-2_wp
a(3,1) = 2.08333333333333333333333333333e-2_wp
a(3,2) = 6.25e-2_wp
a(4,1) = 3.125e-2_wp
a(4,2) = 0.0_wp
a(4,3) = 9.375e-2_wp
a(5,1) = 3.125e-1_wp
a(5,2) = 0.0_wp
a(5,3) = -1.171875_wp
a(5,4) = 1.171875_wp
a(6,1) = 3.75e-2_wp
a(6,2) = 0.0_wp
a(6,3) = 0.0_wp
a(6,4) = 1.875e-1_wp
a(6,5) = 1.5e-1_wp
a(7,1) = 4.79101371111111111111111111111e-2_wp
a(7,2) = 0.0_wp
a(7,3) = 0.0_wp
a(7,4) = 1.12248712777777777777777777778e-1_wp
a(7,5) = -2.55056737777777777777777777778e-2_wp
a(7,6) = 1.28468238888888888888888888889e-2_wp
a(8,1) = 1.6917989787292281181431107136e-2_wp
a(8,2) = 0.0_wp
a(8,3) = 0.0_wp
a(8,4) = 3.87848278486043169526545744159e-1_wp
a(8,5) = 3.59773698515003278967008896348e-2_wp
a(8,6) = 1.96970214215666060156715256072e-1_wp
a(8,7) = -1.72713852340501838761392997002e-1_wp
a(9,1) = 6.90957533591923006485645489846e-2_wp
a(9,2) = 0.0_wp
a(9,3) = 0.0_wp
a(9,4) = -6.34247976728854151882807874972e-1_wp
a(9,5) = -1.61197575224604080366876923982e-1_wp
a(9,6) = 1.38650309458825255419866950133e-1_wp
a(9,7) = 9.4092861403575626972423968413e-1_wp
a(9,8) = 2.11636326481943981855372117132e-1_wp
a(10,1) = 1.83556996839045385489806023537e-1_wp
a(10,2) = 0.0_wp
a(10,3) = 0.0_wp
a(10,4) = -2.46876808431559245274431575997_wp
a(10,5) = -2.91286887816300456388002572804e-1_wp
a(10,6) = -2.6473020233117375688439799466e-2_wp
a(10,7) = 2.84783876419280044916451825422_wp
a(10,8) = 2.81387331469849792539403641827e-1_wp
a(10,9) = 1.23744899863314657627030212664e-1_wp
a(11,1) = -1.21542481739588805916051052503_wp
a(11,2) = 0.0_wp
a(11,3) = 0.0_wp
a(11,4) = 1.66726086659457724322804132886e1_wp
a(11,5) = 9.15741828416817960595718650451e-1_wp
a(11,6) = -6.05660580435747094755450554309_wp
a(11,7) = -1.60035735941561781118417064101e1_wp
a(11,8) = 1.4849303086297662557545391898e1_wp
a(11,9) = -1.33715757352898493182930413962e1_wp
a(11,10) = 5.13418264817963793317325361166_wp
a(12,1) = 2.58860916438264283815730932232e-1_wp
a(12,2) = 0.0_wp
a(12,3) = 0.0_wp
a(12,4) = -4.77448578548920511231011750971_wp
a(12,5) = -4.3509301377703250944070041181e-1_wp
a(12,6) = -3.04948333207224150956051286631_wp
a(12,7) = 5.57792003993609911742367663447_wp
a(12,8) = 6.15583158986104009733868912669_wp
a(12,9) = -5.06210458673693837007740643391_wp
a(12,10) = 2.19392617318067906127491429047_wp
a(12,11) = 1.34627998659334941535726237887e-1_wp
a(13,1) = 8.22427599626507477963168204773e-1_wp
a(13,2) = 0.0_wp
a(13,3) = 0.0_wp
a(13,4) = -1.16586732572776642839765530355e1_wp
a(13,5) = -7.57622116690936195881116154088e-1_wp
a(13,6) = 7.13973588159581527978269282765e-1_wp
a(13,7) = 1.20757749868900567395661704486e1_wp
a(13,8) = -2.12765911392040265639082085897_wp
a(13,9) = 1.99016620704895541832807169835_wp
a(13,10) = -2.34286471544040292660294691857e-1_wp
a(13,11) = 1.7589857770794226507310510589e-1_wp
a(13,12) = 0.0_wp
!
! The coefficients BHAT refer to the formula used to advance the
! integration, here the one of order 8. The coefficients B refer
! to the other formula, here the one of order 7.
!
bhat(1) = 4.17474911415302462220859284685e-2_wp
bhat(2) = 0.0_wp
bhat(3) = 0.0_wp
bhat(4) = 0.0_wp
bhat(5) = 0.0_wp
bhat(6) = -5.54523286112393089615218946547e-2_wp
bhat(7) = 2.39312807201180097046747354249e-1_wp
bhat(8) = 7.0351066940344302305804641089e-1_wp
bhat(9) = -7.59759613814460929884487677085e-1_wp
bhat(10) = 6.60563030922286341461378594838e-1_wp
bhat(11) = 1.58187482510123335529614838601e-1_wp
bhat(12) = -2.38109538752862804471863555306e-1_wp
bhat(13) = 2.5e-1_wp
!
b(1) = 2.9553213676353496981964883112e-2_wp
b(2) = 0.0_wp
b(3) = 0.0_wp
b(4) = 0.0_wp
b(5) = 0.0_wp
b(6) = -8.28606276487797039766805612689e-1_wp
b(7) = 3.11240900051118327929913751627e-1_wp
b(8) = 2.46734519059988698196468570407_wp
b(9) = -2.54694165184190873912738007542_wp
b(10) = 1.44354858367677524030187495069_wp
b(11) = 7.94155958811272872713019541622e-2_wp
b(12) = 4.44444444444444444444444444445e-2_wp
b(13) = 0.0_wp
!
c(1) = 0.0_wp
c(2) = 5.55555555555555555555555555556e-2_wp
c(3) = 8.33333333333333333333333333334e-2_wp
c(4) = 1.25e-1_wp
c(5) = 3.125e-1_wp
c(6) = 3.75e-1_wp
c(7) = 1.475e-1_wp
c(8) = 4.65e-1_wp
c(9) = 5.64865451382259575398358501426e-1_wp
c(10) = 6.5e-1_wp
c(11) = 9.24656277640504446745013574318e-1_wp
c(12) = 1.0_wp
c(13) = c(12)
!
end select
!
! The definitions of all pairs come here for the calculation of
! LAST_STAGE - the position of the last evaluated stage in a method
! RS1, RS2, RS3, RS4 - minimum and maximum rations used is step selection
! COST - the cost of a step
! MAX_STIFF_ITERS - the number of iterations permitted in stiffness detection
! There are at most Q = 3 function calls per iteration. MAX_STIFF_ITERS
! is determined so that Q*MAX_STIFF_ITERS <= 5% of the cost of
! 50 steps and 1 <= MAX_STIFF_ITERS <= 8. This limits the number of
! function calls in each diagnosis of stiffness to 24.
! EXPON - an exponent for use in step slection
! CDIFF - a coefficent used in determining the minimum permissible step
!
last_stage = ptr(no_of_stages)
if (fsal) then
cost = real(no_of_stages-1,kind=wp)
else
cost = real(no_of_stages,kind=wp)
end if
!
max_stiff_iters = min(8,max(1,int(fivepc*cost*fifty)))
!
expon = one/(order+one)
!
! In calculating CDIFF it is assumed that there will be a non-zero
! difference |C(I) - C(J)| less than one. If C(I) = C(J) for any I not
! equal to J, they should be made precisely equal by assignment.
!
cdiff = one
do i = 1, no_of_stages - 1
cdiff = min( cdiff, minval( &
abs((c(i)-c(i+1:no_of_stages))), &
mask=(c(i)-c(i+1:no_of_stages)/=0)) )
end do
!
rs = two; rs1 = one/rs; rs2 = rs**2
rs3 = rs*rs2; rs4 = one/rs3
!
end subroutine method_const
!endb!
!startc!
subroutine setup_r1(comm,t_start,y_start,t_end,tolerance,thresholds, &
method,task,error_assess,h_start,message)
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
real(kind=wp), intent(in) :: t_end, t_start !indep!
real(kind=wp), intent(in) :: tolerance
real(kind=wp), dimension(:), intent(in) :: y_start !dep!
real(kind=wp), dimension(:), intent(in) :: thresholds !shp-dep!
type(rk_comm_real_1d) :: comm
real(kind=wp), intent(in), optional :: h_start !indep!
logical, intent(in), optional :: error_assess, message
character(len=*), intent(in), optional :: task, method
!
character(len=*), parameter :: srname="SETUP"
!
real(kind=wp) :: hmin !indep!
real(kind=wp) :: cdiff
integer :: ier, nrec, tr_dim_of_stages
logical :: legalt
character(len=1) :: task1, method1
!
integer, parameter :: not_ready=-1, fatal=911, just_fine=1
real(kind=wp), parameter :: zero=0.0_wp, pt01=0.01_wp, fivepc=0.05_wp, &
third=1.0_wp/3.0_wp, one=1.0_wp, two=2.0_wp, ten=10.0_wp, fifty=50.0_wp
!
ier = just_fine; nrec = 0
!
! Clear previous state of the suite.
!
call setup_global_stuff
nullify(comm%thresh, comm%err_estimates, comm%weights, comm%y_old, &
comm%scratch, &
comm%y, comm%yp, comm%y_new, comm%yp_old, comm%stages, comm%ge_y, &
comm%ge_yp, comm%ge_err_estimates, comm%ge_assess, comm%ge_y_new, &
comm%ge_stages, comm%v0, comm%v1, comm%v2, comm%v3, comm%vtemp, &
comm%xstage, comm%ytemp, comm%p)
!
! Fetch output channel and machine constants;
!
call machine_const(comm%round_off,comm%sqrrmc,comm%cubrmc,comm%sqtiny, &
comm%outch)
!
body: do
!
! Check for valid shape
if (size(shape(y_start))>0) then
if (any(shape(y_start)==0)) then
ier = fatal; nrec = 2; write (comm%rec,"(a)") &
" ** An extent of Y_START has zero length. This is not permitted."
exit body
end if
end if
!
! Check and process non-trivial optional arguments
if (present(task)) then
task1 = task(1:1); comm%use_range = task1 == "R" .or. task1 == "r"
legalt = comm%use_range .or. task1 == "S" .or. task1 == "s"
if (.not.legalt) then
ier = fatal; nrec = 2; write (comm%rec,"(a,a,a/a)") &
" ** You have set the first character of TASK to be '",TASK1,"'.", &
" ** It must be one of 'R','r','S' or 's'."
exit body
end if
end if
if (present(method)) then
method1 = method(1:1)
select case (method1)
case("L","l"); comm%rk_method = 1
case("M","m"); comm%rk_method = 2
case("H","h"); comm%rk_method = 3
case default
ier = fatal; nrec = 2; write (comm%rec,"(a,a,a/a)") &
" ** You have set the first character of METHOD to be '",METHOD1,"'.", &
" ** It must be one of 'L','l','M','m','H' or 'h'."
exit body
end select
end if
if (present(message)) comm%print_message = message
!
! Check consistency of array arguments
!
if (any(shape(y_start)/=shape(thresholds))) then
ier = fatal; nrec = 1; write (comm%rec,"(a)") &
" ** The shapes of Y_START and THRESHOLDS are not consistent."
exit body
end if
!
! Check and process compulsory arguments
if (t_start == t_end) then
ier = fatal; nrec = 1; write (comm%rec,"(a,e13.5,a)") &
" ** You have set T_START = T_END = ",T_START,"."
exit body
else
comm%t_end = t_end; comm%t_start = t_start
comm%t_old = t_start; comm%t = t_start
comm%dir = sign(one,t_end-t_start)
end if
if ((tolerance > pt01) .or. (tolerance < comm%round_off)) then
ier = fatal; nrec = 2; write (comm%rec,"(a,e13.5,a/a,e13.5,a)") &
" ** You have set TOLERANCE = ",tolerance," which is not permitted. The", &
" ** range of permitted values is (",comm%round_off,",0.01)."
exit body
else
comm%tol = tolerance
end if
if (minval(thresholds) < comm%sqtiny) then !spec-ar!
ier = fatal; nrec = 2; write (comm%rec,"(a/a,e13.5,a)") &
" ** You have set a component of THRESHOLDS to be less than the permitted", &
" ** minimum,",comm%sqtiny,"."
exit body
end if
!
! Set formula definitions and characteristics
call method_const(comm%rk_method, comm%a, comm%b, comm%c, comm%bhat, &
comm%r, comm%e, comm%ptr, comm%no_of_stages, comm%intrp_degree, &
comm%intrp_able, comm%intrp_needs_stages, comm%cost, &
comm%safety, comm%expon, comm%stability_radius, comm%tan_angle, &
comm%rs, comm%rs1, comm%rs2, comm%rs3, comm%rs4, comm%order, &
comm%last_stage, comm%max_stiff_iters, comm%no_of_ge_steps, comm%fsal,&
cdiff)
!
tr_dim_of_stages = maxval(comm%ptr(2:comm%no_of_stages))
comm%toosml = comm%round_off/cdiff
!
! In STEP_INTEGRATE the first step taken will have magnitude H. If
! H_START = ABS(H_START) is not equal to zero, H = H_START. If H_START is
! equal to zero, the code is to find an on-scale initial step size H. To
! start this process, H is set here to an upper bound on the first step
! size that reflects the scale of the independent variable.
! RANGE_INTEGRATE has some additional information, namely the first output
! point, that is used to refine this bound in STEP_INTEGRATE when called
! from RANGE_INTEGRATE. If H_START is not zero, but it is either too big
! or too small, the input H_START is ignored and H_START is set to zero to
! activate the automatic determination of an on-scale initial step size.
!
hmin = max(comm%sqtiny,comm%toosml*max(abs(t_start),abs(t_end)))
if (abs(t_end-t_start) < hmin) then
ier = fatal; nrec = 4; write (comm%rec,"(a/a/a,e13.5,a/a,e13.5,a)") &
" ** You have set values for T_END and T_START that are not clearly", &
" ** distinguishable for the method and the precision of the computer", &
" ** being used. ABS(T_END-T_START) is ",ABS(T_END-T_START)," but should be", &
" ** at least ",hmin,"."
exit body
end if
if (present(h_start)) comm%h_start = abs(h_start)
if (comm%h_start > abs(t_end-t_start) .or. comm%h_start < hmin) &
comm%h_start = zero
if (comm%h_start == zero) then
comm%h = max(abs(t_end-t_start)/comm%rs3,hmin)
else
comm%h = comm%h_start
end if
!
! Allocate a number of arrays using pointers.
!
allocate(comm%thresh(size(y_start,1)), & !alloc!
comm%err_estimates(size(y_start,1)), & !alloc!
comm%weights(size(y_start,1)), & !alloc!
comm%y_old(size(y_start,1)), & !alloc!
comm%scratch(size(y_start,1)), & !alloc!
comm%y(size(y_start,1)), & !alloc!
comm%yp(size(y_start,1)), & !alloc!
comm%stages(size(y_start,1),tr_dim_of_stages), & !alloc!
comm%ymax(size(y_start,1)),stat=ier) !alloc!
if (ier /= 0) then
ier = fatal; nrec = 1 ; write (comm%rec,"(a)") &
" ** Not enough storage available to create workspace required internally."
exit body
else
comm%y = y_start; comm%ymax = abs(y_start)
comm%thresh = thresholds;
comm%y_new => comm%scratch; comm%yp_old => comm%scratch
comm%v0 => comm%err_estimates; comm%vtemp => comm%scratch
comm%v1 => comm%stages(:,1); comm%v2 => comm%stages(:,2)
comm%v3 => comm%stages(:,3)
end if
!
! Pre-allocate storage for interpolation if the TASK = `R' was specified.
!
if (comm%use_range) then
if (comm%intrp_able) then
if (comm%rk_method==1) then
comm%p => comm%stages(:,1:2)
else if (comm%rk_method==2) then
allocate(comm%p(size(y_start,1),5), & !alloc!
comm%ytemp(size(y_start,1)), & !alloc!
comm%xstage(size(y_start,1)),stat=ier) !alloc!
if (ier /= 0) then
ier = fatal; nrec = 1 ; write (comm%rec,"(a)") &
" ** Not enough storage available to create workspace required internally."
exit body
end if
end if
end if
end if
!
! Initialise state and allocate storage for global error assessment
!
comm%t_ge_max_contrib = t_start
if (present(error_assess)) comm%erason = error_assess
if (comm%erason) then
!
! Storage is required for the stages of a secondary integration. The
! stages of the primary intergration can only be overwritten in the
! cases where there is no interpolant or the interpolant does not
! require information about the stages (e.g. METHOD 'H' and METHOD 'L',
! respectively).
if (.not.comm%intrp_needs_stages) then
comm%ge_stages => comm%stages
else
allocate(comm%ge_stages(size(y_start,1),tr_dim_of_stages),stat=ier) !alloc!
if (ier /= 0) then
ier = fatal; nrec = 1 ; write (comm%rec,"(a)") &
" ** Not enough storage available to create workspace required internally."
exit body
end if
end if
allocate(comm%ge_y(size(y_start,1)), & !alloc!
comm%ge_yp(size(y_start,1)), & !alloc!
comm%ge_err_estimates(size(y_start,1)), & !alloc!
comm%ge_assess(size(y_start,1)), & !alloc!
comm%ge_y_new(size(y_start,1)),stat=ier) !alloc!
if (ier /= 0) then
ier = fatal; nrec = 1 ; write (comm%rec,"(a)") &
" ** Not enough storage available to create workspace required internally."
exit body
else
comm%ge_assess = 0.0_wp; comm%ge_y = y_start
end if
end if
exit body
end do body
!
call rkmsg_r1(ier,srname,nrec,comm)
!
contains
subroutine setup_global_stuff
!
comm%h_start=0.0_wp; comm%h_old=0.0_wp
comm%f_count=0; comm%full_f_count=0; comm%step_count=0; comm%bad_step_count=0
comm%at_t_start=.true.; comm%at_t_end=.false.
comm%rk_method=2;
comm%ge_max_contrib=0.0_wp; comm%ge_f_count=0
comm%erason=.false.; comm%erasfl=.false.
comm%print_message=.true.; comm%use_range=.true.
comm%stiff_bad_step_count=0; comm%hit_t_end_count=0
comm%errold=0.0_wp; comm%h_average=0.0_wp
comm%chkeff=.false.; comm%phase2=.true.
comm%save_states(1:7)=not_ready
comm%stop_on_fatal=.true.; comm%saved_fatal_err=.false.
!
end subroutine setup_global_stuff
end subroutine setup_r1
subroutine collect_garbage_r1(comm)
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
! Modified by I.Gladwell (Aug 2002)
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
type(rk_comm_real_1d) :: comm
!
if (associated(comm%thresh)) then
deallocate(comm%thresh); nullify(comm%thresh)
end if
if (associated(comm%y)) then
deallocate(comm%y); nullify(comm%y)
end if
if (associated(comm%yp)) then
deallocate(comm%yp); nullify(comm%yp)
end if
if (associated(comm%ymax)) then
deallocate(comm%ymax); nullify(comm%ymax)
end if
if (associated(comm%scratch)) then
deallocate(comm%scratch); nullify(comm%scratch)
nullify(comm%y_new); nullify(comm%yp_old); nullify(comm%vtemp)
end if
if (associated(comm%weights)) then
deallocate(comm%weights); nullify(comm%weights)
end if
if (associated(comm%ytemp)) then
deallocate(comm%ytemp); nullify(comm%ytemp)
end if
if (associated(comm%y_old)) then
deallocate(comm%y_old); nullify(comm%y_old)
end if
if (associated(comm%err_estimates)) then
deallocate(comm%err_estimates); nullify(comm%err_estimates)
nullify(comm%v0)
end if
if (associated(comm%p,comm%stages(:,1:2))) then
nullify(comm%p)
end if
if (associated(comm%ge_stages,comm%stages)) then
deallocate(comm%stages); nullify(comm%stages); nullify(comm%ge_stages);
nullify(comm%v1,comm%v2,comm%v3)
else if (associated(comm%ge_stages)) then
deallocate(comm%ge_stages); nullify(comm%ge_stages)
end if
if (associated(comm%ge_y_new)) then
deallocate(comm%ge_y_new); nullify(comm%ge_y_new)
end if
if (associated(comm%ge_assess)) then
deallocate(comm%ge_assess); nullify(comm%ge_assess)
end if
if (associated(comm%ge_err_estimates)) then
deallocate(comm%ge_err_estimates); nullify(comm%ge_err_estimates)
end if
if (associated(comm%ge_yp)) then
deallocate(comm%ge_yp); nullify(comm%ge_yp)
end if
if (associated(comm%ge_y)) then
deallocate(comm%ge_y); nullify(comm%ge_y)
end if
if (associated(comm%xstage)) then
deallocate(comm%xstage); nullify(comm%xstage)
end if
if (associated(comm%p)) then
deallocate(comm%p); nullify(comm%p)
end if
if (associated(comm%stages)) then
deallocate(comm%stages); nullify(comm%stages)
end if
!
end subroutine collect_garbage_r1
recursive subroutine range_integrate_r1(comm,f,t_want,t_got,y_got,yderiv_got, &
flag)
!
! Part of rksuite_90 v1.0 (Aug 1994)
! software for initial value problems in ODEs
!
! Authors: R.W. Brankin (NAG Ltd., Oxford, England)
! I. Gladwell (Math Dept., SMU, Dallas, TX, USA)
! see main doc for contact details
!
real(kind=wp), intent(in) :: t_want !indep!
real(kind=wp), intent(out) :: t_got !indep!
real(kind=wp), dimension(:), intent(out) :: y_got, yderiv_got !dep!
integer, intent(out), optional :: flag
type(rk_comm_real_1d), intent(inout) :: comm
!
interface
function f(t,y)
use rksuite_90_prec, only:wp
real(kind=wp), intent(in) :: t !indep!
real(kind=wp), dimension(:), intent(in) :: y !dep!
real(kind=wp), dimension(size(y,1)) :: f !dep!
end function f
end interface
!
character(len=*), parameter :: srname="RANGE_INTEGRATE"
!
real(kind=wp) :: hmin, t_now !indep!
integer :: step_flag, ier, nrec, state
logical :: goback, baderr
!
integer, parameter :: not_ready=-1, usable=-2, fatal=911, catastrophe=912, &
just_fine=1
logical, parameter :: tell=.false., ask=.true.
real(kind=wp), parameter :: zero=0.0_wp
!
ier = just_fine; nrec = 0
goback = .false.; baderr = .false.
body: do
!
! Is it permissible to call RANGE_INTEGRATE?
!
state = get_saved_state_r1("SETUP",comm%save_states)
if (state==fatal) then
ier = catastrophe; nrec = 1; write (comm%rec,"(a)") &
" ** A catastrophic error has already been detected elsewhere."
exit body
end if
if (state==not_ready) then
ier = fatal; nrec = 1; write (comm%rec,"(a)") &
" ** You have not called SETUP, so you cannot use RANGE_INTEGRATE."
exit body
end if
if (.not.comm%use_range) then
ier = fatal; nrec = 2; write (comm%rec,"(a/a)") &
" ** You have called RANGE_INTEGRATE after you specified in SETUP that you",&
" ** were going to use STEP_INTEGRATE. This is not permitted."
exit body
end if
state = get_saved_state_r1(srname,comm%save_states)
if (state==5 .or. state==6) then
ier = fatal; nrec = 1; write (comm%rec,"(a/a)") &
" ** This routine has already returned with a hard failure. You must call",&
" ** SETUP to start another problem."
exit body
end if
state = usable
call set_saved_state_r1(srname,state,comm)
!
if (comm%at_t_start) then
!
! First call.