pyro.mesh package

This is the general mesh module for pyro. It implements everything necessary to work with finite-volume data.

Submodules

pyro.mesh.array_indexer module

An array class that has methods supporting the type of stencil operations we see in finite-difference methods, like i+1, i-1, etc.

class pyro.mesh.array_indexer.ArrayIndexer(d, grid=None)[source]

Bases: ndarray

a class that wraps the data region of a single cell-centered data array (d) and allows us to easily do array operations like d[i+1,j] using the ip() method.

copy(order='C')[source]

make a copy of the array, defined on the same grid

fill_ghost(n=0, bc=None)[source]

Fill the boundary conditions. This operates on a single component, n. We do periodic, reflect-even, reflect-odd, and outflow

We need a BC object to tell us what BC type on each boundary.

ip(shift, buf=0, n=0, s=1)[source]

return a view of the data shifted by shift in the x direction. By default the view is the same size as the valid region, but the buf can specify how many ghost cells on each side to include. The component is n and s is the stride

ip_jp(ishift, jshift, buf=0, n=0, s=1)[source]

return a view of the data shifted by ishift in the x direction and jshift in the y direction. By default the view is the same size as the valid region, but the buf can specify how many ghost cells on each side to include. The component is n and s is the stride

is_asymmetric(nodal=False, tol=1e-14)[source]

return True is the data is left-right asymmetric (to the tolerance tol)—e.g, the sign flips. For node-centered data, set nodal=True

is_symmetric(nodal=False, tol=1e-14, asymmetric=False)[source]

return True is the data is left-right symmetric (to the tolerance tol) For node-centered data, set nodal=True

jp(shift, buf=0, n=0, s=1)[source]

return a view of the data shifted by shift in the y direction. By default the view is the same size as the valid region, but the buf can specify how many ghost cells on each side to include. The component is n and s is the stride

lap(n=0, buf=0)[source]

return the 5-point Laplacian

norm(n=0)[source]

find the norm of the quantity (index n) defined on the same grid, in the domain’s valid region

pretty_print(n=0, fmt=None, show_ghost=True)[source]

Print out a small dataset to the screen with the ghost cells a different color, to make things stand out

v(buf=0, n=0, s=1)[source]

return a view of the valid data region for component n, with stride s, and a buffer of ghost cells given by buf

class pyro.mesh.array_indexer.ArrayIndexerFC(d, idir, grid=None)[source]

Bases: ArrayIndexer

a class that wraps the data region of a single face-centered data array (d) and allows us to easily do array operations like d[i+1,j] using the ip() method.

copy(order='C')[source]

make a copy of the array, defined on the same grid

fill_ghost(n=0, bc=None)[source]

Fill the boundary conditions. This operates on a single component, n. We do periodic, reflect-even, reflect-odd, and outflow

We need a BC object to tell us what BC type on each boundary.

ip_jp(ishift, jshift, buf=0, n=0, s=1)[source]

return a view of the data shifted by ishift in the x direction and jshift in the y direction. By default the view is the same size as the valid region, but the buf can specify how many ghost cells on each side to include. The component is n and s is the stride

is_asymmetric(nodal=False, tol=1e-14)[source]

return True is the data is left-right asymmetric (to the tolerance tol)—e.g, the sign flips. For node-centered data, set nodal=True

is_symmetric(nodal=False, tol=1e-14, asymmetric=False)[source]

return True is the data is left-right symmetric (to the tolerance tol) For node-centered data, set nodal=True

lap(n=0, buf=0)[source]

return the 5-point Laplacian

norm(n=0)[source]

find the norm of the quantity (index n) defined on the same grid, in the domain’s valid region

pretty_print(n=0, fmt=None, show_ghost=True)[source]

Print out a small dataset to the screen with the ghost cells a different color, to make things stand out

pyro.mesh.boundary module

Methods to manage boundary conditions

class pyro.mesh.boundary.BC(xlb='outflow', xrb='outflow', ylb='outflow', yrb='outflow', xl_func=None, xr_func=None, yl_func=None, yr_func=None, grid=None, odd_reflect_dir='')[source]

Bases: object

Boundary condition container – hold the BCs on each boundary for a single variable.

For Neumann and Dirichlet BCs, a function callback can be stored for inhomogeous BCs. This function should provide the value on the physical boundary (not cell center). This is evaluated on the relevant edge when the __init__ routine is called. For this reason, you need to pass in a grid object. Note: this only ensures that the first ghost cells is consistent with the BC value.

