SPEED
MAKE_BUTCHERARRAY.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_butcherarray (order, stages, a_low, b_low, c)
 Makes Butcher array for Runge-Kutta scheme.
 
subroutine make_low_storage_coefficients (stages, a, b, c, a_low, b_low)
 Computes Low-storage memory coefficients for Runge-Kutta scheme.
 

Function/Subroutine Documentation

◆ make_butcherarray()

subroutine make_butcherarray ( integer*4  order,
integer*4  stages,
real*8, dimension(stages)  a_low,
real*8, dimension(stages)  b_low,
real*8, dimension(stages)  c 
)

Makes Butcher array for Runge-Kutta scheme.

Author
Ilario Mazzieri
Note
Compute low storage Butcher array for Runge-Kutta(p,s) time integration scheme
Date
September, 2013
Version
1.0
Parameters
[in]orderp
[in]stagess
[out]A_lowlow storage Butcher array A
[out]B_lowlow storage Butcher array B
[out]cButcher array c

Definition at line 35 of file MAKE_BUTCHERARRAY.f90.

36
37 implicit none
38
39 integer*4:: order, stages
40 real*8, dimension(stages,stages) :: a
41 real*8, dimension(stages) :: b, c, a_low, b_low
42
43 a = 0; b = 0; c = 0; a_low = 0; b_low = 0;
44
45 select case(order)
46
47 case(2)
48 ! RK(2,2)
49 select case(stages)
50 case(2)
51 c(2) = 1
52 b(1) = 0.5; b(2)=0.5
53 a(2,1) = 1
54
55 call make_low_storage_coefficients(stages,a,b,c,a_low,b_low)
56
57 case default
58 write(*,'(A)') 'RK scheme not implemented'
59 end select
60
61 case(3)
62 ! RK(3,3)
63 select case(stages)
64 case(3)
65 c(2) = 0.5; c(3) = 1
66 b(1) = 1.d0/6; b(2) = 2.d0/3; b(3) = 1.d0/6
67 a(2,1) = 0.5;
68 a(3,1) = -1; a(3,2) = 2
69
70 call make_low_storage_coefficients(stages,a,b,c,a_low,b_low)
71
72
73 case default
74 write(*,'(A)') 'RK scheme not implemented'
75 end select
76
77 case(4)
78
79 select case(stages)
80 ! RK(4,4)
81 case(4)
82 c(2) = 0.5; c(3) = 0.5; c(4) = 1
83 b(1) = 1.d0/6; b(2) = 1.d0/3; b(3) = 1.d0/3; b(4) = 1.d0/6
84 a(2,1) = 0.5;
85 a(3,2) = 0.5;
86 a(4,3) = 1;
87
88 call make_low_storage_coefficients(stages,a,b,c,a_low,b_low)
89
90
91 case default
92 write(*,'(A)') 'RK scheme not implemented'
93
94 end select
95
96
97
98
99
100 case default
101 write(*,'(A)') 'RK scheme not implemented'
102
103 end select
104
105
subroutine make_low_storage_coefficients(stages, a, b, c, a_low, b_low)
Computes Low-storage memory coefficients for Runge-Kutta scheme.

References make_low_storage_coefficients().

Here is the call graph for this function:

◆ make_low_storage_coefficients()

subroutine make_low_storage_coefficients ( integer*4  stages,
real*8, dimension(stages,stages)  a,
real*8, dimension(stages)  b,
real*8, dimension(stages)  c,
real*8, dimension(stages)  a_low,
real*8, dimension(stages)  b_low 
)

Computes Low-storage memory coefficients for Runge-Kutta scheme.

Author
Ilario Mazzieri
Note
Compute low storage Butcher array for Runge-Kutta(p,s) time integration scheme where p is the order of the scheme and s is the number of stages
Date
September, 2013
Version
1.0
Parameters
[in]stagess
[in]AButcher array
[in]bButcher array
[in]cButcher array
[out]A_lowlow storage Butcher array A
[out]B_lowlow storage Butcher array B

Definition at line 126 of file MAKE_BUTCHERARRAY.f90.

127
128 implicit none
129
130 integer*4:: stages
131 real*8, dimension(stages,stages) :: a
132 real*8, dimension(stages) :: b, c, a_low, b_low
133 integer*4 :: i
134
135
136 ! Compute Low storage coefficients according to the
137 ! paper Toulorge and Desmet: "Optimal Runge-Kutta schemes for discontinuous
138 ! Galerkin space discretizations applied to wave propagation problems"
139
140 a_low(1) = 0;
141
142 do i = 1, stages -1
143 b_low(i) = a(i+1,i)
144 enddo
145 b_low(stages) = b(stages)
146
147 do i = 2, stages
148 if (b(i) .ne. 0.d0) then
149 a_low(i) = (b(i-1) - b_low(i-1))/b(i);
150 else
151 a_low(i) = (a(i+1,i-1) - c(i))/b_low(i);
152 endif
153 enddo
154
155

Referenced by make_butcherarray().

Here is the caller graph for this function: