compressible package¶
The pyro compressible hydrodynamics solver. This implements the second-order (piecewise-linear), unsplit method of Colella 1990.
Subpackages¶
- compressible.problems package
- Submodules
- compressible.problems.acoustic_pulse module
- compressible.problems.advect module
- compressible.problems.bubble module
- compressible.problems.hse module
- compressible.problems.kh module
- compressible.problems.logo module
- compressible.problems.quad module
- compressible.problems.ramp module
- compressible.problems.rt module
- compressible.problems.rt2 module
- compressible.problems.sedov module
- compressible.problems.sod module
- compressible.problems.test module
Submodules¶
compressible.BC module¶
compressible-specific boundary conditions. Here, in particular, we implement an HSE BC in the vertical direction.
Note: the pyro BC routines operate on a single variable at a time, so some work will necessarily be repeated.
Also note: we may come in here with the aux_data (source terms), so we’ll do a special case for them
-
compressible.BC.
user
(bc_name, bc_edge, variable, ccdata)[source]¶ A hydrostatic boundary. This integrates the equation of HSE into the ghost cells to get the pressure and density under the assumption that the specific internal energy is constant.
Upon exit, the ghost cells for the input variable will be set
Parameters: bc_name : {‘hse’}
The descriptive name for the boundary condition – this allows for pyro to have multiple types of user-supplied boundary conditions. For this module, it needs to be ‘hse’.
bc_edge : {‘ylb’, ‘yrb’}
The boundary to update: ylb = lower y boundary; yrb = upper y boundary.
variable : {‘density’, ‘x-momentum’, ‘y-momentum’, ‘energy’}
The variable whose ghost cells we are filling
ccdata : CellCenterData2d object
The data object
compressible.derives module¶
compressible.eos module¶
This is a gamma-law equation of state: p = rho e (gamma - 1), where gamma is the constant ratio of specific heats.
-
compressible.eos.
dens
(gamma, pres, eint)[source]¶ Given the pressure and the specific internal energy, return the density
Parameters: gamma : float
The ratio of specific heats
pres : float
The pressure
eint : float
The specific internal energy
Returns: out : float
The density
compressible.simulation module¶
-
class
compressible.simulation.
Simulation
(solver_name, problem_name, rp, timers=None, data_class=<class 'mesh.patch.CellCenterData2d'>)[source]¶ Bases:
simulation_null.NullSimulation
The main simulation class for the corner transport upwind compressible hydrodynamics solver
-
initialize
(extra_vars=None, ng=4)[source]¶ Initialize the grid and variables for compressible flow and set the initial conditions for the chosen problem.
-
-
class
compressible.simulation.
Variables
(myd)[source]¶ Bases:
object
a container class for easy access to the different compressible variable by an integer key
compressible.unsplit_fluxes module¶
Implementation of the Colella 2nd order unsplit Godunov scheme. This is a 2-dimensional implementation only. We assume that the grid is uniform, but it is relatively straightforward to relax this assumption.
There are several different options for this solver (they are all discussed in the Colella paper).
- limiter: 0 = no limiting; 1 = 2nd order MC limiter; 2 = 4th order MC limiter
- riemann: HLLC or CGF (for Colella, Glaz, and Freguson solver)
- use_flattening: set to 1 to use the multidimensional flattening at shocks
- delta, z0, z1: flattening parameters (we use Colella 1990 defaults)
The grid indices look like:
j+3/2--+---------+---------+---------+
| | | |
j+1 _| | | |
| | | |
| | | |
j+1/2--+---------XXXXXXXXXXX---------+
| X X |
j _| X X |
| X X |
| X X |
j-1/2--+---------XXXXXXXXXXX---------+
| | | |
j-1 _| | | |
| | | |
| | | |
j-3/2--+---------+---------+---------+
| | | | | | |
i-1 i i+1
i-3/2 i-1/2 i+1/2 i+3/2
We wish to solve
we want U_{i+1/2}^{n+1/2} – the interface values that are input to the Riemann problem through the faces for each zone.
Taylor expanding yields:
n+1/2 dU dU
U = U + 0.5 dx -- + 0.5 dt --
i+1/2,j,L i,j dx dt
dU dF^x dF^y
= U + 0.5 dx -- - 0.5 dt ( ---- + ---- - H )
i,j dx dx dy
dU dF^x dF^y
= U + 0.5 ( dx -- - dt ---- ) - 0.5 dt ---- + 0.5 dt H
i,j dx dx dy
dt dU dF^y
= U + 0.5 dx ( 1 - -- A^x ) -- - 0.5 dt ---- + 0.5 dt H
i,j dx dx dy
dt _ dF^y
= U + 0.5 ( 1 - -- A^x ) DU - 0.5 dt ---- + 0.5 dt H
i,j dx dy
+----------+-----------+ +----+----+ +---+---+
| | |
this is the monotonized this is the source term
central difference term transverse
flux term
There are two components, the central difference in the normal to the interface, and the transverse flux difference. This is done for the left and right sides of all 4 interfaces in a zone, which are then used as input to the Riemann problem, yielding the 1/2 time interface values:
n+1/2
U
i+1/2,j
Then, the zone average values are updated in the usual finite-volume way:
n+1 n dt x n+1/2 x n+1/2
U = U + -- { F (U ) - F (U ) }
i,j i,j dx i-1/2,j i+1/2,j
dt y n+1/2 y n+1/2
+ -- { F (U ) - F (U ) }
dy i,j-1/2 i,j+1/2
Updating U_{i,j}:
- We want to find the state to the left and right (or top and bottom) of each interface, ex. U_{i+1/2,j,[lr]}^{n+1/2}, and use them to solve a Riemann problem across each of the four interfaces.
- U_{i+1/2,j,[lr]}^{n+1/2} is comprised of two parts, the computation
of the monotonized central differences in the normal direction
(eqs. 2.8, 2.10) and the computation of the transverse derivatives,
which requires the solution of a Riemann problem in the transverse
direction (eqs. 2.9, 2.14).
- the monotonized central difference part is computed using the primitive variables.
- We compute the central difference part in both directions before doing the transverse flux differencing, since for the high-order transverse flux implementation, we use these as the input to the transverse Riemann problem.
-
compressible.unsplit_fluxes.
unsplit_fluxes
(my_data, my_aux, rp, ivars, solid, tc, dt)[source]¶ unsplitFluxes returns the fluxes through the x and y interfaces by doing an unsplit reconstruction of the interface values and then solving the Riemann problem through all the interfaces at once
currently we assume a gamma-law EOS
The runtime parameter grav is assumed to be the gravitational acceleration in the y-direction
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
vars : Variables object
The Variables object that tells us which indices refer to which variables
tc : TimerCollection object
The timers we are using to profile
dt : float
The timestep we are advancing through.
Returns: out : ndarray, ndarray
The fluxes on the x- and y-interfaces