-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathmo_ode_solver.f90
2916 lines (2278 loc) · 129 KB
/
mo_ode_solver.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
!> \file mo_ode_solver.f90
!> \brief This module provides a set of iterative methods for the approximation of solutions
!> of Ordinary Differential Equations (ODE).
!> \details
!> It includes the possibilities to integrate a system of Ordinary Differential Equations using Euler,
!> a fourth-order Runge-Kutta with fixed time-steps increments or a fourth-order Runge-Kutta with
!> adaptive stepsize control.
!> \authors Giovanni Dalmasso, Sebastian Mueller
!> \date Mar 2015
module mo_ode_solver
! This module provides a set of iterative methods for the approximation of solutions
! of Ordinary Differential Equations (ODE).
! Written Giovanni Dalmasso, Jul 2012
! Modified Giovanni Dalmasso, Mar 2013 - adapted to JAMS Fortran structure
! - collected together different methods
! - speeded up
! Apr 2013 - added documentation
! Modified Sebastian Mueller, Jan 2015 - added a parameter input for the derivatives to use the solver dynamically
! Mar 2015 - added RBstiff as a solver for stiff ODEs
! License
! -------
! This file is part of the JAMS Fortran package, distributed under the MIT License.
!
! Copyright (c) 2012-2015 Giovanni Dalmasso, Sebastian Mueller
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
! Note on Numerical Recipes License
! ---------------------------------
! Be aware that some code is under the Numerical Recipes License 3rd
! edition <http://www.nr.com/aboutNR3license.html>
!
! The Numerical Recipes Personal Single-User License lets you personally
! use Numerical Recipes code ("the code") on any number of computers,
! but only one computer at a time. You are not permitted to allow anyone
! else to access or use the code. You may, under this license, transfer
! precompiled, executable applications incorporating the code to other,
! unlicensed, persons, providing that (i) the application is
! noncommercial (i.e., does not involve the selling or licensing of the
! application for a fee), and (ii) the application was first developed,
! compiled, and successfully run by you, and (iii) the code is bound
! into the application in such a manner that it cannot be accessed as
! individual routines and cannot practicably be unbound and used in
! other programs. That is, under this license, your application user
! must not be able to use Numerical Recipes code as part of a program
! library or "mix and match" workbench.
!
! Businesses and organizations that purchase the disk or code download,
! and that thus acquire one or more Numerical Recipes Personal
! Single-User Licenses, may permanently assign those licenses, in the
! number acquired, to individual employees. Such an assignment must be
! made before the code is first used and, once made, it is irrevocable
! and can not be transferred.
!
! If you do not hold a Numerical Recipes License, this code is only for
! informational and educational purposes but cannot be used.
implicit none
public :: Euler ! Euler method with equal time-steps increments
public :: RK4 ! Fourth-order Runge-Kutta method with equal time-steps increments
public :: RK4as ! Fourth-order Runge-Kutta method with adaptive stepsize control
public :: RBstiff ! Rosenbrock Method for stiff ODEs with adaptive stepsize control
! ------------------------------------------------------------------
! NAME
! Euler
! PURPOSE
! Integration of Ordinary Differential Equations using Euler method.
! Starting from N initial values vstart known at x1, use Euler to advance nstep equal increments to x2.
! Results are stored in the variables xout and yout.
!> \brief Euler.
!> \details Starting from $N$ initial values $v_{start}$ known at $x_1$,
!> use Euler to advance $n$-steps equal increments to $x_2$.
!> Results are stored in the variables $x_{out}$ and $y_{out}$.
!> If you use "para" as parameters for derivs, the interface for derivs has to look like:
!>
!> interface
!> subroutine derivs( x, y, para, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(in) :: para ! parameter for the derivatives
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
!>
!> Elsewise "para" must be left out:
!>
!> interface
!> subroutine derivs( x, y, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
! INTENT(IN)
!> \param[in] "real(sp/dp), dimension(:) :: vstart" initial conditions.
!> $N$ inital values known at the time $x_1$
!> (given $N$ ODEs)
!> \param[in] "real(sp/dp) :: x1" initial time.
!> \param[in] "real(sp/dp) :: x2" final time.
!> \param[in] "real(sp/dp) :: h" size of the incremental time step (fixed).
!> \param[in] "interface :: derivs_sp/dp" returns derivatives $dydx$ of $y$ at $x$.
! INTENT(INOUT)
! None
! INTENT(OUT)
!> \param[out] "real(sp/dp), dimension(:), allocatable :: xout" storage for outputs (time).
!> \param[out] "real(sp/dp), dimension(:), allocatable :: yout" storage for outputs
!> (incremented variables).
!> dim 1 = function evaluation
!> at any time point.
!> dim 2 = number of equations.
! INTENT(IN), OPTIONAL
!> \param[in] "real(sp/dp), dimension(:) :: para" parameter for the derivatives $dydx$
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
!> \note The user has to supply the subroutine derivs(x,y,dydx), which returns derivatives $dydx$ at $x$.
! EXAMPLE
! call Euler( vstart, x1, x2, h, derivs, xout, yout )
! --> see example in test directory --> test_mo_ode_solver
! LITERATURE
! http://en.wikipedia.org/wiki/Euler_method
! HISTORY
!> \author Giovanni Dalmasso
!> \date Jul 2012
! Modified, Giovanni Dalmasso, Mar 2013
! Modified, Sebastian Mueller, Jan 2015
interface Euler
module procedure Euler_sp, Euler_dp, Euler_para_sp, Euler_para_dp
end interface Euler
! ------------------------------------------------------------------
! NAME
! RK4
! PURPOSE
! Integration of Ordinary Differential Equations using a fourth-order Runge-Kutta method
! with fixed time-steps increments.
! Starting from N initial values vstart known at x1, use Euler to advance nstep equal increments to x2.
! Results are stored in the variables xout and yout.
!> \brief fourth-order Runge-Kutta.
!> \details Starting from $N$ initial values $v_{start}$ known at $x_1$,
!> use a fourth-order Runge-Kutta with fixed time-steps increments
!> to advance $n$-steps equal increments to $x_2$.
!> Results are stored in the variables $x_{out}$ and $y_{out}$.
!> If you use "para" as parameters for derivs, the interface for derivs has to look like:
!>
!> interface
!> subroutine derivs( x, y, para, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(in) :: para ! parameter for the derivatives
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
!>
!> Elsewise "para" must be left out:
!>
!> interface
!> subroutine derivs( x, y, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
! INTENT(IN)
!> \param[in] "real(sp/dp), dimension(:) :: vstart" initial conditions.
!> $N$ inital values known at the time $x_1$
!> (given $N$ ODEs)
!> \param[in] "real(sp/dp) :: x1" initial time.
!> \param[in] "real(sp/dp) :: x2" final time.
!> \param[in] "real(sp/dp) :: h" size of the incremental time step (fixed).
!> \param[in] "interface :: derivs_sp/dp" returns derivatives $dydx$ of $y$ at $x$.
! INTENT(INOUT)
! None
! INTENT(OUT)
!> \param[out] "real(sp/dp), dimension(:), allocatable :: xout" storage for outputs (time).
!> \param[out] "real(sp/dp), dimension(:), allocatable :: yout" storage for outputs
!> (incremented variables).
!> dim 1 = function evaluation
!> at any time point.
!> dim 2 = number of equations.
! INTENT(IN), OPTIONAL
!> \param[in] "real(sp/dp), dimension(:) :: para" parameter for the derivatives $dydx$
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
!> \note The user has to supply the subroutine derivs(x,y,dydx), which returns derivatives $dydx$ at $x$.
! EXAMPLE
! call RK4( vstart, x1, x2, h, derivs, xout, yout )
! --> see example in test directory --> test_mo_ode_solver
! LITERATURE
! 1) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 77 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 1 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1992
! 2) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
!> \author Giovanni Dalmasso
!> \date Jul 2012
! Modified, Giovanni Dalmasso, Mar 2013
! Modified, Sebastian Mueller, Jan 2015
interface RK4
module procedure RK4_sp, RK4_dp, RK4_para_sp, RK4_para_dp
end interface RK4
! ------------------------------------------------------------------
! NAME
! RK4as
! PURPOSE
! Integration of Ordinary Differential Equations using a fourth-order Runge-Kutta method
! with adaptive stepsize control.
! Integrate the array of starting values ystart from x1 to x2 with accuracy eps, storing intermediate results
! in the module variables. h1 should be set as a guessed first stepsize,
! hmin as the minimum allowed stepsize (can be zero).
! On output ystart is replaced by values at the end of the integration interval.
!> \brief fourth-order Runge-Kutta with adaptive stepsize control.
!> \details Integrate the array of starting values $y_{start} from $x_1$ to $x_2$ with accuracy $\varepsilon$,
!> storing intermediate results in the module variables.
!> $h_1$ should be set as a guessed first stepsize, $h_{min} as the minimum allowed stepsize (can be zero).
!> On output $y_{start}$ is replaced by values at the end of the integration interval.
!> If you use "para" as parameters for derivs, the interface for derivs has to look like:
!>
!> interface
!> subroutine derivs( x, y, para, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(in) :: para ! parameter for the derivatives
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
!>
!> Elsewise "para" must be left out:
!>
!> interface
!> subroutine derivs( x, y, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> end interface
! INTENT(IN)
!> \param[in] "real(sp/dp), dimension(:) :: vstart" initial conditions.
!> $N$ inital values known at the time $x_1$
!> (given $N$ ODEs)
!> \param[in] "real(sp/dp) :: x1" initial time.
!> \param[in] "real(sp/dp) :: x2" final time.
!> \param[in] "real(sp/dp) :: h" guessed first stepsize.
!> \param[in] "interface :: derivs_sp/dp" derivatives $dydx$ of $y$ at $x$.
! INTENT(INOUT)
! None
! INTENT(OUT)
!> \param[out] "real(sp/dp), dimension(:), allocatable :: xout" storage for outputs (time).
!> \param[out] "real(sp/dp), dimension(:), allocatable :: yout" storage for outputs
!> (incremented variables).
!> dim 1 = function evaluations
!> at any time point.
!> dim 2 = number of equations.
! INTENT(IN), OPTIONAL
!> \param[in] "real(sp/dp), optional :: hmin" minimum allowed stepsize (can be 0.)
!> DEFAULT: 0.0
!> \param[in] "real(sp/dp), optional :: eps" accuracy (overall tolerance level)
!> DEFAULT: 1E-6
!> \param[in] "real(sp/dp), dimension(:) :: para" parameter for the derivatives $dydx$
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
!> \note The user has to supply the subroutine derivs(x,y,dydx), which returns derivatives $dydx$ at $x$.
! EXAMPLE
! call RK4as( ystart, x1, x2, h, derivs, xout, yout, hmin, eps )
! --> see example in test directory --> test_mo_ode_solver
! LITERATURE
! 1) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 77 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 1 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1992
! 2) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
!> \author Giovanni Dalmasso
!> \date Jul 2012
! Modified, Giovanni Dalmasso, Mar 2013
! Modified, Sebastian Mueller, Jan 2015
interface RK4as
module procedure RK4as_sp, RK4as_dp, RK4as_para_sp, RK4as_para_dp
end interface RK4as
! ------------------------------------------------------------------
! NAME
! RBstiff
! PURPOSE
! Integration of Ordinary Differential Equations using a Rosenbrock Method
! with adaptive stepsize control.
! Integrate the array of starting values ystart from x1 to x2 with accuracy eps, storing intermediate results
! in the module variables. h1 should be set as a guessed first stepsize,
! hmin as the minimum allowed stepsize (can be zero).
! On output ystart is replaced by values at the end of the integration interval.
!> \brief Rosenbrock Method with adaptive stepsize control.
!> \details Integrate the array of starting values $y_{start} from $x_1$ to $x_2$ with accuracy $\varepsilon$,
!> storing intermediate results in the module variables.
!> $h_1$ should be set as a guessed first stepsize, $h_{min} as the minimum allowed stepsize (can be zero).
!> On output $y_{start}$ is replaced by values at the end of the integration interval.
!> If you use "para" as parameters for derivs/jacobn, the interface for derivs has to look like:
!>
!> interface
!> subroutine derivs( x, y, para, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(in) :: para ! parameter for the derivatives
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> subroutine jacobn(x, y, para, dfdx, dfdy)
!> use mo_kind, only: sp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp),dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp),dimension(:), intent(in) :: para ! parameters for derivs
!> real(sp/dp),dimension(:), intent(out) :: dfdx ! derivatives of f for x
!> real(sp/dp),dimension(:,:), intent(out) :: dfdy ! jacobi-matrix df/dy
!> end subroutine jacobn
!> end interface
!>
!> Elsewise "para" must be left out:
!>
!> interface
!> subroutine derivs( x, y, dydx )
!> use mo_kind, only: sp/dp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp), dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp), dimension(:), intent(out) :: dydx ! derivatives of y
!> end subroutine derivs
!> subroutine jacobn(x, y, dfdx, dfdy)
!> use mo_kind, only: sp
!> implicit none
!> real(sp/dp), intent(in) :: x ! time
!> real(sp/dp),dimension(:), intent(in) :: y ! unknowns of the equations
!> real(sp/dp),dimension(:), intent(out) :: dfdx ! derivatives of f for x
!> real(sp/dp),dimension(:,:), intent(out) :: dfdy ! jacobi-matrix df/dy
!> end subroutine jacobn
!> end interface
! INTENT(IN)
!> \param[in] "real(sp/dp), dimension(:) :: vstart" initial conditions.
!> $N$ inital values known at the time $x_1$
!> (given $N$ ODEs)
!> \param[in] "real(sp/dp) :: x1" initial time.
!> \param[in] "real(sp/dp) :: x2" final time.
!> \param[in] "real(sp/dp) :: h" guessed first stepsize.
!> \param[in] "interface :: derivs_sp/dp" derivatives $dydx$ of $y$ at $x$.
!> \param[in] "interface :: jacobn_sp/dp" derivatives of the RHS: $dfdx$ and $dfdy$
! INTENT(INOUT)
! None
! INTENT(OUT)
!> \param[out] "real(sp/dp), dimension(:), allocatable :: xout" storage for outputs (time).
!> \param[out] "real(sp/dp), dimension(:), allocatable :: yout" storage for outputs
!> (incremented variables).
!> dim 1 = function evaluations
!> at any time point.
!> dim 2 = number of equations.
! INTENT(IN), OPTIONAL
!> \param[in] "real(sp/dp), optional :: hmin" minimum allowed stepsize (can be 0.)
!> DEFAULT: 0.0
!> \param[in] "real(sp/dp), optional :: eps" accuracy (overall tolerance level)
!> DEFAULT: 1E-6
!> \param[in] "real(sp/dp), dimension(:) :: para" parameter for the derivatives
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
!> \note The user has to supply the subroutine derivs(x,y,dydx), which returns derivatives $dydx$ at $x$.
! EXAMPLE
! call RK4as( ystart, x1, x2, h, derivs, jacobn, xout, yout, hmin, eps )
! --> see example in test directory --> test_mo_ode_solver
! LITERATURE
! 1) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 77 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 1 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1992
! 2) Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
!> \author Sebastian Mueller
!> \date Mar 2015
interface RBstiff
module procedure RBstiff_sp, RBstiff_dp, RBstiff_para_sp, RBstiff_para_dp
end interface RBstiff
private
! Private method
interface CashKarpRK
module procedure CashKarpRK_sp, CashKarpRK_dp, CashKarpRK_para_sp, CashKarpRK_para_dp
end interface CashKarpRK
interface RBstep
module procedure RBstep_sp, RBstep_dp, RBstep_para_sp, RBstep_para_dp
end interface RBstep
! ------------------------------------------------------------------
contains
! ------------------------------------------------------------------
! SINGLE PRECISION Euler
subroutine Euler_sp( ystart, x1, x2, h, derivs, xout, yout ) ! all obligatory
use mo_kind, only: i4, sp
implicit none
! Intent IN
real(sp), dimension(:), intent(in) :: ystart ! initial conditions
real(sp), intent(in) :: x1, x2 ! initial and final time
real(sp), intent(in) :: h ! step size
! Intent OUT
real(sp), dimension(:), allocatable, intent(out) :: xout
real(sp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: sp
implicit none
real(sp), intent(in) :: x ! time
real(sp), dimension(:), intent(in) :: y ! unknowns of the equations
real(sp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: j ! counter
integer(i4) :: nstep ! number of steps
real(sp) :: x
real(sp), dimension( size(ystart) ) :: dy, y
y(:) = ystart(:) ! load starting values
nstep = nint( (x2-x1)/h, i4 ) ! find number of steps
if ( allocated(xout) ) deallocate(xout) ! clear out old stored variables if necessary
if ( allocated(yout) ) deallocate(yout)
allocate( xout(nstep+1_i4) ) ! allocate storage for saved values
allocate( yout(nstep+1_i4, size(ystart)) )
yout(1,:) = y(:)
xout(1) = x1
x = x1
do j=1, nstep ! take nstep steps
call derivs( x, y, dy )
y = y + h*dy ! step EULER
if ( h .lt. tiny(1._sp) ) stop 'Euler_sp --> stepsize not significant!'
x = x+h
xout(j+1_i4) = x ! store intermediate steps
yout(j+1_i4,:) = y(:)
end do
end subroutine Euler_sp
! DOUBLE PRECISION Euler
subroutine Euler_dp( ystart, x1, x2, h, derivs, xout, yout ) ! all obligatory
use mo_kind, only: i4, dp
implicit none
! Intent IN
real(dp), dimension(:), intent(in) :: ystart ! initial conditions
real(dp), intent(in) :: x1, x2 ! initial and final time
real(dp), intent(in) :: h ! step size
! Intent OUT
real(dp), dimension(:), allocatable, intent(out) :: xout
real(dp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: dp
implicit none
real(dp), intent(in) :: x ! time
real(dp), dimension(:), intent(in) :: y ! unknowns of the equations
real(dp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: j ! counter
integer(i4) :: nstep ! nuber of steps
real(dp) :: x
real(dp), dimension( size(ystart) ) :: dy, y
y(:) = ystart(:) ! load starting values
nstep = nint( (x2-x1)/h, i4 ) ! find number of steps
if ( allocated(xout) ) deallocate(xout) ! clear out old stored variables if necessary
if ( allocated(yout) ) deallocate(yout)
allocate( xout(nstep+1_i4) ) ! allocate storage for saved values
allocate( yout(nstep+1_i4, size(ystart)) )
yout(1,:) = y(:)
xout(1) = x1
x = x1
do j=1, nstep ! take nstep steps
call derivs( x, y, dy )
y = y + h*dy ! step EULER
if ( h .lt. tiny(1._dp) ) stop 'Euler_dp --> stepsize not significant!'
x = x+h
xout(j+1_i4) = x ! store intermediate steps
yout(j+1_i4,:) = y(:)
end do
end subroutine Euler_dp
! ------------------------------------------------------------------
! SINGLE PRECISION 4th order RUNGE-KUTTA
subroutine RK4_sp( ystart, x1, x2, h, derivs, xout, yout ) ! all obligatory
use mo_kind, only: i4, sp
implicit none
! Intent IN
real(sp), dimension(:), intent(in) :: ystart ! initial conditions
real(sp), intent(in) :: x1, x2 ! initial and final time
real(sp), intent(in) :: h ! step size
! Intent OUT
real(sp), dimension(:), allocatable, intent(out) :: xout
real(sp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: sp
implicit none
real(sp), intent(in) :: x ! time
real(sp), dimension(:), intent(in) :: y ! unknowns of the equations
real(sp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: j ! counter
integer(i4) :: nstep ! nuber of steps
real(sp) :: x
real(sp) :: hh, h6, xh
real(sp), dimension( size(ystart) ) :: dy, dyt, dym, y, yt
y(:) = ystart(:) ! load starting values
nstep = nint( (x2-x1)/h, i4 ) ! find number of steps
if ( allocated(xout) ) deallocate(xout) ! clear out old stored variables if necessary
if ( allocated(yout) ) deallocate(yout)
allocate( xout(nstep+1_i4) ) ! allocate storage for saved values
allocate( yout(nstep+1_i4, size(ystart)) )
yout(1,:) = y(:)
xout(1) = x1
x = x1
hh = h*0.5_sp
h6 = h/6.0_sp
xh = x+hh
do j=1, nstep ! take nstep steps
call derivs( x, y, dy )
yt = y + hh*dy ! first step
call derivs( xh, yt, dyt ) ! second step
yt = y + hh*dyt
call derivs( xh, yt, dym ) ! third step
yt = y + h*dym
dym = dyt + dym
call derivs( x+h, yt, dyt ) ! fourth step
y = y + h6*( dy+dyt+2.0_sp*dym ) ! accumulate increments with proper weights
if ( h .lt. tiny(1._sp) ) stop 'RK4_sp --> stepsize not significant!'
x = x+h
xout(j+1_i4) = x ! store intermediate steps
yout(j+1_i4,:) = y(:)
end do
end subroutine RK4_sp
! DOUBLE PRECISION 4th order RUNGE-KUTTA
subroutine RK4_dp( ystart, x1, x2, h, derivs, xout, yout ) ! all obligatory
use mo_kind, only: i4, dp
implicit none
! Intent IN
real(dp), dimension(:), intent(in) :: ystart ! initial conditions
real(dp), intent(in) :: x1, x2 ! initial and final time
real(dp), intent(in) :: h ! step size
! Intent OUT
real(dp), dimension(:), allocatable, intent(out) :: xout
real(dp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: dp
implicit none
real(dp), intent(in) :: x ! time
real(dp), dimension(:), intent(in) :: y ! unknowns of the equations
real(dp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: j ! counter
integer(i4) :: nstep ! nuber of steps
real(dp) :: x
real(dp) :: hh, h6, xh
real(dp), dimension( size(ystart) ) :: dy, dyt, dym, y, yt
y(:) = ystart(:) ! load starting values
nstep = nint( (x2-x1)/h, i4 ) ! find number of steps
if ( allocated(xout) ) deallocate(xout) ! clear out old stored variables if necessary
if ( allocated(yout) ) deallocate(yout)
allocate( xout(nstep+1_i4) ) ! allocate storage for saved values
allocate( yout(nstep+1_i4, size(ystart)) )
yout(1,:) = y(:)
xout(1) = x1
x = x1
hh = h*0.5_dp
h6 = h/6.0_dp
xh = x+hh
do j=1, nstep ! take nstep steps
call derivs( x, y, dy )
yt = y + hh*dy ! first step
call derivs( xh, yt, dyt ) ! second step
yt = y + hh*dyt
call derivs( xh, yt, dym ) ! third step
yt = y + h*dym
dym = dyt + dym
call derivs( x+h, yt, dyt ) ! fourth step
y = y + h6*( dy+dyt+2.0_dp*dym ) ! accumulate increments with proper weights
if ( h .lt. tiny(1._dp) ) stop 'RK4_sp --> stepsize not significant!'
x = x+h
xout(j+1_i4) = x ! store intermediate steps
yout(j+1_i4,:) = y(:)
end do
end subroutine RK4_dp
! ------------------------------------------------------------------
! SINGLE PRECISION 4th order RUNGE-KUTTA Adaptive Step-size
subroutine RK4as_sp( ystart, x1, x2, h, derivs, xout, yout, & ! obligatory
hmin, eps ) ! optional
use mo_kind, only: i4, sp
use mo_nrutil, only: reallocate
implicit none
! Intent IN
real(sp), intent(in) :: x1, x2 ! initial and final time
real(sp), intent(in) :: h ! guessed step size
real(sp), dimension(:), intent(in) :: ystart ! initial conditions
real(sp), optional, intent(in) :: hmin
real(sp), optional, intent(in) :: eps
! Intent OUT
real(sp), dimension(:), allocatable, intent(out) :: xout
real(sp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: sp
implicit none
real(sp), intent(in) :: x ! time
real(sp), dimension(:), intent(in) :: y ! unknowns of the equations
real(sp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: nstep ! nuber of steps
integer(i4) :: nok, nbad ! number of good and bad (but retried and fixed) steps taken
integer(i4) :: kount ! total number of saved steps
real(sp) :: hIN, hdid, hnext, x, hminIN, epsIN, xsav, dxsav, errmax, htemp, xnew
real(sp), dimension(size(ystart)) :: dydx, y, yscal, yerr, ytemp
! Pointers
real(sp), dimension(:), pointer :: xp
real(sp), dimension(:,:), pointer :: yp
! parameters
integer(i4), parameter :: MAXstp = 1000000_i4 ! max number of steps
real(sp), parameter :: safety = 0.9_sp, pgrow = -0.2_sp, pshrnk = -0.25_sp
real(sp) :: errcon
errcon = exp( (1._sp/pgrow)*log(5._sp/safety) )
if( present(hmin) ) then
hminIN = hmin
else
hminIN = 0.0_sp
end if
if( present(eps) ) then
epsIN = eps
else
epsIN = 1e-6_sp
end if
x = x1
hIN = sign( h, x2-x1 )
nok = 0_i4
nbad = 0_i4
kount = 0_i4
y(:) = ystart(:)
dxsav = tiny(1._sp)
xsav = x-2.0_sp*dxsav ! assures storage of first step
nullify( xp, yp )
allocate( xp(256) )
allocate( yp(size(xp), size(ystart)) )
call save_a_step ! save initial step
do nstep=1, MAXstp ! take at most MAXstp steps
call derivs( x, y, dydx )
yscal(:) = abs( y(:) ) + abs( hIN*dydx(:) ) + tiny(1._sp) ! scaling used to monitor accuracy
! --> CAN BE MODIFIED...
if ( abs(x-xsav) .gt. abs(dxsav) ) call save_a_step ! store intermediate results
if ( (x+hIN-x2)*(x+hIN-x1) .gt. 0.0_sp ) hIN = x2-x ! if stepsize can overshoot, decrease
do
call CashKarpRK( y, dydx, x, hIN, ytemp, yerr, derivs ) ! take a step
errmax = maxval( abs(yerr(:)/yscal(:)) )/epsIN ! evaluate accuracy
if ( errmax .lt. 1.0_sp ) exit ! step succeeded
htemp = safety*hIN*(errmax**pshrnk) ! truncation error too large, reduce stepsize
hIN = sign( max( abs(htemp), 0.1_sp*abs(hIN) ), hIN ) ! no more than a factor of 10
xnew = x+hIN
if ( abs(xnew-x) .lt. epsilon(1.0_sp) ) stop 'RK4as_sp --> hey !!! stepsize underflow!'
end do
if ( errmax .gt. errcon ) then ! compute size of next step
hnext = safety*hIN*(errmax**pgrow)
else ! no more than a factor of 5 increase
hnext = 5.0_sp*hIN
end if
hdid = hIN
x = x+hIN
y(:) = ytemp(:)
if ( abs(hdid-hIN) .gt. epsilon(1.0_sp) ) then
nbad = nbad+1_i4
else
nok = nok+1_i4
end if
if ( (x-x2)*(x2-x1) .ge. 0.0_sp ) then ! are we done?!?!
call save_a_step ! save final step
allocate( xout(kount) ) ! allocate storage for outputs
xout(:) = xp(1:kount)
allocate( yout(kount, size(yp,2)) ) ! allocate storage for outputs
yout(:,:) = yp(1:kount, :)
deallocate( xp, yp ) ! clear out old stored variables
return ! normal exit
end if
if ( abs(hnext-hminIN) .lt. epsilon(1.0_sp) ) stop 'RK4as_sp --> WTF! ...stepsize smaller than minimum!!!'
hIN = hnext
end do
stop 'RK4as_sp --> dude, too many steps!!!'
contains
subroutine save_a_step ! --> like a macro in C
kount = kount+1
if ( kount .gt. size(xp) ) then
xp=>reallocate( xp, 2*size(xp) )
yp=>reallocate( yp, size(xp), size(yp,2) )
end if
xp(kount) = x
yp(kount, :) = y(:)
xsav = x
end subroutine save_a_step
end subroutine RK4as_sp
! DOUBLE PRECISION 4th order RUNGE-KUTTA Adaptive Step-size
subroutine RK4as_dp( ystart, x1, x2, h, derivs, xout, yout, & ! obligatory
hmin, eps ) ! optional
use mo_kind, only: i4, dp
use mo_nrutil, only: reallocate
implicit none
! Intent IN
real(dp), intent(in) :: x1, x2 ! initial and final time
real(dp), intent(in) :: h ! guessed step size
real(dp), dimension(:), intent(in) :: ystart ! initial conditions
real(dp), optional, intent(in) :: hmin
real(dp), optional, intent(in) :: eps
! Intent OUT
real(dp), dimension(:), allocatable, intent(out) :: xout
real(dp), dimension(:,:), allocatable, intent(out) :: yout
interface
subroutine derivs( x, y, dydx )
use mo_kind, only: dp
implicit none
real(dp), intent(in) :: x ! time
real(dp), dimension(:), intent(in) :: y ! unknowns of the equations
real(dp), dimension(:), intent(out) :: dydx ! derivatives of y
end subroutine derivs
end interface
! Internal variables
integer(i4) :: nstep ! nuber of steps
integer(i4) :: nok, nbad ! number of good and bad (but retried and fixed) steps taken
integer(i4) :: kount ! total number of saved steps
real(dp) :: hIN, hdid, hnext, x, hminIN, epsIN, xsav, dxsav, errmax, htemp, xnew
real(dp), dimension(size(ystart)) :: dydx, y, yscal, yerr, ytemp
! Pointers
real(dp), dimension(:), pointer :: xp
real(dp), dimension(:,:), pointer :: yp
! parameters
integer(i4), parameter :: MAXstp = 1000000_i4 ! max nuber of steps
real(dp), parameter :: safety = 0.9_dp, pgrow = -0.2_dp, pshrnk = -0.25_dp
real(dp) :: errcon
errcon = exp( (1._dp/pgrow)*log(5._dp/safety) )
if( present(hmin) ) then
hminIN = hmin
else
hminIN = 0.0_dp
end if
if( present(eps) ) then
epsIN = eps
else
epsIN = 1e-6_dp
end if
x = x1
hIN = sign( h, x2-x1 )
nok = 0_i4
nbad = 0_i4
kount = 0_i4
y(:) = ystart(:)
dxsav = tiny(1._dp)
xsav = x-2.0_dp*dxsav ! assures storage of first step
nullify( xp, yp )
allocate( xp(256) )
allocate( yp(size(xp), size(ystart)) )
call save_a_step ! save initial step
do nstep=1, MAXstp ! take at most MAXstp steps
call derivs( x, y, dydx )
yscal(:) = abs( y(:) ) + abs( hIN*dydx(:) ) + tiny(1._dp) ! scaling used to monitor accuracy --> CAN BE MODIFIED...
if ( abs(x-xsav) .gt. abs(dxsav) ) call save_a_step ! store intermediate results
if ( (x+hIN-x2)*(x+hIN-x1) .gt. 0.0_dp ) hIN = x2-x ! if stepsize can overshoot, decrease
do
call CashKarpRK( y, dydx, x, hIN, ytemp, yerr, derivs ) ! take a step
errmax = maxval( abs(yerr(:)/yscal(:)) )/epsIN ! evaluate accuracy
if ( errmax .lt. 1.0_dp ) exit ! step succeeded
htemp = safety*hIN*(errmax**pshrnk) ! truncation error too large, reduce stepsize
hIN = sign( max( abs(htemp), 0.1_dp*abs(hIN) ), hIN ) ! no more than a factor of 10
xnew = x+hIN
if ( abs(xnew-x) .lt. epsilon(1._dp) ) stop 'RK4as_dp --> hey!!! stepsize underflow!'
end do
if ( errmax .gt. errcon ) then ! compute size of next step
hnext = safety*hIN*(errmax**pgrow)
else ! no more than a factor of 5 increase
hnext = 5.0_dp*hIN
end if
hdid = hIN