pynbody.analysis.ionfrac.IonFractionTable#

class pynbody.analysis.ionfrac.IonFractionTable(redshift_values, log_temp_values, log_den_values, tables)[source]#

Bases: IonFractionTableBase

Class for calculating ion fractions from a grid of cloudy models.

Rather than use this class directly, for many uses it is simpler to use calculate() which will use the currently loaded or default table.

Methods

calculate(simulation[, ion])

Calculate the ionisation fraction for the gas particles of a given simulation

from_cloudy(cloudy_path[, table, ...])

Generate a table by running cloudy with the specified ionising radiation table and parameters.

load(filename)

Load a table from a numpy .npz file, generated using save()

plot([ion, redshift])

Use matplotlib to plot the ion fraction table for a given ion at a given redshift

save(filename)

Save the table to a numpy .npz file

__init__(redshift_values, log_temp_values, log_den_values, tables)[source]#

Initialise an ion fraction table from raw data.

Most users will instead want to use load() to load a pre-existing table, or from_cloudy() to generate a new table by executing cloudy.

Parameters:
  • redshift_values (array_like) – Redshift values at which the tables are defined

  • log_temp_values (array_like) – Log10 temperature values at which the tables are defined

  • log_den_values (array_like) – Log10 density values at which the tables are defined

  • tables (dict) – Dictionary of tables, with keys being ion names and values being 3D numpy arrays of ion fraction values. The shape of each array should be (len(redshift_values), len(log_temp_values), len(log_den_values)).

calculate(simulation, ion='ovi')[source]#

Calculate the ionisation fraction for the gas particles of a given simulation

Values are interpolated from the (rho,T) table at the appropriate redshift, and any values outside the range of the table are clamped to the nearest valid value.

Parameters:
  • simulation (pynbody.snapshot.SimSnap) – The simulation snapshot to calculate the ion fractions for. The gas particles must have ‘rho’ and ‘temp’ fields.

  • ion (str) – The name of the ion to calculate the fraction for, e.g. HI, MgII, OVI etc. The only molecular fraction available is H2. The default is ‘ovi’. Case insensitive.

Returns:

The ion fraction for each gas particle in the simulation, according to the table.

Return type:

array-like

classmethod from_cloudy(cloudy_path, table='hm12', redshift_range=(0, 15), num_redshifts=10, log_temp_range=(2.0, 8.0), num_temps=10, log_den_range=(-8.0, 2.0), num_dens=10)[source]#

Generate a table by running cloudy with the specified ionising radiation table and parameters.

This can take a long time, but the resulting table can then be saved using the save() method and reused by calling load(). The grid is computed in parallel using the default number of processors detected by multiprocessing.Pool. A progress bar is displayed using tqdm.

Parameters:
  • cloudy_path (str) – Path to the cloudy executable

  • table (str) – Name of the cloudy radiation table to use. The default is ‘hm12’.

  • redshift_range (tuple) – Minimum and maximum redshift values to use

  • num_redshifts (int) – Number of redshift values to use. These are spaced equally in log(1+z).

  • log_temp_range (tuple) – Minimum and maximum log10 temperature values to use

  • num_temps (int) – Number of temperature values to use, spaced equally in log space

  • log_den_range (tuple) – Minimum and maximum log10 density values to use

  • num_dens (int) – Number of density values to use, spaced equally in log space

classmethod load(filename)[source]#

Load a table from a numpy .npz file, generated using save()

If the file is not found, it is assumed to be a pynbody-provided table. If such a built-in table exists, the path is modified automatically to point at it. If it does not exist, an attempt is made to download it from a zenodo repository.

plot(ion='ovi', redshift=0.0)[source]#

Use matplotlib to plot the ion fraction table for a given ion at a given redshift

save(filename)[source]#

Save the table to a numpy .npz file