Skip to content

Commit 58682a8

Browse files
author
Daniel Ruprecht
committed
add split-explicit method by Knoth et al.
1 parent 079ff2c commit 58682a8

File tree

3 files changed

+173
-6
lines changed

3 files changed

+173
-6
lines changed

examples/boussinesq_2d_imex/plotgmrescounter.py

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -14,7 +14,9 @@
1414
udirk = np.load('dirk.npy')
1515
uimex = np.load('rkimex.npy')
1616
uref = np.load('uref.npy')
17+
usplit = np.load('split.npy')
1718

19+
print("Estimated discretisation error split explicit: %5.3e" % ( np.linalg.norm(usplit.flatten() - uref.flatten(), np.inf)/np.linalg.norm(uref.flatten(),np.inf) ))
1820
print("Estimated discretisation error of DIRK: %5.3e" % ( np.linalg.norm(udirk.flatten() - uref.flatten(), np.inf)/np.linalg.norm(uref.flatten(),np.inf) ))
1921
print("Estimated discretisation error of RK-IMEX: %5.3e" % ( np.linalg.norm(uimex.flatten() - uref.flatten(), np.inf)/np.linalg.norm(uref.flatten(),np.inf) ))
2022
print("Estimated discretisation error of SDC: %5.3e" % ( np.linalg.norm(uend.flatten() - uref.flatten(), np.inf)/np.linalg.norm(uref.flatten(),np.inf) ))

examples/boussinesq_2d_imex/rungmrescounter.py

Lines changed: 28 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -21,7 +21,7 @@
2121

2222
from unflatten import unflatten
2323

24-
from standard_integrators import dirk, bdf2, trapezoidal, rk_imex
24+
from standard_integrators import dirk, bdf2, trapezoidal, rk_imex, SplitExplicit
2525

2626
if __name__ == "__main__":
2727

@@ -46,8 +46,8 @@
4646

4747
# setup parameters "in time"
4848
t0 = 0
49-
Tend = 3000
50-
Nsteps = 500
49+
Tend = 3000
50+
Nsteps = 500
5151
dt = Tend/float(Nsteps)
5252

5353
# This comes as read-in for the problem class
@@ -96,19 +96,33 @@
9696
print("CFL number of acoustics (horizontal): %4.2f" % cfl_acoustic_hor)
9797
print("CFL number of acoustics (vertical): %4.2f" % cfl_acoustic_ver)
9898

99-
dirkp = dirk(P, np.min([4,dirk_order]))
99+
method_split = 'MIS4_4'
100+
# method_split = 'RK3'
101+
splitp = SplitExplicit(P, method_split, pparams)
100102
u0 = uinit.values.flatten()
103+
usplit = np.copy(u0)
104+
# for i in range(0,Nsteps):
105+
print(np.linalg.norm(usplit))
106+
for i in range(0,2*Nsteps):
107+
usplit = splitp.timestep(usplit, dt/2)
108+
print(np.linalg.norm(usplit))
109+
110+
dirkp = dirk(P, np.min([4,dirk_order]))
101111
udirk = np.copy(u0)
102112
print("Running DIRK ....")
113+
print(np.linalg.norm(udirk))
103114
for i in range(0,Nsteps):
104115
udirk = dirkp.timestep(udirk, dt)
116+
print np.linalg.norm(udirk)
105117

106118
rkimex = rk_imex(P, dirk_order)
107119
uimex = np.copy(u0)
108120
dt_imex = dt
109121
print("Running RK-IMEX ....")
110122
for i in range(0,Nsteps):
123+
# print("Running RK-IMEWX Step: %4.2f" % dt_imex)
111124
uimex = rkimex.timestep(uimex, dt_imex)
125+
print(np.linalg.norm(uimex))
112126

113127
# call main function to get things done...
114128
print("Running SDC...")
@@ -123,6 +137,7 @@
123137
for i in range(0,10*Nsteps):
124138
uref = rkimexref.timestep(uref, dt_ref)
125139

