pyro.burgers package#

The pyro inviscid burgers solver. This implements a second-order, unsplit method for inviscid burgers equations based on the Colella 1990 paper.

Subpackages#

Submodules#

pyro.burgers.burgers_interface module#

pyro.burgers.burgers_interface.apply_transverse_corrections(grid, dt, u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr)[source]#
Parameters:
gridGrid2d

The grid object

dtfloat

The timestep

u_xl, u_xrndarray ndarray

left and right states of x-velocity in x-interface.

u_yl, u_yrndarray ndarray

left and right states of x-velocity in y-interface.

v_xl, v_xrndarray ndarray

left and right states of y-velocity in x-interface.

v_yl, u_yrndarray ndarray

left and right states of y-velocity in y-interface.

Returns
——-
outndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray

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 transverse terms.

pyro.burgers.burgers_interface.construct_unsplit_fluxes(grid, u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr)[source]#

Construct the interface fluxes for the burgers equation:

\[u_t + u u_x + v u_y = 0 v_t + u v_x + v v_y = 0\]
Parameters:
gridGrid2d

The grid object

u_xl, u_xrndarray ndarray

left and right states of x-velocity in x-interface.

u_yl, u_yrndarray ndarray

left and right states of x-velocity in y-interface.

v_xl, v_xrndarray ndarray

left and right states of y-velocity in x-interface.

v_yl, u_yrndarray ndarray

left and right states of y-velocity in y-interface.

——-
Returns
——-
outndarray, ndarray, ndarray, ndarray

The u,v fluxes on the x- and y-interfaces

pyro.burgers.burgers_interface.get_interface_states(grid, dt, u, v, ldelta_ux, ldelta_vx, ldelta_uy, ldelta_vy)[source]#

Construct the interface states for the burgers equation:

\[u_t + u u_x + v u_y = 0 v_t + u v_x + v v_y = 0\]

Compute the unsplit predictions of u and v on both the x- and y-interfaces.

Parameters:
gridGrid2d

The grid object

dtfloat

The timestep

u, vndarray

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

Returns:
outndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray, ndarray

get the predictions of the left and right states of u and v on both the x- and y-interfaces interface states without the transverse terms.

pyro.burgers.burgers_interface.riemann(grid, q_l, q_r)[source]#

Solve the Burger’s Riemann problem given the input left and right states and return the state on the interface.

This uses the expressions from Almgren, Bell, and Szymczak 1996.

Parameters:
gridGrid2d

The grid object

q_l, q_rndarray

left and right states

Returns:
outndarray

Interface state

pyro.burgers.burgers_interface.riemann_and_upwind(grid, q_l, q_r)[source]#

First solve the Riemann problem given q_l and q_r to give the velocity on the interface and: use this velocity to upwind to determine the state (q_l, q_r, or a mix) on the interface).

This differs from upwind, above, in that we don’t take in a velocity to upwind with).

Parameters:
gridGrid2d

The grid object

q_l, q_rndarray

left and right states

Returns:
outndarray

Upwinded state

pyro.burgers.burgers_interface.upwind(grid, q_l, q_r, s)[source]#

upwind the left and right states based on the specified input velocity, s. The resulting interface state is q_int

Parameters:
gridGrid2d

The grid object

q_l, q_rndarray

left and right states

sndarray

velocity

Returns:
outndarray

Upwinded state

pyro.burgers.simulation module#

class pyro.burgers.simulation.Simulation(solver_name, problem_name, problem_func, rp, *, problem_finalize_func=None, problem_source_func=None, timers=None, data_class=<class 'pyro.mesh.patch.CellCenterData2d'>)[source]#

Bases: NullSimulation

dovis()[source]#

Do runtime visualization

evolve()[source]#

Evolve the burgers equation through one timestep.

initialize()[source]#

Initialize the grid and variables for advection and set the initial conditions for the chosen problem.

method_compute_timestep()[source]#

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.