synthesizer.particle.stars

A module for working with arrays of stellar particles.

Contains the Stars class for use with particle based systems. This houses all the data detailing collections of stellar particles. Each property is stored in (N_star, ) shaped arrays for efficiency.

We also provide functions for creating “fake” stellar distributions by sampling a SFZH.

In both cases a myriad of extra optional properties can be set by providing them as keyword arguments.

Example usages:

stars = Stars(initial_masses, ages, metallicities,

redshift=redshift, current_masses=current_masses, …)

stars = sample_sfhz(sfzh, n, total_initial_mass,

smoothing_lengths=smoothing_lengths, tau_v=tau_vs, coordinates=coordinates, …)

Functions

synthesizer.particle.stars.sample_sfhz(sfzh, log10ages, log10metallicities, nstar, initial_mass=1, **kwargs)[source]

Create “fake” stellar particles by sampling a SFZH.

Parameters:
  • sfhz (array-like, float) – The Star Formation Metallicity History grid (from parametric.Stars).

  • log10ages (array-like, float) – The log of the SFZH age axis.

  • log10metallicities (array-like, float) – The log of the SFZH metallicities axis.

  • nstar (int) – The number of stellar particles to produce.

  • intial_mass (int) – The intial mass of the fake stellar particles.

Returns:

stars (Stars)

An instance of Stars containing the fake stellar particles.

Examples using synthesizer.particle.stars.sample_sfhz

Create sampled SED

Create sampled SED

Compare parametric and particle SEDs

Compare parametric and particle SEDs

Compare SPS grid assignment methods

Compare SPS grid assignment methods

Plot line of sight diagnostics

Plot line of sight diagnostics

Plot line of sight diagnostics

Plot line of sight diagnostics

Classes

class synthesizer.particle.stars.Stars(initial_masses, ages, metallicities, redshift=None, tau_v=None, alpha_enhancement=None, coordinates=None, velocities=None, current_masses=None, smoothing_lengths=None, s_oxygen=None, s_hydrogen=None, softening_length=None, centre=None, metallicity_floor=1e-05)[source]

The base Stars class. This contains all data a collection of stars could contain. It inherits from the base Particles class holding attributes and methods common to all particle types.

The Stars class can be handed to methods elsewhere to pass information about the stars needed in other computations. For example a Galaxy object can be passed a stars object for use with any of the Galaxy helper methods.

Note that due to the many possible operations, this class has a large number of optional attributes which are set to None if not provided.

initial_masses

The intial stellar mass of each particle in Msun.

Type:

array-like, float

ages

The age of each stellar particle in yrs.

Type:

array-like, float

metallicities

The metallicity of each stellar particle.

Type:

array-like, float

tau_v

V-band dust optical depth of each stellar particle.

Type:

array-like, float

alpha_enhancement

The alpha enhancement [alpha/Fe] of each stellar particle.

Type:

array-like, float

log10ages

Convenience attribute containing log10(age in yr).

Type:

array-like, float

log10metallicities

Convenience attribute containing log10(metallicity).

Type:

array-like, float

resampled

Flag for whether the young particles have been resampled.

Type:

bool

current_masses

The current mass of each stellar particle in Msun.

Type:

array-like, float

smoothing_lengths

The smoothing lengths (describing the sph kernel) of each stellar particle in simulation length units.

Type:

array-like, float

s_oxygen

fractional oxygen abundance.

Type:

array-like, float

s_hydrogen

fractional hydrogen abundance.

Type:

array-like, float

imf_hmass_slope

The slope of high mass end of the initial mass function (WIP).

Type:

float

nstars

The number of stellar particles in the object.

Type:

int

generate_line(grid, line_id, fesc, mask=None, method='cic', nthreads=0)[source]

Calculate rest frame line luminosity and continuum from an SPS Grid.

This is a flexible base method which extracts the rest frame line luminosity of this stellar population from the SPS grid based on the passed arguments.

Parameters:
  • grid (Grid) – A Grid object.

  • line_id (list/str) – A list of line_ids or a str denoting a single line. Doublets can be specified as a nested list or using a comma (e.g. ‘OIII4363,OIII4959’).

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • mask (array) – A mask to apply to the particles (only applicable to particle)

  • method (str) – The method to use for the interpolation. Options are: ‘cic’ - Cloud in cell ‘ngp’ - Nearest grid point

  • nthreads (int) – The number of threads to use in the C extension. If -1 then all available threads are used.

Returns:

Line

An instance of Line contain this lines wavelenth, luminosity, and continuum.

generate_lnu(grid, spectra_name, fesc=0.0, young=None, old=None, verbose=False, do_grid_check=False, grid_assignment_method='cic', parametric_young_stars=None, parametric_sfh='constant', aperture=None, nthreads=0)[source]

