pyro.advection package

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

Subpackages

Submodules

pyro.advection.advective_fluxes module

pyro.advection.advective_fluxes.unsplit_fluxes(my_data, rp, dt, scalar_name, interface)[source]

Construct the fluxes through the interfaces for the linear advection equation:

\[a_t + u a_x + v a_y = 0\]

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

In the pure advection case, there is no Riemann problem we need to solve – we just simply do upwinding. So there is only one ‘state’ at each interface, and the zone the information comes from depends on the sign of the velocity.

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:
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

Returns:
outndarray, ndarray

The fluxes on the x- and y-interfaces

pyro.advection.simulation module

class pyro.advection.simulation.Simulation(solver_name, problem_name, rp, timers=None, data_class=<class 'pyro.mesh.patch.CellCenterData2d'>)[source]

Bases: NullSimulation

dovis()[source]

Do runtime visualization.

evolve()[source]

Evolve the linear advection equation through one timestep. We only consider the “density” variable in the CellCenterData2d object that is part of the Simulation.

initialize()[source]

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

method_compute_timestep()[source]

Compute the advective timestep (CFL) constraint. We use the driver.cfl parameter to control what fraction of the CFL step we actually take.