API Documentation

class lightcone_io.ShellArray(basedir, basename, remote_dir=None)

Bases: Sequence, LocalOrRemoteFile

Sequence-like container for a set of lightcone shells. This class is the recommended way to read lightcone HEALPix maps.

A ShellArray can be indexed with an integer to return a Shell.

Parameters:
  • basedir (str) – directory containing the lightcone outputs

  • basename (str) – name of the subdirectory for this lightcone (e.g. lightcone0)

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

property have_shell_dir

Return True if the map files are in shell_X subdirectories.

Returns:

True if the shell_0 directory exists, False otherwise

Return type:

bool

property nr_shells

Return the number of shells in this lightcone output by counting files.

Returns:

the number of shells.

Return type:

int

class lightcone_io.Shell(basedir, basename, shell_nr, remote_dir=None, have_shell_dir=None)

Bases: Mapping, LocalOrRemoteFile

Dict-like container for all of the HealpixMap instances associated with a lightcone shell. Subscripting a Shell with a map name returns a HealpixMap.

Parameters:
  • basedir (str) – directory containing the lightcone outputs

  • basename (str) – name of the subdirectory for this lightcone (e.g. lightcone0)

  • shell_nr (int) – index of the lightcone shell to read

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

  • have_shell_dir (bool or None) – True if each shell is in a shell_X subdirectory, False if not, None if not known

static filename(basedir, basename, shell_nr, file_nr, have_shell_dir)

Return the filename for one of the shell files.

Parameters:
  • basedir (str) – directory containing the lightcone outputs

  • basename (str) – name of the subdirectory for this lightcone (e.g. lightcone0)

  • shell_nr (int) – index of the lightcone shell to read

  • file_nr (int) – index of the lightcone shell file to read

  • have_shell_dir (bool) – True if each shell is in a shell_X subdirectory, False if not

have_shell_dir(basedir, basename, shell_nr)
class lightcone_io.HealpixMap(filenames, map_name, remote_dir=None)

Bases: Sequence, LocalOrRemoteFile

Class used to read lightcone HEALPix maps written out by SWIFT. Numpy-style indexing can be used to read pixels from the map.

Parameters:
  • filenames (list of str) – names of the HDF5 files containing map data

  • map_name (str) – name of the map to read

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

property nside

Return the nside (resolution) parameter of this HEALPix map

property dtype

Return the numpy data type of this HEALPix map

property units

Return the unyt units of this HEALPix map

read_pixels(start, stop)

Read the specified range of pixels from the map and return a numpy array.

Parameters:
  • start (int) – index of the first pixel to read

  • stop (int) – index just after the last pixel to read

Returns:

a numpy array with the specified pixel data

Return type:

numpy.ndarray

class lightcone_io.ParticleLightcone(fname, comm=None, extra_filename=None, remote_dir=None)

Bases: Mapping, LocalOrRemoteFile

Class used to read particle lightcones. This is a dict-like container for the IndexedLightconeParticleType instances which represent the different particle types in this lightcone (gas, DM, stars etc).

Particle lightcones are opened by specifying the name of any one of the particle data files.

If an MPI communicator is specified, then reads are parallelized over MPI ranks in the communicator.

Parameters:
  • filename (str) – name of one of the lightcone particle files

  • comm (mpi4py.MPI.Comm, or None) – MPI communicator

  • extra_filename (str, or None) – name of a file with extra particle properties

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

class lightcone_io.IndexedLightconeParticleType(type_name, metadata, index, units, filenames, extra_filenames=None, remote_dir=None)

Bases: LocalOrRemoteFile

Class to read a single particle type from a lightcone. This assumes that the lightcone particles have been sorted by redshift and healpix pixel index. We can then efficiently extract the particles in particular redshift ranges and regions of the sky.

Parameters:
  • type_name (str) – HDF5 group name for this particle type

  • metadata (dict) – lightcone metadata read by ParticleLightcone

  • units (dict) – system of units used in this output

  • filenames (list of str) – names of the HDF5 files containing particle data

  • extra_filenames (list of str or None) – names of the HDF5 files containing additional data

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

set_mpi_mode(comm)

