Skip to content

Commit 29065dc

Browse files
BenjaminRodenbergMakisHIshaanDesai
authored
Fix checkpointing and subcycling for solid.py (#475)
* Fix checkpointing. Requires precice/fenics-adapter#170. --------- Co-authored-by: Gerasimos Chourdakis <[email protected]> Co-authored-by: Ishaan Desai <[email protected]>
1 parent e0b33fa commit 29065dc

File tree

1 file changed

+14
-8
lines changed

1 file changed

+14
-8
lines changed

perpendicular-flap/solid-fenics/solid.py

Lines changed: 14 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -59,7 +59,6 @@ def neumann_boundary(x, on_boundary):
5959
v = TestFunction(V)
6060

6161
u_np1 = Function(V)
62-
saved_u_old = Function(V)
6362

6463
# function known from previous timestep
6564
u_n = Function(V)
@@ -80,15 +79,18 @@ def neumann_boundary(x, on_boundary):
8079

8180
precice_dt = precice.get_max_time_step_size()
8281
fenics_dt = precice_dt # if fenics_dt == precice_dt, no subcycling is applied
83-
# fenics_dt = 0.02 # if fenics_dt < precice_dt, subcycling is applied
82+
# n_substeps = 5 # number of substeps per window
83+
# fenics_dt = precice_dt / n_substeps # if fenics_dt < precice_dt, subcycling is applied
8484
dt = Constant(np.min([precice_dt, fenics_dt]))
8585

8686
# clamp the beam at the bottom
8787
bc = DirichletBC(V, Constant((0, 0)), fixed_boundary)
8888

8989
# alpha method parameters
90-
alpha_m = Constant(0)
91-
alpha_f = Constant(0)
90+
alpha_m = Constant(0.2)
91+
alpha_f = Constant(0.4)
92+
# alpha_m = Constant(0)
93+
# alpha_f = Constant(0)
9294

9395
"""
9496
Check requirements for alpha_m and alpha_f from
@@ -196,13 +198,14 @@ def avg(x_old, x_new, alpha):
196198
while precice.is_coupling_ongoing():
197199

198200
if precice.requires_writing_checkpoint(): # write checkpoint
199-
precice.store_checkpoint(u_n, t, n)
201+
precice.store_checkpoint((u_n, v_n, a_n), t, n)
200202

201203
precice_dt = precice.get_max_time_step_size()
202204
dt = Constant(np.min([precice_dt, fenics_dt]))
203205

204206
# read data from preCICE and get a new coupling expression
205-
read_data = precice.read_data(dt)
207+
# sample force F at $F(t_{n+1-\alpha_f})$ (see generalized alpha paper)
208+
read_data = precice.read_data((1 - float(alpha_f)) * dt)
206209

207210
# Update the point sources on the coupling boundary with the new read data
208211
Forces_x, Forces_y = precice.get_point_sources(read_data)
@@ -227,17 +230,20 @@ def avg(x_old, x_new, alpha):
227230

228231
# Either revert to old step if timestep has not converged or move to next timestep
229232
if precice.requires_reading_checkpoint(): # roll back to checkpoint
230-
u_cp, t_cp, n_cp = precice.retrieve_checkpoint()
233+
uva_cp, t_cp, n_cp = precice.retrieve_checkpoint()
234+
u_cp, v_cp, a_cp = uva_cp
231235
u_n.assign(u_cp)
236+
v_n.assign(v_cp)
237+
a_n.assign(a_cp)
232238
t = t_cp
233239
n = n_cp
234240
else:
241+
update_fields(u_np1, u_n, v_n, a_n)
235242
u_n.assign(u_np1)
236243
t += float(dt)
237244
n += 1
238245

239246
if precice.is_time_window_complete():
240-
update_fields(u_np1, saved_u_old, v_n, a_n)
241247
if n % 10 == 0:
242248
displacement_out << (u_n, t)
243249

0 commit comments

Comments
 (0)