Generate the integrated rest frame spectra for a given grid key spectra for all stars. Can optionally apply masks.

Parameters:
  • grid (Grid) – The spectral grid object.

  • spectra_name (string) – The name of the target spectra inside the grid file.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • young (bool/float) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (bool/float) – If not None, specifies age in Myr at which to filter for old star particles.

  • verbose (bool) – Flag for verbose output.

  • do_grid_check (bool) –

    Whether to check how many particles lie outside the grid. This is True by default and provides a vital sanity check. There are instances when you may want to turn this off: - You know particles will lie outside the grid and want

    this behaviour. In this case the check is redundant.

    • You know your particle lie within the grid but don’t want to waste compute checking. This case is useful when working with large particle counts.

  • grid_assignment_method (string) – The type of method used to assign particles to a SPS grid point. Allowed methods are cic (cloud in cell) or nearest grid point (ngp) or there uppercase equivalents (CIC, NGP). Defaults to cic.

  • parametric_young_stars (bool/float) – If not None, specifies age in Myr below which we replace individual star particles with a parametric SFH.

  • parametric_sfh (string) – Form of the parametric SFH to use for young stars. Currently two are supported, Constant and TruncatedExponential, selected using the keyword arguments constant and exponential.

  • aperture (float) – If not None, specifies the radius of a spherical aperture to apply to the particles.

  • nthreads (int) – The number of threads to use in the C extension. If -1 then all available threads are used.

Returns:

Numpy array of integrated spectra in units of (erg / s / Hz).

generate_particle_line(grid, line_id, fesc, mask=None, method='cic', nthreads=0)[source]

Calculate rest frame line luminosity and continuum from an SPS Grid.

This is a flexible base method which extracts the rest frame line luminosity of this stellar population from the SPS grid based on the passed arguments and calculate the luminosity and continuum for each individual particle.

Parameters:
  • grid (Grid) – A Grid object.

  • line_id (list/str) – A list of line_ids or a str denoting a single line. Doublets can be specified as a nested list or using a comma (e.g. ‘OIII4363,OIII4959’).

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • mask (array) – A mask to apply to the particles (only applicable to particle)

  • method (str) – The method to use for the interpolation. Options are: ‘cic’ - Cloud in cell ‘ngp’ - Nearest grid point

  • nthreads (int) – The number of threads to use in the C extension. If -1 then all available threads are used.

Returns:

Line

An instance of Line contain this lines wavelenth, luminosity, and continuum.

generate_particle_lnu(grid, spectra_name, fesc=0.0, young=None, old=None, verbose=False, do_grid_check=False, grid_assignment_method='cic', nthreads=0)[source]

Generate the particle rest frame spectra for a given grid key spectra for all stars. Can optionally apply masks.

Parameters:
  • grid (Grid) – The spectral grid object.

  • spectra_name (string) – The name of the target spectra inside the grid file.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • young (bool/float) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (bool/float) – If not None, specifies age in Myr at which to filter for old star particles.

  • verbose (bool) – Flag for verbose output. By default False.

  • do_grid_check (bool) –

    Whether to check how many particles lie outside the grid. This is True by default and provides a vital sanity check. There are instances when you may want to turn this off:

    • You know particles will lie outside the grid and want this behaviour. In this case the check is redundant.

    • You know your particle lie within the grid but don’t want to waste compute checking. This case is useful when working with large particle counts.

  • grid_assignment_method (string) – The type of method used to assign particles to a SPS grid point. Allowed methods are cic (cloud in cell) or nearest grid point (ngp) or there uppercase equivalents (CIC, NGP). Defaults to cic.

  • nthreads (int) – The number of threads to use in the C extension. If -1 then all available threads are used.

Returns:

Numpy array of integrated spectra in units of (erg / s / Hz).

get_particle_line_attenuated(grid, line_ids, fesc=0.0, tau_v_nebular=None, tau_v_stellar=None, dust_curve_nebular=<synthesizer.dust.attenuation.PowerLaw object>, dust_curve_stellar=<synthesizer.dust.attenuation.PowerLaw object>, mask=None, method='cic', label='')[source]

Get a LineCollection containing attenuated lines for each particle.

Calculates attenuated properties (luminosity, continuum, EW) for a set of lines. Allows the nebular and stellar attenuation to be set separately.