140+
usplit = unflatten(usplit, 4, P.N[0], P.N[1])
126141
udirk = unflatten(udirk, 4, P.N[0], P.N[1])
127142
uimex = unflatten(uimex, 4, P.N[0], P.N[1])
128143
uref = unflatten(uref, 4, P.N[0], P.N[1])
@@ -131,8 +146,17 @@
131146
np.save('sdc', uend.values)
132147
np.save('dirk', udirk)
133148
np.save('rkimex', uimex)
149+
np.save('split', usplit)
134150
np.save('uref', uref)
151+
152+
print "diff split ",np.linalg.norm(uref-usplit)
153+
print "diff dirk ",np.linalg.norm(uref-udirk)
154+
print "diff rkimex ",np.linalg.norm(uref-uimex)
155+
print "diff sdc ",np.linalg.norm(uref-uend.values)
135156

157+
print(" #### Logging report for Split #### " )
158+
print("Total number of matrix multiplcations: %5i" % splitp.logger.nsmall)
159+
136160
print(" #### Logging report for DIRK-%1i #### " % dirkp.order)
137161
print("Number of calls to implicit solver: %5i" % dirkp.logger.solver_calls)
138162
print("Total number of GMRES iterations: %5i" % dirkp.logger.iterations)

examples/boussinesq_2d_imex/standard_integrators.py

Lines changed: 143 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -270,8 +270,149 @@ def f_solve(self, b, alpha, u0):
270270
return sol
271271

