advection_weno package¶
The pyro advection solver. This implements a finite difference Lax-Friedrichs flux split method with WENO reconstruction based on Shu’s review from 1998 (https://ntrs.nasa.gov/archive/nasa/casi.ntrs.nasa.gov/19980007543.pdf) although the notation more follows Gerolymos et al (https://doi.org/10.1016/j.jcp.2009.07.039).
Most of the code is taken from advection_rk and toy-conslaw.
The general flow of the solver when invoked through pyro.py is:
create grid
initial conditions
main loop
- fill ghost cells
- compute dt
- compute fluxes
- conservative update
- output
Subpackages¶
Submodules¶
advection_weno.fluxes module¶
-
advection_weno.fluxes.
fluxes
(my_data, rp, dt)[source]¶ Construct the fluxes through the interfaces for the linear advection equation
\[a_t + u a_x + v a_y = 0\]We use a high-order flux split WENO method to construct the interface fluxes. No Riemann problems are solved. The Lax-Friedrichs flux split will probably make it diffusive; the lack of a transverse solver also will not help.
Parameters: - my_data : CellCenterData2d object
The data object containing the grid and advective scalar that we are advecting.
- rp : RuntimeParameters object
The runtime parameters for the simulation
- dt : float
The timestep we are advancing through.
- scalar_name : str
The name of the variable contained in my_data that we are advecting
Returns: - out : ndarray, ndarray
The fluxes on the x- and y-interfaces
-
advection_weno.fluxes.
fvs
(q, order, u, alpha)[source]¶ Perform Flux-Vector-Split (LF) finite differencing using WENO in 1d.
Parameters: - q : np array
input data with at least order+1 ghost zones
- order : int
WENO order (k)
- u : float
Advection velocity in this direction
- alpha : float
Maximum characteristic speed
Returns: - f : np array
flux
advection_weno.simulation module¶
-
class
advection_weno.simulation.
Simulation
(solver_name, problem_name, rp, timers=None, data_class=<class 'mesh.patch.CellCenterData2d'>)[source]¶ Bases:
advection.simulation.Simulation
-
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.
-