Enable MPI mode for this output. In MPI mode each MPI rank reads a subset of the selected particles.

Parameters:

comm (mpi4py.MPI.Comm) – MPI communicator

get_redshift_bins_in_range(redshift_min, redshift_max)

Return indexes of all redshift bins overlapping specified range

Parameters:
  • redshift_min (float) – minimum redshift

  • redshift_max (float) – maximum redshift

Returns:

a numpy array of ints with the bin indexes

Return type:

numpy.ndarray

get_pixels_in_radius(vector, radius)

Return indexes of all healpix pixels within radius of vector

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

Returns:

a numpy array of ints with the pixel indexes

Return type:

numpy.ndarray

get_cell_indexes_from_bins(redshift_bins, healpix_bins)

Given arrays of redshift and healpix bins to read, return indexes of all cells which must be read.

Parameters:
  • redshift_bins (numpy.ndarray) – redshift bin indexes as array of ints

  • healpix_bins (numpy.ndarray) – healpix pixel indexes as array of ints

Returns:

a numpy array of ints with the cell indexes

Return type:

numpy.ndarray

read_cells(property_names, cells_to_read)

Given an array of cells to read, read the specified properties from the lightcone particle files.

Parameters:
  • property_names (list of str) – list of properties to read

  • cells_to_read (numpy.ndarray) – a numpy array of ints with the cell indexes

Returns:

dict of ndarrays with the requested properties

Return type:

dict

get_cell_indexes_from_vector_radius_redshift(vector, radius, redshift_range)

Given a line of sight vector, and angular radius and a redshift range, compute which cells we need to read.

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • redshift_range (sequence of two floats [z_min, z_max]) – redshift range to read

Returns:

a numpy array of ints with the cell indexes

Return type:

numpy.ndarray

count_particles(vector=None, radius=None, redshift_range=None)

Return the number of particles in a region defined by a line of sight vector, an angular radius and a redshift range.

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • redshift_range (sequence of two floats [z_min, z_max]) – redshift range to read

Returns:

the number of particles selected

Return type:

int

iterate_chunks(property_names, vector=None, radius=None, redshift_range=None, max_particles=1048576)

Iterate over chunks of particles in a region defined by a line of sight vector, an angular radius and a redshift range. This can be used to avoid reading all of the particles into memory at once. Yields dicts of arrays where the keys are the requested particle property names and the values are the particle property values.

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • redshift_range (sequence of two floats [z_min, z_max]) – redshift range to read

  • max_particles (int) – maximum number of particles to return in one chunk

Yield:

dict of ndarray with the particle properties for each chunk

Return type:

dict

read(property_names, vector=None, radius=None, redshift_range=None)

Read all particles in a region defined by a line of sight vector, an angular radius and a redshift range. Returns a dict of arrays where the keys are the requested particle property names and the values are the particle property values. Always returns complete cells, so it may return particles outside the requested region.

Parameters:
  • property_names (list of str) – list of properties to read

  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • redshift_range (sequence of two floats [z_min, z_max]) – redshift range to read

Returns:

dict of ndarray with the particle properties

Return type:

dict

read_exact(property_names, vector=None, radius=None, redshift_range=None)

Read all particles in a region defined by a line of sight vector, an angular radius and a redshift range. Returns a dict of arrays where the keys are the requested particle property names and the values are the particle property values. This version filters out any particles outside the requested region.

Parameters:
  • property_names (list of str) – list of properties to read

  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • redshift_range (sequence of two floats [z_min, z_max]) – redshift range to read

Returns:

dict of ndarray with the particle properties

Return type:

dict

class lightcone_io.HaloLightconeFile(filename, soap_filename=None, remote_dir=None)

Bases: LocalOrRemoteFile

Class used to read a single file from the halo lightcone. Each file contains halos which were identified in a single snapshot of the simulation. Halos are sorted according to their pixel index in a low resolution HEALPix map so that regions on the sky can be extracted efficiently.

If the name of the SOAP file which was used to generate the lightcone halo catalogue is supplied, then it will be possible to look up extra halo properties from SOAP.

