output

Helper functions for reading and modifying flash output files. For more infromation on using this module to make figures with matplotlib and yt, please refer to the following presentation: http://flash.uchicago.edu/site/flashcode/user_support/yt_flash_informal.pdf

flash.output.slice(axis, coord, field, pf, bounds=None, resolution=600, method='nearest', **kwargs)[source]

Grabs the slice of a certain field (parameter, data) from the flash output file pf. This slice is performed at the coord along the axis.

This function requires both scipy and yt and returns x, y, and z data that is suitable for plotting with matplotlib.imshow().

Parameters :

axis : int

The axis along which to slice. Can be 0, 1, or 2 for x, y, z.

coord : float

The coordinate along the axis at which to slice. This is in “domain” coordinates.

field : str

A field to retrieve, e.g. ‘dens’.

pf : str or yt.data_objects.static_output

The checkpoint or plot file, may be either a path to the file or the yt object returned by pf.mods.load().

bounds : len-4 sequence of floats or None, optional

This defines the area within the domain that should be sliced. By default, the whole domain is taken. Applies the bounds in axes order, ie if axis=2 then bounds=[xmin, xmax, ymin, ymax].

resolution : int or len-2 sequence of ints, optional

If a list of intergers, this is the number of points (pixels) along each dimension. If a single integer, then this is the resolution the longest axis. The resolution on the other axis is then calculated to preserve the aspect ratio.

method : str, optional

Interpolation method flag, passed directly down to the function scipy.interpolate.griddata().

kwargs : optional

All other keyword arguments are passes to the slice() method on the pf.h object from yt.

Returns :

xdat : 2d numpy array

Meshed data for x.

ydat : 2d numpy array

Meshed data for y.

zdat : 2d numpy array

Meshed data for z.

flash.output.slice_gradient(axis, coord, field, pf, bounds=None, resolution=600, method='nearest', **kwargs)[source]

Grabs the gradient of a slice of a certain field (parameter, data) from the flash output file pf. This slice is performed at the coord along the axis.

This function requires both scipy and yt and returns x, y, and z data that is suitable for plotting with matplotlib.imshow().

Parameters :

axis : int

The axis along which to slice. Can be 0, 1, or 2 for x, y, z.

coord : float

The coordinate along the axis at which to slice. This is in “domain” coordinates.

field : str

A field to retrieve, e.g. ‘dens’.

pf : str or yt.data_objects.static_output

The checkpoint or plot file, may be either a path to the file or the yt object returned by pf.mods.load().

bounds : len-4 sequence of floats or None, optional

This defines the area within the domain that should be sliced. By default, the whole domain is taken. Applies the bounds in axes order, ie if axis=2 then bounds=[xmin, xmax, ymin, ymax].

resolution : int or len-2 sequence of ints, optional

If a list of intergers, this is the number of points (pixels) along each dimension. If a single integer, then this is the resolution the longest axis. The resolution on the other axis is then calculated to preserve the aspect ratio.

method : str, optional

Interpolation method flag, passed directly down to the function scipy.interpolate.griddata().

kwargs : optional

All other keyword arguments are passes to the slice() method on the pf.h object from yt.

Returns :

adat : 2d numpy array

Meshed data for the first non-axis axis, eg x when axis=z.

bdat : 2d numpy array

Meshed data for the second non-axis axis, eg y when axis=z.

dfdadat : 2d numpy array

Gradient of field along a.

dfdbdat : 2d numpy array

Gradient of field along b.

magdat : 2d numpy array

Magnitude of the gradient of field.

flash.output.lineout(p1, p2, field, pf, **kwargs)[source]

Grabs a line out (ray) of a certain field (parameter, data) from the flash output file pf.

This function requires both scipy and yt and returns x, y, z and value data that is suitable for plotting with matplotlib.plot().

Parameters :

p1 : three-tuple of floats

The first point in the line-out.

p2 : three-tuple of floats

The second point in the line-out.

field : str

A field to retrieve, e.g. ‘dens’.

pf : str or yt.data_objects.static_output

The checkpoint or plot file, may be either a path to the file or the yt object returned by pf.mods.load().

kwargs : optional

All other keyword arguments are passes to the ray() method on the pf.h object from yt.

Returns :

x : 1D numpy array

Interpolated x data.

y : 1D numpy array

Interpolated y data.

z : 1D numpy array

Interpolated z data.

v : 1D numpy array

Interpolated field values along ray.

flash.output.shock_on_lineout(p1, p2, field, pf, threshold=1e-06, min_threshold=1e-36, **kwargs)[source]

Finds the shock of a certain field (parameter, data) along a line out (ray) in a flash output file pf. Currently this assumes that p1 is closer to the center of the shock than p2. Moreover the geometry from the FLASH simulation must be cartesian or 2D cylindrical (rz).

Parameters :

p1 : three-tuple of floats

The first point in the line-out.

p2 : three-tuple of floats

The second point in the line-out.

field : str

A field to retrieve, e.g. ‘dens’.

pf : str or yt.data_objects.static_output

The checkpoint or plot file, may be either a path to the file or the yt object returned by pf.mods.load().

threshold : float, optional

The value above which gradients must exceed to be considered a shock. Required for ignoring low-level noise.

min_threshold : float, optional

If a shock is not found at a given threshold level, the threshold value is reduced by an order of magnitude until a shock is found. This continues until a minimum threshold is reached.

kwargs : optional

All other keyword arguments are passes to the ray() method on the pf.h object from yt.

Returns :

shock_p : three-tuple of floats

Shock position.

shock_v : field type (float)

Shock peak value.

See also

flash.output.lineout
Used to get the lineout from two points and a data file.
flash.analysis.shock_detect
Used to find the shock itself.
flash.output.load_laser_dat(filename)[source]

Reads in a LaserEnergyProfile.dat file as a numpy structured array.

Parameters :

filename : str or file-like

The string path to the file or an object which implements the file interface

Returns :

dat : ndarray

This structured array has the following dtype and units:

Field

dtype

units

Step

int32

unitless

Time

float64

s

dt

float64

s

Energy In

float64

erg

Energy Out

float64

erg

dE In

float64

erg

dE Out

float64

erg

Previous topic

analysis

Next topic

domain specific languages

This Page