synthesizer.load_data.load_eagle

A module for interfacing with the outputs of the EAGLE hydrodynamical simulations suite.

The EAGLE hdf5 data loading scripts have been taken from the eagle_io package (https://github.com/flaresimulations/eagle_IO/). This has been to make the call as self-contained as possible.

The functions here used in https://github.com/flaresimulations/synthesizer-pipeline

Functions

synthesizer.load_data.load_eagle.apply_CGSUnits_conversion(filename, dataset, dat, verbose=True)[source]
Parameters:
  • filename (str) – filename to read from

  • dataset (str) – dataset to read attribute

  • dat (array) – dataset array to apply conversion

Return type:

ndarray[Any, dtype[Any]]

Returns:

Numpy array of the dataset converted to CGS units

synthesizer.load_data.load_eagle.apply_PhotUnits(filename, dataset, dat, verbose=True)[source]

” :type filename: str :param filename: filename to read from :type filename: str :type dataset: str :param dataset: dataset to read attribute :type dataset: str :type dat: ndarray[Any, dtype[Any]] :param dat: dataset array to apply conversion :type dat: array

Return type:

unyt_array[unyt_quantity]

Returns:

Unyt array of the dataset in the photometric units

synthesizer.load_data.load_eagle.apply_hfreeUnits_conversion(filename, dataset, dat, verbose=True)[source]
Parameters:
  • filename (str) – filename to read from

  • dataset (str) – dataset to read attribute

  • dat (array) – dataset array to apply conversion

Return type:

ndarray[Any, dtype[Any]]

Returns:

Numpy array of the dataset converted to h free units

synthesizer.load_data.load_eagle.apply_physicalUnits_conversion(filename, dataset, dat, verbose=True)[source]
Parameters:
  • filename (str) – filename to read from

  • dataset (str) – dataset to read attribute

  • dat (array) – dataset array to apply conversion

Return type:

ndarray[Any, dtype[Any]]

Returns:

Numpy array of the dataset converted to physical units

synthesizer.load_data.load_eagle.assign_galaxy_prop(ii, zed, boxl, aperture, grpno, sgrpno, cop, s_grpno, s_sgrpno, s_imasses, s_masses, s_ages, s_Zsmooth, s_coords, s_hsml, s_oxygen, s_hydrogen, g_grpno, g_sgrpno, g_masses, g_Zsmooth, g_sfr, g_coords, g_hsml, verbose, s_kwargs={}, g_kwargs={})[source]

Load stellar and gas particle data into synthesizer galaxy object.

Parameters:
  • ii (int) – galaxy number

  • zed (float) – redshift

  • boxl (float) – simulation box side length

  • aperture (float) – aperture to use from centre of potential

  • grpno (array) – Group numbers in this chunk

  • sgrpno (array) – Subgroup numbers in this chunk

  • cop (array) – centre of potential of subhalos in this chunk

  • s_grpno (array) – Stellar particle group numbers

  • s_sgrpno (array) – Stellar particle subgroup numbers

  • s_imasses (array) – Stellar particle initial masses in Msun

  • s_masses (array) – Stellar particle current masses in Msun

  • s_ages (array) – Stellar particle ages in Gyr

  • s_Zsmooth (array) – Stellar particle smoothed metallicity

  • s_coords (array) – Stellar particle coordinates in pMpc

  • s_hsml (array) – Stellar particle smoothing length in pMpc

  • s_oxygen (array) – Stellar particle abundance in oxygen

  • s_hydrogen (array) – Stellar particle abundance in hydrogen

  • g_grpno (array) – Gas particle group number

  • g_sgrpno (array) – Gas particle subgroup number

  • g_masses (array) – Gas particle masses in Msun

  • g_Zsmooth (array) – Gas particle smoothed metallicity

  • g_sfr (array) – Gas particle instantaneous SFR in Msun/yr

  • g_coords (ndarray[Any, dtype[float32]]) – (array) Gas particle coordinates in pMpc

  • g_hsml (array) – Gas particle smoothing length in pMpc

  • verbose (bool) – Are we talking?

  • s_kwargs (dictionary) – kwargs for stars

  • g_kwargs (dictionary) – kwargs for gas

Return type:

Galaxy

Returns:

Galaxy

synthesizer galaxy object

synthesizer.load_data.load_eagle.get_age(scale_factors, z, numThreads=4)[source]

Function to convert scale factor to z

Parameters:
  • scale_factors (array) – scale factor of the star particle

  • z (float) – redshift

  • numThreads (int) – number of threads to use

Return type:

ndarray[Any, dtype[Any]]

Returns:

array of ages of the star particles in years

synthesizer.load_data.load_eagle.get_files(fileType, directory, tag)[source]

Fetch filename for the different eagle outputs for reading

Parameters:
  • fileType (string) – type of file in the eagle outputs

  • directory (string) – eagle data file location

  • tag (string) – snapshot tag to load

Return type:

List[str]

Returns:

array of files to read

synthesizer.load_data.load_eagle.get_star_formation_time(scale_factor)[source]

Function to convert scale factor to z

Parameters:

scale_factor (float) – scale factor of the star particle

Return type:

float

Returns:

age of the star particle in years

synthesizer.load_data.load_eagle.load_EAGLE(fileloc, args, tot_chunks=1535, verbose=False)[source]

Load EAGLE required EAGLE galaxy properties for generating their SEDs Most useful for running on high-z snaps or individual chunks

Parameters:
  • fileloc (string) – eagle data file location

  • args (namedtuple) – parser arguments passed on to this job

  • tot_chunks (int) – total number of files to process

  • verbose (bool) – Are we talking?

Return type:

List[Union[Galaxy, Never]]

Returns:

a dictionary of Galaxy objects with stars and gas components

synthesizer.load_data.load_eagle.load_EAGLE_shm(chunk, fileloc, tag, s_len, g_len, args, numThreads=1, dtype='float32', tot_chunks=1535, verbose=False)[source]

Load EAGLE required EAGLE galaxy properties for generating their SEDs using numpy memmap. Most useful for running on high-z snaps

Parameters:
  • chunk (int) – file number to process

  • fileloc (string) – eagle data file location

  • tag (string) – snapshot tag to load

  • s_len (int) – total number of star particles

  • g_len (int) – total number of gas particles

  • args (namedtuple) – parser arguments passed on to this job

  • numThreads (int) – number of threads to use

  • dtype (numpy object) – data type of the array in memory

  • tot_chunks (int) – total number of files to process

  • verbose (bool) – Are we talking?

Return type:

List[Union[Galaxy, Never]]

Returns:

a dictionary of Galaxy objects with stars and gas components

synthesizer.load_data.load_eagle.np_shm_read(name, array_shape, dtype)[source]

Read the required numpy memmap array

Parameters:
  • name (str) – location of the memmap

  • array_shape (tuple) – shape of the memmap array

  • dtype (str) – data type of the memmap array

Return type:

ndarray[Any, dtype[Any]]

Returns:

Required eagle array from numpy memmap location

synthesizer.load_data.load_eagle.read_array(ftype, directory, tag, dataset, numThreads=1, noH=False, physicalUnits=False, CGS=False, photUnits=False, verbose=True)[source]
Parameters:
  • ftype (str) – eagle file type to read

  • directory (str) – location of the eagle simulation directory

  • tag (str) – snapshot tag to read

  • dataset (str) – name of the dataset to read

  • numThreads (int) – number threads to use

  • noH (bool) – remove any reduced Hubble factors

  • physicalUnits (bool) – return in physical units

  • CGS (bool) – return in CGS units

  • verbose (bool) – verbose condition

Return type:

Union[ndarray[Any, dtype[Any]], unyt_array[unyt_quantity]]

Returns:

Numpy/unyt array from eagle hdf5 filetype

synthesizer.load_data.load_eagle.read_hdf5(filename, dataset)[source]

Read the required dataset from the eagle hdf5 file

Parameters:
  • filename (str) – name of the hdf5 file

  • dataset (str) – name of the dataset to extract

Return type:

ndarray[Any, dtype[Any]]

Returns:

Numpy array of the required hdf5 dataset

synthesizer.load_data.load_eagle.read_header(ftype, directory, tag, dataset)[source]
ftype (str)

eagle file type to read

directory (str)

location of the eagle simulation directory

tag (str)

snapshot tag to read

dataset (str)

name of the dataset to read

Return type:

ndarray[Any, dtype[Any]]

Returns:

Reads in required eagle hdf5 header data