class pyro.mesh.boundary.BCProp(xl_prop, xr_prop, yl_prop, yr_prop)[source]

Bases: object

A simple container to hold properties of the boundary conditions.

pyro.mesh.boundary.bc_is_solid(bc)[source]

return a container class indicating which boundaries are solid walls

pyro.mesh.boundary.define_bc(bc_type, function, is_solid=False)[source]

use this to extend the types of boundary conditions supported on a solver-by-solver basis. Here we pass in the reference to a function that can be called with the data that needs to be filled. is_solid indicates whether it should be interpreted as a solid wall (no flux through the BC)”

pyro.mesh.fv module

This implements support for 4th-order accurate finite-volume data by adding support for converting between cell averages and centers.

class pyro.mesh.fv.FV2d(grid, dtype=<class 'numpy.float64'>)[source]

Bases: CellCenterData2d

this is a finite-volume grid. We expect the data to represent cell-averages, and do operations to 4th order. This assumes dx = dy

from_centers(name)[source]

treat the stored data as if it lives at cell-centers and convert it to an average

to_centers(name)[source]

convert variable name from an average to cell-centers

pyro.mesh.integration module

A generic Runge-Kutta type integrator for integrating CellCenterData2d. We support a generic Butcher tableau for explicit the Runge-Kutta update:

0   |
c_2 | a_21
c_3 | a_31 a_32
:   |  :        .
:   |  :          .
c_s | a_s1 a_s2 ... a_s,s-1
----+---------------------------
    | b_1  b_2  ... b_{s-1}  b_s

the update is:

y_{n+1} = y_n + dt sum_{i=1}^s {b_i k_i}

and the s increment is:

k_s = f(t + c_s dt, y_n + dt (a_s1 k1 + a_s2 k2 + ... + a_s,s-1 k_{s-1})
class pyro.mesh.integration.RKIntegrator(t, dt, method='RK4')[source]

Bases: object

the integration class for CellCenterData2d, supporting RK integration

compute_final_update()[source]

this constructs the final t + dt update, overwriting the initial data

get_stage_start(istage)[source]

get the starting conditions (a CellCenterData2d object) for stage istage

nstages()[source]

return the number of stages

set_start(start)[source]

store the starting conditions (should be a CellCenterData2d object)

store_increment(istage, k_stage)[source]

store the increment for stage istage – this should not have a dt weighting

pyro.mesh.patch module

The patch module defines the classes necessary to describe finite-volume data and the grid that it lives on.

Typical usage:

  • create the grid:

    grid = Grid2d(nx, ny)
    
  • create the data that lives on that grid:

    data = CellCenterData2d(grid)
    
    bc = BC(xlb="reflect", xrb="reflect",
           ylb="outflow", yrb="outflow")
    data.register_var("density", bc)
    ...
    
    data.create()
    
  • initialize some data:

    dens = data.get_var("density")
    dens[:, :] = ...
    
  • fill the ghost cells:

    data.fill_BC("density")
    
class pyro.mesh.patch.CellCenterData2d(grid, dtype=<class 'numpy.float64'>)[source]

Bases: object

A class to define cell-centered data that lives on a grid. A CellCenterData2d object is built in a multi-step process before it can be used.

  • Create the object. We pass in a grid object to describe where the data lives:

    my_data = patch.CellCenterData2d(myGrid)
    
  • Register any variables that we expect to live on this patch. Here BC describes the boundary conditions for that variable:

    my_data.register_var('density', BC)
    my_data.register_var('x-momentum', BC)
    ...
    
  • Register any auxiliary data – these are any parameters that are needed to interpret the data outside of the simulation (for example, the gamma for the equation of state):

    my_data.set_aux(keyword, value)
    
  • Finish the initialization of the patch:

    my_data.create()
    

This last step actually allocates the storage for the state variables. Once this is done, the patch is considered to be locked. New variables cannot be added.

add_derived(func)[source]

Register a function to compute derived variable

Parameters:
funcfunction

A function to call to derive the variable. This function should take two arguments, a CellCenterData2d object and a string variable name (or list of variables)

add_ivars(ivars)[source]

Add ivars

create()[source]

Called after all the variables are registered and allocates the storage for the state data.

fill_BC(name)[source]

Fill the boundary conditions. This operates on a single state variable at a time, to allow for maximum flexibility.

We do periodic, reflect-even, reflect-odd, and outflow

Each variable name has a corresponding BC stored in the CellCenterData2d object – we refer to this to figure out the action to take at each boundary.

Parameters:
namestr

The name of the variable for which to fill the BCs.

fill_BC_all()[source]

Fill boundary conditions on all variables.

get_aux(keyword)[source]

Get the auxiliary data associated with keyword

Parameters:
keywordstr

The name of the auxiliary data to access

Returns:
outvariable type

The value corresponding to the keyword

get_var(name)[source]

Return a data array for the variable described by name. Stored variables will be checked first, and then any derived variables will be checked.

For a stored variable, changes made to this are automatically reflected in the CellCenterData2d object.

Parameters:
namestr

The name of the variable to access

Returns:
outndarray

The array of data corresponding to the variable name

get_var_by_index(n)[source]

Return a data array for the variable with index n in the data array. Any changes made to this are automatically reflected in the CellCenterData2d object.

Parameters:
nint

The index of the variable to access

Returns:
outndarray

The array of data corresponding to the index

get_vars()[source]

Return the entire data array. Any changes made to this are automatically reflected in the CellCenterData2d object.

Returns:
outndarray

The array of data

max(name, ng=0)[source]

return the maximum of the variable name in the domain’s valid region

min(name, ng=0)[source]

return the minimum of the variable name in the domain’s valid region

pretty_print(var, fmt=None)[source]

print out the contents of the data array with pretty formatting indicating where ghost cells are.

prolong(varname)[source]

Prolong the data in the current (coarse) grid to a finer (factor of 2 finer) grid. Return an array with the resulting data (and same number of ghostcells). Only the data for the variable varname will be operated upon.

We will reconstruct the data in the zone from the zone-averaged variables using the same limited slopes as in the advection routine. Getting a good multidimensional reconstruction polynomial is hard – we want it to be bilinear and monotonic – we settle for having each slope be independently monotonic:

          (x)         (y)
f(x,y) = m    x/dx + m    y/dy + <f>

where the m’s are the limited differences in each direction. When averaged over the parent cell, this reproduces <f>.

Each zone’s reconstrution will be averaged over 4 children:

+-----------+     +-----+-----+
|           |     |     |     |
|           |     |  3  |  4  |
|    <f>    | --> +-----+-----+
|           |     |     |     |
|           |     |  1  |  2  |
+-----------+     +-----+-----+

We will fill each of the finer resolution zones by filling all the 1’s together, using a stride 2 into the fine array. Then the 2’s and …, this allows us to operate in a vector fashion. All operations will use the same slopes for their respective parents.

register_var(name, bc)[source]

Register a variable with CellCenterData2d object.

Parameters:
namestr

The variable name

bcBC object

The boundary conditions that describe the actions to take for this variable at the physical domain boundaries.

restrict(varname, N=2)[source]

Restrict the variable varname to a coarser grid (factor of 2 coarser) and return an array with the resulting data (and same number of ghostcells)

set_aux(keyword, value)[source]

Set any auxiliary (scalar) data. This data is simply carried along with the CellCenterData2d object

Parameters:
keywordstr

The name of the datum

valueany time

The value to associate with the keyword

write(filename)[source]

create an output file in HDF5 format and write out our data and grid.

write_data(f)[source]

write the data out to an hdf5 file – here, f is an h5py File pbject

zero(name)[source]

Zero out the data array associated with variable name.

Parameters:
namestr

The name of the variable to zero

class pyro.mesh.patch.FaceCenterData2d(grid, idir, dtype=<class 'numpy.float64'>)[source]

Bases: CellCenterData2d

A class to define face-centered data that lives on a grid. Data can be face-centered in x or y. This is built in the same multistep process as a CellCenterData2d object

add_derived(func)[source]

Register a function to compute derived variable

Parameters:
funcfunction

A function to call to derive the variable. This function should take two arguments, a CellCenterData2d object and a string variable name (or list of variables)

create()[source]

Called after all the variables are registered and allocates the storage for the state data. For face-centered data, we have one more zone in the face-centered direction.

fill_BC(name)[source]

Fill the boundary conditions. This operates on a single state variable at a time, to allow for maximum flexibility.

We do periodic, reflect-even, reflect-odd, and outflow

Each variable name has a corresponding BC stored in the CellCenterData2d object – we refer to this to figure out the action to take at each boundary.

Parameters:
namestr

The name of the variable for which to fill the BCs.

get_var_by_index(n)[source]

Return a data array for the variable with index n in the data array. Any changes made to this are automatically reflected in the CellCenterData2d object.

Parameters:
nint

The index of the variable to access

Returns:
outndarray

The array of data corresponding to the index

get_vars()[source]

Return the entire data array. Any changes made to this are automatically reflected in the CellCenterData2d object.

Returns:
outndarray

The array of data

prolong(varname)[source]

Prolong the data in the current (coarse) grid to a finer (factor of 2 finer) grid. Return an array with the resulting data (and same number of ghostcells). Only the data for the variable varname will be operated upon.

We will reconstruct the data in the zone from the zone-averaged variables using the same limited slopes as in the advection routine. Getting a good multidimensional reconstruction polynomial is hard – we want it to be bilinear and monotonic – we settle for having each slope be independently monotonic:

          (x)         (y)
f(x,y) = m    x/dx + m    y/dy + <f>

where the m’s are the limited differences in each direction. When averaged over the parent cell, this reproduces <f>.

Each zone’s reconstrution will be averaged over 4 children:

+-----------+     +-----+-----+
|           |     |     |     |
|           |     |  3  |  4  |
|    <f>    | --> +-----+-----+
|           |     |     |     |
|           |     |  1  |  2  |
+-----------+     +-----+-----+

We will fill each of the finer resolution zones by filling all the 1’s together, using a stride 2 into the fine array. Then the 2’s and …, this allows us to operate in a vector fashion. All operations will use the same slopes for their respective parents.

restrict(varname, N=2)[source]

Restrict the variable varname to a coarser grid (factor of 2 coarser) and return an array with the resulting data (and same number of ghostcells)

write_data(f)[source]

write the data out to an hdf5 file – here, f is an h5py File pbject

class pyro.mesh.patch.Grid2d(nx, ny, ng=1, xmin=0.0, xmax=1.0, ymin=0.0, ymax=1.0)[source]

Bases: object

the 2-d grid class. The grid object will contain the coordinate information (at various centerings).

A basic (1-d) representation of the layout is:

|     |      |     X     |     |      |     |     X     |      |     |
+--*--+- // -+--*--X--*--+--*--+- // -+--*--+--*--X--*--+- // -+--*--+
   0          ng-1    ng   ng+1         ... ng+nx-1 ng+nx      2ng+nx-1

                     ilo                      ihi

|<- ng guardcells->|<---- nx interior zones ----->|<- ng guardcells->|

The ‘*’ marks the data locations.

coarse_like(N)[source]

return a new grid object coarsened by a factor n, but with all the other properties the same

fine_like(N)[source]

return a new grid object finer by a factor n, but with all the other properties the same

scratch_array(nvar=1)[source]

return a standard numpy array dimensioned to have the size and number of ghostcells as the parent grid

pyro.mesh.patch.cell_center_data_clone(old)[source]

Create a new CellCenterData2d object that is a copy of an existing one

Parameters:
oldCellCenterData2d object

The CellCenterData2d object we wish to copy

pyro.mesh.patch.do_demo()[source]

show examples of the patch methods / classes

pyro.mesh.reconstruction module

Support for computing limited differences needed in reconstruction of slopes in constructing interface states.

pyro.mesh.reconstruction.flatten(myg, q, idir, ivars, rp)[source]

compute the 1-d flattening coefficients

pyro.mesh.reconstruction.flatten_multid(myg, q, xi_x, xi_y, ivars)[source]

compute the multidimensional flattening coefficient

pyro.mesh.reconstruction.limit(data, myg, idir, limiter)[source]

a single driver that calls the different limiters based on the value of the limiter input variable.

pyro.mesh.reconstruction.limit2(a, myg, idir)[source]

2nd order monotonized central difference limiter

pyro.mesh.reconstruction.limit4(a, myg, idir)[source]

4th order monotonized central difference limiter

pyro.mesh.reconstruction.nolimit(a, myg, idir)[source]

just a centered difference without any limiting

pyro.mesh.reconstruction.well_balance(q, myg, limiter, iv, grav)[source]

subtract off the hydrostatic pressure before limiting. Note, this only considers the y direction.

pyro.mesh.reconstruction.weno(q, order)[source]

Perform WENO reconstruction

Parameters:
qnp array

input data with 3 ghost zones

orderint

WENO order (k)

Returns:
q_plus, q_minusnp array

data reconstructed to the right / left respectively

pyro.mesh.reconstruction.weno_upwind(q, order)[source]

Perform upwinded (left biased) WENO reconstruction

Parameters:
qnp array

input data

orderint

WENO order (k)

Returns:
q_plusnp array

data reconstructed to the right