Parameters:
  • grid (Grid) – The Grid object.

  • line_ids (list/str) – A list of line_ids or a str denoting a single line. Doublets can be specified as a nested list or using a comma (e.g. ‘OIII4363,OIII4959’).

  • fesc (float) – The Lyman continuum escaped fraction, the fraction of ionising photons that entirely escaped.

  • tau_v_BS (float) – V-band optical depth of the nebular emission.

  • tau_v_stellar (float) – V-band optical depth of the stellar emission.

  • dust_curve_nebular (dust_curve) – A dust_curve object specifying the dust curve. for the nebular emission

  • dust_curve_stellar (dust_curve) – A dust_curve object specifying the dust curve for the stellar emission.

  • mask (array) – A mask to apply to the particles (only applicable to particle)

  • method (str) – The method to use for the interpolation. Options are: ‘cic’ - Cloud in cell ‘ngp’ - Nearest grid point

Returns:

LineCollection

A dictionary like object containing line objects.

get_particle_line_intrinsic(grid, line_ids, fesc=0.0, mask=None, method='cic', label='')[source]

Get a LineCollection containing intrinsic lines for each particle.

The resulting LineCollection contains the intrinsic lines for each individual particle.

Parameters:
  • grid (Grid) – A Grid object.

  • line_ids (list/str) – A list of line_ids or a str denoting a single line. Doublets can be specified as a nested list or using a comma (e.g. ‘OIII4363,OIII4959’).

  • fesc (float) – The Lyman continuum escaped fraction, the fraction of ionising photons that entirely escaped.

  • mask (array) – A mask to apply to the particles (only applicable to particle)

  • method (str) – The method to use for the interpolation. Options are: ‘cic’ - Cloud in cell ‘ngp’ - Nearest grid point

Returns:

LineCollection

A dictionary like object containing line objects.

get_particle_line_screen(grid, line_ids, fesc=0.0, tau_v=None, dust_curve=<synthesizer.dust.attenuation.PowerLaw object>, mask=None, method='cic', label='')[source]

Get a LineCollection with screen attenuated lines for each particle.

Calculates attenuated properties (luminosity, continuum, EW) for a set of lines assuming a simple dust screen (i.e. both nebular and stellar emission feels the same dust attenuation). This is a wrapper around the more general method above.

Parameters:
  • grid (Grid) – The Grid object.

  • line_ids (list/str) – A list of line_ids or a str denoting a single line. Doublets can be specified as a nested list or using a comma (e.g. ‘OIII4363,OIII4959’).

  • fesc (float) – The Lyman continuum escaped fraction, the fraction of ionising photons that entirely escaped.

  • tau_v (float) – V-band optical depth.

  • dust_curve (dust_curve) – A dust_curve object specifying the dust curve.

  • mask (array) – A mask to apply to the particles (only applicable to particle)

  • method (str) – The method to use for the interpolation. Options are: ‘cic’ - Cloud in cell ‘ngp’ - Nearest grid point

Returns:

LineCollection

A dictionary like object containing line objects.

get_particle_spectra_incident(grid, young=None, old=None, label='', **kwargs)[source]

Generate the incident (equivalent to pure stellar for stars) spectra using the provided Grid.

Parameters:
  • grid (obj) – Spectral grid object.

  • young (unyt_quantity) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (unyt_quantity) – If not None, specifies age in Myr at which to filter for old star particles.

  • label (string) – A modifier for the spectra dictionary key such that the key is label + “_incident”.

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Returns:

Sed

An Sed object containing the stellar spectra.

get_particle_spectra_linecont(grid, fesc=0.0, fesc_LyA=1.0, young=None, old=None, **kwargs)[source]

Generate the line contribution spectra. This is only invoked if fesc_LyA < 1.

Parameters:
  • grid (obj) – Spectral grid object.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • fesc_LyA (float) – Fraction of Lyman-alpha emission that can escape unimpeded by the ISM/IGM.

  • young (unyt_quantity) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (unyt_quantity) – If not None, specifies age in Myr at which to filter for old star particles.

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Returns:

numpy.ndarray

The line contribution spectra.

get_particle_spectra_nebular(grid, fesc=0.0, young=None, old=None, label='', **kwargs)[source]

Generate nebular spectra from a grid object and star particles. The grid object must contain a nebular component.

Parameters:
  • grid (obj) – Spectral grid object.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • young (unyt_quantity) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (unyt_quantity) – If not None, specifies age in Myr at which to filter for old star particles.

  • label (string) – A modifier for the spectra dictionary key such that the key is label + “_nebular”.

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Returns:

Sed

An Sed object containing the nebular spectra.

get_particle_spectra_reprocessed(grid, fesc=0.0, fesc_LyA=1.0, young=None, old=None, label='', **kwargs)[source]

Generates the intrinsic spectra, this is the sum of the escaping radiation (if fesc>0), the transmitted emission, and the nebular emission. The transmitted emission is the emission that is transmitted through the gas. In addition to returning the intrinsic spectra this saves the incident, nebular, and escaped spectra if update is set to True.