Parameters:
  • filename (str) – Name of the halo lightcone file to open

  • soap_filename (str) – Name of the SOAP output with extra halo properties

  • remote_dir (hdfstream.RemoteDirectory, or None for local files) – remote directory containing the file, or None

property properties

Return a list of all halo property names in this file. These are the values which are valid to pass to the properties parameter of read_halos_in_pixels() and read_halos_in_radius().

Note that this does not include any additional properties which we may be able to read from SOAP.

get_pixels_in_radius(vector, radius)

Return indexes of all healpix pixels within radius of vector

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

Returns:

a numpy array of ints with the pixel indexes

Return type:

numpy.ndarray

read_halos_in_pixels(pixels, properties)

Read halos in the specified HEALPix pixels and return the requested halo properties in a dict of unyt arrays.

If a SOAP filename was specifed then we’ll check the SOAP file for any halo properties which are not present in the halo lightcone file.

Parameters:
  • pixels (numpy.ndarray, or None) – array of HEALPix pixel indexes to read, or None to read all

  • properties (list of str) – list of halo properties (i.e. HDF5 dataset names) to read

Returns:

dict of arrays with the halo properties

Return type:

dict of unyt.unyt_array

read_halos_in_radius(vector, radius, properties)

Read halos in an angular radius around the specified line of sight vector and return the requested halo properties in a dict of unyt arrays. May also return some halos slightly outside the radius because if we read a pixel we read all of the halos in it.

Parameters:
  • vector (numpy.ndarray) – direction vector as an array of 3 floats

  • radius (float) – angular radius in radians

  • properties (list of str) – list of halo properties (i.e. HDF5 dataset names) to read

Returns:

dict of arrays with the halo properties

Return type:

dict of unyt.unyt_array

read_halos(properties)

Read all halos in this file and return the requested halo properties in a dict of unyt arrays.

Parameters:

properties (list of str) – list of halo properties (i.e. HDF5 dataset names) to read

Returns:

dict of arrays with the halo properties

Return type:

dict of unyt.unyt_array

class lightcone_io.XrayCalculator_LC(redshifts, table, bands, observing_types)

Bases: object

check_for_restframe(observing_types, bands)

remove the restframe X-ray observation types and corresponding bands

load_all_tables(redshifts, table, bands, observing_types)

Load the x-ray tables for the specified bands and observing-types

Params:

redshifts: of the particles table: either a path to the table to be read or a dictionary containing the tables themselves bands, observing_types: the bands and observation types, within the band to add to tables.

static get_index_1d(bins, subdata)

Find the closest bin index below the specified value, and the relative offset compared to that bin.

static get_index_1d_irregular(bins, subdata)

Find the closest bin index below the specified value, and the relative offset compared to that bin. Unlike get_index_1d, this allows for irregular bin spacing

static get_table_interp(idx_z, idx_he, idx_T, idx_n, t_z, d_z, t_T, d_T, t_nH, d_nH, t_He, d_He, X_Ray, abundance_to_solar)

4D interpolate the x-ray table for each traced metal Scale the metals with their respective relative solar abundances Add the metals to the background case without metals

find_indices(densities, temperatures, element_mass_fractions, masses, redshifts, fill_value=0)

Check interpolation table bounds Compute all interpolation bin indices, and the offsets from those bins Compute all the indices for the flattened x-ray table

interpolate_X_Ray(idx_z, idx_he, idx_T, idx_n, t_z, d_z, t_T, d_T, t_nH, d_nH, t_He, d_He, abundance_to_solar, joint_mask, volumes, data_n, bands=None, observing_types=None, fill_value=None)

Main function Calculate the particle emissivities through interpolation Convert to luminosity using the particle volume

class lightcone_io.Snapshot_Cosmology_For_Lightcone(snapshot_root_dir, search_snapshot_pattern=None)

Bases: object

get_snapshot_filename(rootdir, re_search_snapshot_pattern=None)
Relies on FLAMINGO snapshot virtual files to be named as:

“/snapshots/flamingo_4DigitInt/flamingo_4DigitInt.hdf5”

Params:
rootdir: path to directory with flamingo snapshots

“./path/to/flamingo/runs/LBoxsizeNresolution/SimName/snapshots”

