pyro.burgers_viscous package#

The pyro viscous burgers solver. This implements a second-order, unsplit method for viscous burgers equations.

Subpackages#

Submodules#

pyro.burgers_viscous.interface module#

pyro.burgers_viscous.interface.apply_diffusion_corrections(grid, dt, eps, u, v, u_xl, u_xr, u_yl, u_yr, v_xl, v_xr, v_yl, v_yr)[source]#

Apply diffusion correction term to the interface state

\[u_t + u u_x + v u_y = eps (u_xx + u_yy) v_t + u v_x + v v_y = eps (v_xx + v_yy)\]

We use a second-order (piecewise linear) unsplit Godunov method (following Colella 1990).

Our convection is that the fluxes are going to be defined on the left edge of the computational zones:

 |             |             |             |
 |             |             |             |
-+------+------+------+------+------+------+--
 |     i-1     |      i      |     i+1     |

          a_l,i  a_r,i   a_l,i+1

a_r,i and a_l,i+1 are computed using the information in zone i,j.

Parameters:
gridGrid2d

The grid object

dtfloat

The timestep

epsfloat

The viscosity

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

unsplit predictions of the left and right states of u and v on both the x- and y-interfaces along with diffusion correction terms.

pyro.burgers_viscous.interface.diffuse(my_data, rp, dt, scalar_name, A)[source]#

A routine to solve the Helmhotlz equation with constant coefficient and update the state.

(a + b lap) phi = f

using Crank-Nicolson discretization with multigrid V-cycle.

Parameters:
my_dataCellCenterData2d object

The data object containing the grid and advective scalar that we are advecting.

rpRuntimeParameters object

The runtime parameters for the simulation

dtfloat

The timestep we are advancing through.

scalar_namestr

The name of the variable contained in my_data that we are advecting

A: ndarray

The advective source term for diffusion

Returns:
outndarray (solution of the Helmholtz equation)
pyro.burgers_viscous.interface.get_lap(grid, a)[source]#
Parameters:
grid: grid object
a: ndarray

the variable that we want to find the laplacian of

Returns:
outndarray (laplacian of state a)

pyro.burgers_viscous.simulation module#

class pyro.burgers_viscous.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: Simulation

evolve()[source]#

Evolve the viscous burgers equation through one timestep.