Skip to content

Commit 4b7b925

Browse files
move history update to before fluid solve
Formerly, history values were updated prior to post-processing, so that all evaluations happened with lhs==lhs0 etc. This patch fixes this, improving VTK and matplotlib results slightly.
1 parent 628a1e6 commit 4b7b925

File tree

1 file changed

+14
-16
lines changed

1 file changed

+14
-16
lines changed

perpendicular-flap/fluid-nutils/fluid.py

Lines changed: 14 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -69,8 +69,8 @@ def main(inflow: 'inflow velocity' = 10,
6969
meshdofs0 = meshdofs
7070
meshdofs00 = meshdofs
7171
meshdofs000 = meshdofs
72-
meshdofs0000 = meshdofs
73-
lhs0 = numpy.zeros(len(ns.ubasis))
72+
lhs = numpy.zeros(len(ns.ubasis))
73+
lhs0 = lhs
7474

7575
# for visualization
7676
bezier = domain.sample('bezier', 2)
@@ -98,7 +98,6 @@ def main(inflow: 'inflow velocity' = 10,
9898

9999
# initialize preCICE
100100
precice_dt = interface.initialize()
101-
dt = min(precice_dt, timestepsize)
102101

103102
# boundary conditions for fluid equations
104103
sqr = domain.boundary['wall,flap'].integral('urel_k urel_k d:x0' @ ns, degree=4)
@@ -132,7 +131,6 @@ def main(inflow: 'inflow velocity' = 10,
132131
# better initial guess: start from Stokes solution, comment out for comparison with other solvers
133132
#res_stokes = domain.integral('(ubasis_ni,j ((u_i,j + u_j,i) rho nu - p δ_ij) + pbasis_n u_k,k) d:x' @ ns, degree=4)
134133
#lhs0 = solver.solve_linear('lhs', res_stokes, constrain=cons, arguments=dict(meshdofs=meshdofs, meshdofs0=meshdofs0, meshdofs00=meshdofs00, meshdofs000=meshdofs000, dt=dt))
135-
lhs00 = lhs0
136134

137135
timestep = 0
138136

@@ -149,9 +147,19 @@ def main(inflow: 'inflow velocity' = 10,
149147

150148
# save checkpoint
151149
if interface.is_action_required(precice.action_write_iteration_checkpoint()):
152-
checkpoint = lhs0, lhs00, timestep, meshdofs0, meshdofs00, meshdofs000, meshdofs0000
150+
checkpoint = timestep, lhs, lhs0, meshdofs, meshdofs0, meshdofs00, meshdofs000
153151
interface.mark_action_fulfilled(precice.action_write_iteration_checkpoint())
154152

153+
# advance variables
154+
timestep += 1
155+
lhs00 = lhs0
156+
lhs0 = lhs
157+
meshdofs0000 = meshdofs000
158+
meshdofs000 = meshdofs00
159+
meshdofs00 = meshdofs0
160+
meshdofs0 = meshdofs
161+
dt = min(precice_dt, timestepsize)
162+
155163
# solve fluid equations
156164
lhs = solver.newton('lhs', res, lhs0=lhs0, constrain=cons,
157165
arguments=dict(lhs0=lhs0, dt=dt, meshdofs=meshdofs, meshdofs0=meshdofs0,
@@ -171,20 +179,10 @@ def main(inflow: 'inflow velocity' = 10,
171179

172180
# do the coupling
173181
precice_dt = interface.advance(dt)
174-
dt = min(precice_dt, timestepsize)
175-
176-
# advance variables
177-
timestep += 1
178-
lhs00 = lhs0
179-
lhs0 = lhs
180-
meshdofs0000 = meshdofs000
181-
meshdofs000 = meshdofs00
182-
meshdofs00 = meshdofs0
183-
meshdofs0 = meshdofs
184182

185183
# read checkpoint if required
186184
if interface.is_action_required(precice.action_read_iteration_checkpoint()):
187-
lhs0, lhs00, timestep, meshdofs0, meshdofs00, meshdofs000, meshdofs0000 = checkpoint
185+
timestep, lhs, lhs0, meshdofs, meshdofs0, meshdofs00, meshdofs000 = checkpoint
188186
interface.mark_action_fulfilled(precice.action_read_iteration_checkpoint())
189187

190188
if interface.is_time_window_complete():

0 commit comments

Comments
 (0)