Parameters:
  • grid (obj) – Spectral grid object.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • fesc_LyA (float) – Fraction of Lyman-alpha emission that can escape unimpeded by the ISM/IGM.

  • young (unyt_quantity) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (unyt_quantity) – If not None, specifies age in Myr at which to filter for old star particles.

  • label (string) – A modifier for the spectra dictionary key such that the key is label + “_transmitted”.

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Updates:

incident: transmitted nebular reprocessed intrinsic

if fesc>0:

escaped

Returns:

Sed

An Sed object containing the intrinsic spectra.

get_particle_spectra_screen(grid=None, fesc=0.0, tau_v=None, dust_curve=<synthesizer.dust.attenuation.PowerLaw object>, young=None, old=None, label='')[source]

Generates the dust attenuated spectra. First generates the intrinsic spectra if this hasn’t already been calculated.

Parameters:
  • grid (obj) – Spectral grid object.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Updates:

incident: transmitted nebular reprocessed intrinsic attenuated

if fesc>0:

escaped

Returns:

Sed

An Sed object containing the attenuated spectra.

get_particle_spectra_transmitted(grid, fesc=0.0, young=None, old=None, label='', **kwargs)[source]

Generate the transmitted spectra using the provided Grid. This is the emission which is transmitted through the gas as calculated by the photoionisation code.

Parameters:
  • grid (obj) – Spectral grid object.

  • fesc (float/array-like, float) – Fraction of stellar emission that escapes unattenuated from the birth cloud. Can either be a single value or an value per star (defaults to 0.0).

  • young (unyt_quantity) – If not None, specifies age in Myr at which to filter for young star particles.

  • old (unyt_quantity) – If not None, specifies age in Myr at which to filter for old star particles.

  • label (string) – A modifier for the spectra dictionary key such that the key is label + “_transmitted”.

  • kwargs – Any keyword arguments which can be passed to generate_particle_lnu.

Returns:

Sed

An Sed object containing the transmitted spectra.

get_sfzh(grid, grid_assignment_method='cic', nthreads=0)[source]

Generate the binned SFZH history of this collection of particles.

The binned SFZH produced by this method is equivalent to the weights used to extract spectra from the grid.

Parameters:
  • grid (Grid) – The spectral grid object.

  • grid_assignment_method (string) – The type of method used to assign particles to a SPS grid point. Allowed methods are cic (cloud in cell) or nearest grid point (ngp) or there uppercase equivalents (CIC, NGP). Defaults to cic.

  • nthreads (int) – The number of threads to use in the computation. If set to -1 all available threads will be used.

Returns:

Numpy array of containing the SFZH.

property log10ages

Return stellar particle ages in log (base 10)

Returns:

log10ages (array)

log10 stellar ages

property log10metallicities

Return stellar particle ages in log (base 10). Zero valued metallicities are set to metallicity_floor, which is set on initialisation of this stars object. To check it, run stars.metallicity_floor.

Returns:

log10metallicities (array)

log10 stellar metallicities.

plot_sfzh(grid, grid_assignment_method='cic', show=True)[source]

Plot the binned SZFH.

Parameters:

show (bool) – Should we invoke plt.show()?

Returns:

fig

The Figure object contain the plot axes.

ax

The Axes object containing the plotted data.

renormalise_mass(stellar_mass)[source]

Renormalises and overwrites the initial masses. Useful when rescaling the mass of the system of stellar particles.

Parameters:

stellar_mass (array-like, float) – The stellar mass array to be renormalised.

resample_young_stars(min_age=100000000.0, min_mass=700, max_mass=1000000.0, power_law_index=-1.3, n_samples=1000.0, force_resample=False, verbose=False)[source]

Resample young stellar particles into individual HII regions, with a power law distribution of masses. A young stellar particle is a stellar particle with an age < min_age (defined in Myr?).

This function overwrites the properties stored in attributes with the resampled properties.

Note: Resampling and imaging are not supported. If attempted an error

is thrown.

Parameters:
  • min_age (float) – The age below which stars will be resampled, in yrs.

  • min_mass (float) – The lower bound of the mass interval used in the power law sampling, in Msun.

  • max_mass (float) – The upper bound of the mass interval used in the power law sampling, in Msun.

  • power_law_index (float) – The index of the power law from which to sample stellar

  • n_samples (int) – The number of samples to generate for each stellar particles younger than min_age.

  • force_resample (bool) – A flag for whether resampling should be forced. Only applicable if trying to resample and already resampled Stars object.

  • verbose (bool) – Are we talking?