re_search_snapshot_pattern: file naming pattern in regex terms.

Returns:

path to lowest redshift snapshot file.

cosmology_from_hdf5(snapshot_path)
get_internal_units(snapshot_path)
z2Myr(z)

Return a unyt_array [Myrs] with age of the universe at the input redshift(s).

redshift_with_time_offset(z, dt)
Params:

z: redshift of observer dt: length of time [Myr]

Returns:

The redshift equal to the input redshift - length of time (dt)

z2r(z)

Returns the comoving radius for a given redshift in cMpc

shell_comoving_raddii(redshift_filename)
Params:

redshift_filename: path to .txt file with shell redshift ranges

Returns:

N x 2 array of comoving radii, per shell returns [inner, outer]

get_AGN_delta_T_K(snapshot_path)

Return the Delta_T_AGN param required to identify recently heated particles

lightcone_volume(comoving_inner_radius, comoving_outer_radius, use_redshift=False, phi=None)

Compute the comoving volume of the lightcone given some inner and outer radius

Params:

comoving_inner_radius, comoving_outer_radius: the radius of the lightcone in Mpc use_redshift [boolean]: if False, then radii are in comoving distances. if True, then

input radii are redshift values.

phi: 1/2 angle of the cones apature [radian]

class lightcone_io.Xray_Filter(particle_data, numb_particles, remove_SF=True, log_T_min_K=5, remove_rhp=True, cosmo=None)

Bases: object

filter_by_StarFormationRates(sfr)

remove all star forming gas particles & update mask

Params

sfr: gas particle star formation rates

filter_by_Temperatures(temperatures, log_T_min_K)

keep only `hot’ gas particles & update mask

Params

temperatures: temperatures of gas particles log_T_min_K: minimum temperature to be considered a hot gas particle

time_from_last_AGN_feedback(scalefactors, last_agn_feedback_scalefactors, f_z2t, handle_negatives=0.0)

Compute the difference in time from current scalefactor and the scalefactor when last heated by AGN feedback. Due to compression on LastAGNFeedbackScaleFactors values, sometimes LastAGNFeedbackScaleFactors > scalefactor, therefore we must allow negative values.

Params

scalefactors: scalefactor of gas particles last_agn_feedback_scalefactor: scalefactor that particle was last

directly heated by AGN feedback and -1 if it has never been heated by AGN feedback

f_z2t: function to give age of universe for a given redshift handle_negatives: time difference given to where

last_agn_feedback_scalefactor > scalefactors

identify_recently_heated(cosmo, scalefactors, last_agn_feedback_scalefactors, temperatures, densities, RHP_filter_max_time_Myr, RHP_filter_log_density_cm3=-2.25)

Identify recently heated particles. Due to the BFloat16 lossy compression filter applied to LastAGNFeedbackScaleFactors we cannot reliably determine if a particle was last heated by AGN feedback within the last 15 Myr. Note the minimum density criterion is included to compensate for the uncertanty about LastAGNFeedbackScaleFactors values due to the lossy compression filter used for the particle lightcones. To reproduce snapshot selection of recently heated particles:

RHP_filter_max_time_Myr=15 RHP_filter_log_density=None

Params

cosmo: cosmology object, i.e. Snapshot_Cosmology_For_Lightcone() scalefactors: scalefactor of gas particles last_agn_feedback_scalefactor: scalefactor that particle was last

directly heated by AGN feedback and -1 if it has never been heated by AGN feedback

temperatures (unyt_array): temperatures of gas particles [K] densities (unyt_array): densities of gas particles [g/cm^3] RHP_filter_max_time_Myr: maximum time in Myrs that has elapsed from the

LastAGNFeedbackScaleFactors to be considered a recently heated particle. If == None then adopt 70 Myr as compromise for uncertanty.

RHP_filter_log_density: log10(minimum gas density [1/cm^3]), the minimum density

to be considered a recently heated particle. If == None, then do not use density to select recently heated particles. -2.25 selected as default for McDonald+26 based on HYDRO_FIDUCIAL L1_m9 snapshots & ROSAT photon intrinsic maps

Returns

Boolean mask, True if RHP, otherwise False