-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtimestepping_methods.py
162 lines (140 loc) · 4.14 KB
/
timestepping_methods.py
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
#timestepping options for averaged_sw_explicit.py
#rk2/rk4/heuns/ssprk3/leapfrog
def rk2(U, USlow_in, USlow_out, DU, V, W,
expt, ensemble, cheby, cheby2, SlowSolver, wt, dt):
#Average the nonlinearity
cheby.apply(U, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V)
#Step forward V
V.assign(U + 0.5*V)
#transform forwards to U^{n+1/2}
cheby2.apply(V, DU, dt/2)
#Average the nonlinearity
cheby.apply(DU, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V)
#Advance U
cheby2.apply(U, DU, dt/2)
V.assign(DU + V)
#transform forwards to next timestep
cheby2.apply(V, U, dt/2)
def rk4(U, USlow_in, USlow_out, DU, U1, U2, U3, V1, V2, V3, V, W,
expt, ensemble, cheby, cheby2, SlowSolver, wt, dt):
#Average the nonlinearity
cheby.apply(U, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V1)
#Step forward U1
U1.assign(U + 0.5*V1)
#Average the nonlinearity
cheby2.apply(U1, DU, dt/2)
cheby.apply(DU, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V2)
#Step forward U2
cheby2.apply(U, V, dt/2)
U2.assign(V + 0.5*V2)
#Average the nonlinearity
cheby.apply(U2, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V3)
#Step forward U1
U3.assign(V + V3)
#Average the nonlinearity
U1.assign(V2 + V3)
cheby2.apply(U1, U2, dt/2)
cheby2.apply(U3, DU, dt/2)
cheby.apply(DU, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U3)
cheby2.apply(V1, U1, dt)
cheby2.apply(U, V, dt)
U.assign(V + 1/6*U1 + 1/3*U2 + 1/6*U3)
def heuns(U, USlow_in, USlow_out, DU, U1, U2, W,
expt, ensemble, cheby, cheby2, SlowSolver, wt, dt):
#Average the nonlinearity
cheby.apply(U, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U1)
#Step forward U1
U1.assign(U + U1)
#Average the nonlinearity
cheby2.apply(U1, DU, dt)
cheby.apply(DU, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U2)
#Step forward U2
cheby2.apply(U, DU, dt)
U2.assign(DU + U2)
#transform forwards to next timestep
cheby2.apply(U1, U, dt)
U.assign(0.5*U + 0.5*U2)
def ssprk3(U, USlow_in, USlow_out, DU, U1, U2, W,
expt, ensemble, cheby, cheby2, SlowSolver, wt, dt):
#Average the nonlinearity
cheby.apply(U, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U1)
#Step forward U1
DU.assign(U + U1)
cheby2.apply(DU, U1, dt)
#Average the nonlinearity
cheby.apply(U1, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U2)
#Step forward U2
DU.assign(U1 + U2)
cheby2.apply(DU, U2, -dt/2)
cheby2.apply(U, U1, dt/2)
U2.assign(0.75*U1 + 0.25*U2)
#Average the nonlinearity
cheby.apply(U2, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, U1)
#Advance U
DU.assign(U2 + U1)
cheby2.apply(DU, U2, dt/2)
cheby2.apply(U, U1, dt)
U.assign(1/3*U1 + 2/3*U2)
def leapfrog(U, USlow_in, USlow_out, U_old, U_new, DU, U1, U2, V, W,
expt, ensemble, cheby, cheby2, SlowSolver, wt, dt, asselin):
#Average the nonlinearity
cheby.apply(U, USlow_in, expt)
SlowSolver.solve()
cheby.apply(USlow_out, DU, -expt)
DU *= wt
ensemble.allreduce(DU, V)
#Step forward V
cheby2.apply(U_old, DU, dt)
V.assign(DU + 2*V)
cheby2.apply(V, U_new, dt)
#Asselin filter
cheby2.apply(U_old, U1, dt)
cheby2.apply(U_new, U2, -dt)
V.assign((U1+U2)*0.5 - U)
U_old.assign(U + asselin*V)
#Advance U
U.assign(U_new)