Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

simplify incompressible solver based on burgers solver #169

Merged
merged 41 commits into from
Aug 8, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
Show all changes
41 commits
Select commit Hold shift + click to select a range
715f6ad
initial commit
zhichen3 Apr 20, 2023
fb29ed7
fix burgers
zhichen3 Apr 21, 2023
15fa38c
simplify incomp_interface.py
zhichen3 Apr 21, 2023
c354dc8
change visual plottings to output both u and v
zhichen3 Apr 21, 2023
af56db8
change plotting
zhichen3 Apr 21, 2023
28a97bf
change figure
zhichen3 Apr 21, 2023
51293a7
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 21, 2023
16b5b5d
fix pylint and flake8
zhichen3 Apr 21, 2023
034f4de
let incompressible inherit from burgers
zhichen3 Apr 21, 2023
0d48b7d
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 21, 2023
02492c1
fix pylint flake8
zhichen3 Apr 21, 2023
00b3cb9
change import
zhichen3 Apr 21, 2023
9702078
move riemann and upwind to burgers_interface to avoid circular import
zhichen3 Apr 21, 2023
5fef42e
move riemann and upwinds functions in incomp_interface to burgers_int…
zhichen3 Apr 21, 2023
0101c94
fix function call
zhichen3 Apr 21, 2023
2a20b56
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 21, 2023
272fef8
change import name
zhichen3 Apr 21, 2023
58cfda2
change riemann to riemann_and_upwind
zhichen3 Apr 22, 2023
47d7af1
change buf=1 to buf=2
zhichen3 Apr 23, 2023
0053b53
split transverse correction to a separate function
zhichen3 Apr 23, 2023
78198f7
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 23, 2023
8e51e5c
apply gradp correction first instead of transverse correction
zhichen3 Apr 23, 2023
dbf9cb5
fix pylint and flake8
zhichen3 Apr 23, 2023
25d0603
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 23, 2023
f165466
fix pylint
zhichen3 Apr 24, 2023
be549b4
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 24, 2023
e2d97c0
fix pylint
zhichen3 Apr 24, 2023
c569a36
fix pylint
zhichen3 Apr 24, 2023
31da292
Merge branch 'main' into burgers-fix
zhichen3 Apr 24, 2023
bf37342
update docs
zhichen3 Apr 24, 2023
de1f108
revert timestep calculation
zhichen3 Apr 24, 2023
da501fd
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 24, 2023
e721f4a
Merge branch 'main' into burgers-fix
zingale Apr 26, 2023
57b73da
Merge branch 'main' into incompressible_simplify
zhichen3 Apr 27, 2023
047ad09
revert to apply transverse terms first
zhichen3 Apr 27, 2023
8074edb
more revert
zhichen3 Apr 27, 2023
104363b
break up lines
zhichen3 Apr 27, 2023
10d00c6
Merge branch 'burgers-fix' of github.com:zhichen3/pyro2 into burgers-fix
zhichen3 Apr 27, 2023
b2d1d5b
Merge branch 'burgers-fix' into incompressible_simplify
zhichen3 Apr 27, 2023
d68599c
break up lines
zhichen3 Apr 27, 2023
bc17d46
Merge branch 'main' into incompressible_simplify
zhichen3 Aug 4, 2023
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
196 changes: 76 additions & 120 deletions pyro/incompressible/incomp_interface.py
Original file line number Diff line number Diff line change
Expand Up @@ -32,12 +32,29 @@ def mac_vels(grid, dt,

# get the full u and v left and right states (including transverse
# terms) on both the x- and y-interfaces
# pylint: disable-next=unused-variable
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = get_interface_states(grid, dt,
u, v,
ldelta_ux, ldelta_vx,
ldelta_uy, ldelta_vy,
gradp_x, gradp_y)

u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
burgers_interface.get_interface_states(grid, dt,
u, v,
ldelta_ux, ldelta_vx,
ldelta_uy, ldelta_vy)

# apply transverse terms on both x- and y- interfaces
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
burgers_interface.apply_transverse_corrections(grid, dt,
u_xl, u_xr,
u_yl, u_yr,
v_xl, v_xr,
v_yl, v_yr)

# apply pressure gradient correction terms to the interface state
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
apply_gradp_corrections(dt,
u_xl, u_xr,
u_yl, u_yr,
v_xl, v_xr,
v_yl, v_yr,
gradp_x, gradp_y)

# Riemann problem -- this follows Burger's equation. We don't use
# any input velocity for the upwinding. Also, we only care about
Expand Down Expand Up @@ -82,13 +99,29 @@ def states(grid, dt,
x and y velocities predicted to the interfaces
"""

# get the full u and v left and right states (including transverse
# terms) on both the x- and y-interfaces
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = get_interface_states(grid, dt,
u, v,
ldelta_ux, ldelta_vx,
ldelta_uy, ldelta_vy,
gradp_x, gradp_y)
# get the full u and v left and right states without transverse terms and gradp
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
burgers_interface.get_interface_states(grid, dt,
u, v,
ldelta_ux, ldelta_vx,
ldelta_uy, ldelta_vy)

# apply transverse terms on both x- and y- interfaces
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
burgers_interface.apply_transverse_corrections(grid, dt,
u_xl, u_xr,
u_yl, u_yr,
v_xl, v_xr,
v_yl, v_yr)

# apply pressure gradient correction terms to the interface state
u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr = \
apply_gradp_corrections(dt,
u_xl, u_xr,
u_yl, u_yr,
v_xl, v_xr,
v_yl, v_yr,
gradp_x, gradp_y)

# upwind using the MAC velocity to determine which state exists on
# the interface
Expand All @@ -100,128 +133,51 @@ def states(grid, dt,
return u_xint, v_xint, u_yint, v_yint


def get_interface_states(grid, dt,
u, v,
ldelta_ux, ldelta_vx,
ldelta_uy, ldelta_vy,
gradp_x, gradp_y):
def apply_gradp_corrections(dt,
u_xl, u_xr,
u_yl, u_yr,
v_xl, v_xr,
v_yl, v_yr,
gradp_x, gradp_y):
r"""
Compute the unsplit predictions of u and v on both the x- and
y-interfaces. This includes the transverse terms.

Parameters
----------
grid : Grid2d
The grid object
dt : float
The timestep
u, v : ndarray
x-velocity and y-velocity
ldelta_ux, ldelta_uy: ndarray
Limited slopes of the x-velocity in the x and y directions
ldelta_vx, ldelta_vy: ndarray
Limited slopes of the y-velocity in the x and y directions
gradp_x, gradp_y : ndarray
Pressure gradients in the x and y directions

u_xl, u_xr : ndarray ndarray
left and right states of x-velocity in x-interface.
u_yl, u_yr : ndarray ndarray
left and right states of x-velocity in y-interface.
v_xl, v_xr : ndarray ndarray
left and right states of y-velocity in x-interface.
v_yl, u_yr : ndarray ndarray
left and right states of y-velocity in y-interface.
Returns
-------
out : ndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray
unsplit predictions of u and v on both the x- and
y-interfaces
correct the interface states of the left and right states of u and v on
both the x- and y-interfaces interface states with the pressure gradient
terms.
"""

u_xl = grid.scratch_array()
u_xr = grid.scratch_array()
u_yl = grid.scratch_array()
u_yr = grid.scratch_array()

v_xl = grid.scratch_array()
v_xr = grid.scratch_array()
v_yl = grid.scratch_array()
v_yr = grid.scratch_array()

# first predict u and v to both interfaces, considering only the normal
# part of the predictor. These are the 'hat' states.

dtdx = dt / grid.dx
dtdy = dt / grid.dy

# u on x-edges
u_xl.ip(1, buf=2)[:, :] = u.v(buf=2) + \
0.5 * (1.0 - dtdx * u.v(buf=2)) * ldelta_ux.v(buf=2)
u_xr.v(buf=2)[:, :] = u.v(buf=2) - \
0.5 * (1.0 + dtdx * u.v(buf=2)) * ldelta_ux.v(buf=2)

# v on x-edges
v_xl.ip(1, buf=2)[:, :] = v.v(buf=2) + \
0.5 * (1.0 - dtdx * u.v(buf=2)) * ldelta_vx.v(buf=2)
v_xr.v(buf=2)[:, :] = v.v(buf=2) - \
0.5 * (1.0 + dtdx * u.v(buf=2)) * ldelta_vx.v(buf=2)

# u on y-edges
u_yl.jp(1, buf=2)[:, :] = u.v(buf=2) + \
0.5 * (1.0 - dtdy * v.v(buf=2)) * ldelta_uy.v(buf=2)
u_yr.v(buf=2)[:, :] = u.v(buf=2) - \
0.5 * (1.0 + dtdy * v.v(buf=2)) * ldelta_uy.v(buf=2)

# v on y-edges
v_yl.jp(1, buf=2)[:, :] = v.v(buf=2) + \
0.5 * (1.0 - dtdy * v.v(buf=2)) * ldelta_vy.v(buf=2)
v_yr.v(buf=2)[:, :] = v.v(buf=2) - \
0.5 * (1.0 + dtdy * v.v(buf=2)) * ldelta_vy.v(buf=2)

# now get the normal advective velocities on the interfaces by solving
# the Riemann problem.
uhat_adv = burgers_interface.riemann(grid, u_xl, u_xr)
vhat_adv = burgers_interface.riemann(grid, v_yl, v_yr)

# now that we have the advective velocities, upwind the left and right
# states using the appropriate advective velocity.

# on the x-interfaces, we upwind based on uhat_adv
u_xint = burgers_interface.upwind(grid, u_xl, u_xr, uhat_adv)
v_xint = burgers_interface.upwind(grid, v_xl, v_xr, uhat_adv)

# on the y-interfaces, we upwind based on vhat_adv
u_yint = burgers_interface.upwind(grid, u_yl, u_yr, vhat_adv)
v_yint = burgers_interface.upwind(grid, v_yl, v_yr, vhat_adv)

# at this point, these states are the `hat' states -- they only
# considered the normal to the interface portion of the predictor.

ubar = grid.scratch_array()
vbar = grid.scratch_array()

ubar.v(buf=2)[:, :] = 0.5 * (uhat_adv.v(buf=2) + uhat_adv.ip(1, buf=2))
vbar.v(buf=2)[:, :] = 0.5 * (vhat_adv.v(buf=2) + vhat_adv.jp(1, buf=2))

# v du/dy is the transerse term for the u states on x-interfaces
vu_y = grid.scratch_array()
vu_y.v(buf=2)[:, :] = vbar.v(buf=2) * (u_yint.jp(1, buf=2) - u_yint.v(buf=2))

u_xl.ip(1, buf=2)[:, :] += -0.5 * dtdy * vu_y.v(buf=2) - 0.5 * dt * gradp_x.v(buf=2)
u_xr.v(buf=2)[:, :] += -0.5 * dtdy * vu_y.v(buf=2) - 0.5 * dt * gradp_x.v(buf=2)

# v dv/dy is the transverse term for the v states on x-interfaces
vv_y = grid.scratch_array()
vv_y.v(buf=2)[:, :] = vbar.v(buf=2) * (v_yint.jp(1, buf=2) - v_yint.v(buf=2))

v_xl.ip(1, buf=2)[:, :] += -0.5 * dtdy * vv_y.v(buf=2) - 0.5 * dt * gradp_y.v(buf=2)
v_xr.v(buf=2)[:, :] += -0.5 * dtdy * vv_y.v(buf=2) - 0.5 * dt * gradp_y.v(buf=2)
# Apply pressure gradient correction terms

# u dv/dx is the transverse term for the v states on y-interfaces
uv_x = grid.scratch_array()
uv_x.v(buf=2)[:, :] = ubar.v(buf=2) * (v_xint.ip(1, buf=2) - v_xint.v(buf=2))
# transverse term for the u states on x-interfaces
u_xl.ip(1, buf=2)[:, :] += - 0.5 * dt * gradp_x.v(buf=2)
u_xr.v(buf=2)[:, :] += - 0.5 * dt * gradp_x.v(buf=2)

v_yl.jp(1, buf=2)[:, :] += -0.5 * dtdx * uv_x.v(buf=2) - 0.5 * dt * gradp_y.v(buf=2)
v_yr.v(buf=2)[:, :] += -0.5 * dtdx * uv_x.v(buf=2) - 0.5 * dt * gradp_y.v(buf=2)
# transverse term for the v states on x-interfaces
v_xl.ip(1, buf=2)[:, :] += - 0.5 * dt * gradp_y.v(buf=2)
v_xr.v(buf=2)[:, :] += - 0.5 * dt * gradp_y.v(buf=2)

# u du/dx is the transverse term for the u states on y-interfaces
uu_x = grid.scratch_array()
uu_x.v(buf=2)[:, :] = ubar.v(buf=2) * (u_xint.ip(1, buf=2) - u_xint.v(buf=2))
# transverse term for the v states on y-interfaces
v_yl.jp(1, buf=2)[:, :] += - 0.5 * dt * gradp_y.v(buf=2)
v_yr.v(buf=2)[:, :] += - 0.5 * dt * gradp_y.v(buf=2)

u_yl.jp(1, buf=2)[:, :] += -0.5 * dtdx * uu_x.v(buf=2) - 0.5 * dt * gradp_x.v(buf=2)
u_yr.v(buf=2)[:, :] += -0.5 * dtdx * uu_x.v(buf=2) - 0.5 * dt * gradp_x.v(buf=2)
# transverse term for the u states on y-interfaces
u_yl.jp(1, buf=2)[:, :] += - 0.5 * dt * gradp_x.v(buf=2)
u_yr.v(buf=2)[:, :] += - 0.5 * dt * gradp_x.v(buf=2)

return u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr
26 changes: 3 additions & 23 deletions pyro/incompressible/simulation.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,14 +4,15 @@
import numpy as np

import pyro.mesh.array_indexer as ai
from pyro.burgers import Simulation as burgers_simulation
from pyro.incompressible import incomp_interface
from pyro.mesh import patch, reconstruction
from pyro.multigrid import MG
from pyro.particles import particles
from pyro.simulation_null import NullSimulation, bc_setup, grid_setup
from pyro.simulation_null import bc_setup, grid_setup


class Simulation(NullSimulation):
class Simulation(burgers_simulation):

def initialize(self):
"""
Expand Down Expand Up @@ -51,27 +52,6 @@ def initialize(self):
problem = importlib.import_module(f"pyro.incompressible.problems.{self.problem_name}")
problem.init_data(self.cc_data, self.rp)

def method_compute_timestep(self):
"""
The timestep() function computes the advective timestep
(CFL) constraint. The CFL constraint says that information
cannot propagate further than one zone per timestep.

We use the driver.cfl parameter to control what fraction of the CFL
step we actually take.
"""

cfl = self.rp.get_param("driver.cfl")

u = self.cc_data.get_var("x-velocity")
v = self.cc_data.get_var("y-velocity")

# the timestep is min(dx/|u|, dy|v|)
xtmp = self.cc_data.grid.dx/(abs(u))
ytmp = self.cc_data.grid.dy/(abs(v))

self.dt = cfl*float(min(xtmp.min(), ytmp.min()))

def preevolve(self):
"""
preevolve is called before we being the timestepping loop. For
Expand Down