272272
#
273-
# A diagonally implicit Runge-Kutta method of order 2, 3 or 4
273+
# Split-Explicit method
274274
#
275+
276+
class SplitExplicit:
277+
278+
def __init__(self, problem, method, pparams):
279+
280+
assert isinstance(problem, boussinesq_2d_imex), "problem is wrong type of object"
281+
self.Ndof = np.shape(problem.M)[0]
282+
self.method = method
283+
self.logger = logging()
284+
self.problem = problem
285+
self.pparams = pparams
286+
self.NdofTher = 2*problem.N[0]*problem.N[1]
287+
self.NdofMom = 2*problem.N[0]*problem.N[1]
288+
289+
print "dx ",problem.h[0]
290+
print "dz ",problem.h[1]
291+
292+
assert self.method in ["MIS4_4","RK3"], 'Method must be MIS4_4'
293+
294+
if (self.method=='RK3'):
295+
self.nstages=3
296+
self.aRunge=np.zeros((4,4))
297+
self.aRunge[0,0]= 1./3.
298+
self.aRunge[1,1]= 1./2.
299+
self.aRunge[2,2]= 1.
300+
self.dRunge=np.zeros((4,4))
301+
self.gRunge=np.zeros((4,4))
302+
if (self.method=='MIS4_4'):
303+
self.nstages=4
304+
self.aRunge=np.zeros((4,4))
305+
self.aRunge[0,0]= 0.38758444641450318
306+
self.aRunge[1,0]= -2.5318448354142823e-002
307+
self.aRunge[1,1]= 0.38668943087310403
308+
self.aRunge[2,0]= 0.20899983523553325
309+
self.aRunge[2,1]= -0.45856648476371231
310+
self.aRunge[2,2]= 0.43423187573425748
311+
self.aRunge[3,0]= -0.10048822195663100
312+
self.aRunge[3,1]= -0.46186171956333327
313+
self.aRunge[3,2]= 0.83045062122462809
314+
self.aRunge[3,3]= 0.27014914900250392
315+
self.dRunge=np.zeros((4,4))
316+
self.dRunge[1,1]= 0.52349249922385610
317+
self.dRunge[2,1]= 1.1683374366893629
318+
self.dRunge[2,2]= -0.75762080241712637
319+
self.dRunge[3,1]= -3.6477233846797109e-002
320+
self.dRunge[3,2]= 0.56936148730740477
321+
self.dRunge[3,3]= 0.47746263002599681
322+
self.gRunge=np.zeros((4,4))
323+
self.gRunge[1,1]= 0.13145089796226542
324+
self.gRunge[2,1]= -0.36855857648747881
325+
self.gRunge[2,2]= 0.33159232636600550
326+
self.gRunge[3,1]= -6.5767130537473045E-002
327+
self.gRunge[3,2]= 4.0591093109036858E-002
328+
self.gRunge[3,3]= 6.4902111640806712E-002
329+
self.dtRunge=np.zeros(self.nstages);
330+
for i in range(0,self.nstages):
331+
self.dtRunge[i]=0
332+
temp=1.
333+
for j in range(0,i+1):
334+
self.dtRunge[i]=self.dtRunge[i]+self.aRunge[i,j]
335+
temp=temp-self.dRunge[i,j]
336+
self.dRunge[i,0]=temp
337+
for j in range(0,i+1):
338+
self.aRunge[i,j]=self.aRunge[i,j]/self.dtRunge[i]
339+
self.gRunge[i,j]=self.gRunge[i,j]/self.dtRunge[i]
340+
341+
self.U = np.zeros((self.Ndof,self.nstages+1))
342+
self.F = np.zeros((self.Ndof,self.nstages))
343+
self.FSlow = np.zeros((self.Ndof))
344+
self.nsMin = 8
345+
self.logger.nsmall = 0
346+
347+
def NumSmallTimeSteps(self , dx, dz, dt):
348+
349+
cs = self.pparams['c_s']
350+
ns = dt / (.9 / np.sqrt(1/(dx*dx)+1/(dz*dz)) / cs)
351+
ns = max(np.int(np.ceil(ns)),self.nsMin)
352+
# print " dx " , dx
353+
# print " dz " , dz
354+
# print " cs " , cs
355+
# print 1/np.sqrt(1/(dx*dx)+1/(dz*dz))/cs
356+
# print " dt ",dt," ns ",ns
357+
358+
return ns
359+
360+
361+
def timestep(self, u0, dt):
362+
363+
self.U[:,0] = u0
364+
365+
self.ns = self.NumSmallTimeSteps(self.problem.h[0], self.problem.h[1], dt)
366+
367+
for i in range(0,self.nstages):
368+
self.F[:,i] = self.f_slow(self.U[:,i])
369+
self.FSlow[:] = 0.
370+
for j in range(0,i+1):
371+
self.FSlow += (self.aRunge[i,j]*self.F[:,j] + self.gRunge[i,j]/dt*(self.U[:,j] - u0))
372+
self.U[:,i+1] = 0
373+
for j in range(0,i+1):
374+
self.U[:,i+1] += self.dRunge[i,j]*self.U[:,j]
375+
nsLoc = np.int(np.ceil(self.ns*self.dtRunge[i]))
376+
self.logger.nsmall += nsLoc
377+
dtLoc = dt*self.dtRunge[i]
378+
dTau = dtLoc/nsLoc
379+
self.U[:,i+1] = self.VerletLin(self.U[:,i+1], self.FSlow, nsLoc, dTau)
380+
u0 = self.U[:,self.nstages]
381+
return u0
382+
383+
384+
def VerletLin(self, u0, FSlow, ns, dTau):
385+
for i in range(0,ns):
386+
u0[0:self.NdofMom] += dTau * (self.f_fastMom(u0) + FSlow[0:self.NdofMom])
387+
u0[self.NdofMom:self.Ndof] += dTau * (self.f_fastTher(u0) + FSlow[self.NdofMom:self.Ndof])
388+
389+
return u0
390+
391+
def RK3Lin(self, u0, FSlow, ns, dTau):
392+
393+
u = u0
394+
for i in range(0,ns):
395+
u = u0 + dTau/3.*(self.f_fast(u) + FSlow)
396+
u = u0 + dTau/2.*(self.f_fast(u) + FSlow)
397+
u = u0 + dTau*(self.f_fast(u) + FSlow)
398+
u0 = u
399+
400+
return u0
401+
402+
def f_slow(self, u):
403+
return self.problem.D_upwind.dot(u)
404+
405+
def f_fast(self, u):
406+
return self.problem.M.dot(u)
407+
408+
def f_fastMom(self, u):
409+
return self.problem.M[0:self.NdofMom,self.NdofMom:self.Ndof]\
410+
.dot(u[self.NdofMom:self.Ndof])
411+
412+
def f_fastTher(self, u):
413+
return self.problem.M[self.NdofMom:self.Ndof,0:self.NdofMom]\
414+
.dot(u[0:self.NdofMom])
415+
275416
class dirk:
276417

277418
def __init__(self, problem, order):
@@ -283,7 +424,7 @@ def __init__(self, problem, order):
283424
self.problem = problem
284425

285426
assert self.order in [2,22,3,4], 'Order must be 2,22,3,4'
286-
427+
287428
if (self.order==2):
288429
self.nstages = 1
289430
self.A = np.zeros((1,1))

0 commit comments

Comments
 (0)