galapy.Galaxy

Implements galaxy classes wrapping up the interplay between the different components contribution to the overall emission. By instantiating an object of class GXY or derived it is possible to set-up the different components, compute global and component-specific emission, compute the derived quantities and modify the parameter setting.

Classes

GXY([age, redshift, cosmo, lstep, do_Xray, ...])

Base galaxy class.

PhotoGXY(*args[, pms])

Galaxy class including photometric system.

SpectralGXY([age, redshift, cosmo, lstep, ...])

Not implemented yet

class galapy.Galaxy.GXY(age=None, redshift=None, cosmo='Planck18', lstep=None, do_Xray=False, do_Radio=False, do_AGN=False, do_IGM=False, sfh=None, csp=None, ism=None, agn=None, nff=None, syn=None)

Base galaxy class. All other galaxy models derive from this object. All the arguments are optional. Whether a particular component is included or not (e.g. radio, x-ray) the wavelength grid will always span from 1 to \(10^{10}\) angstroms. Note that by default, the only components that will be built are

  1. star formation history (SFH): modelling the stellar mass growth

  2. composite stellar populations (CSP): modelling unattenuated stellar emission

  3. inter-stellar medium (ISM): modelling the absorption/re-radiation due to inter-stellar dust. This is divided in two components: molecular clouds (MC) and diffuse dust (DD)

All the other components are not strictly necessay to generate an SED and, for the sake of performances, they should be built only if necessary.

Parameters:
  • age (float) – (None defaults to 1.e+6) age in years of the galaxy at the time it is observed

  • redshift (float) – (None defaults to 0.0) redshift of the galaxy

  • cosmo (str or dict) – (default = 'Planck18') cosmological model of choice, takes the same arguments an object of type galapy.Cosmology.CSM would take. If a string is passed, it should name one of the pre-computed cosmologies available in the database. Available cosmologies are 'WMAP7', 'WMAP9', 'Planck15', 'Planck18'. If a dictionary is passed, the class expects to find 3 key-value couples: key = 'redshift', value = an array of redshift values; key = 'luminosity_distance', value = an array of luminosity distances corresponding to the redshift values of the first key-value couple; key = 'age', value = an array of ages of the Universe corresponding to the redshift values of the first key-value couple. All these arrays should have the same length.

  • lstep (scalar int or boolean mask) – (default = None) Reduces the granularity of the wavelength grid. If the input is scalar, select one any lstep wavelengths. If the input is a boolean mask, only select the indexes of the array masked with lstep.

  • do_Xray (bool) – (default = False) build the components modelling X-Ray emission, i.e. X-Ray binaries (low & high mass) and, if do_AGN=True, the X-Ray part of the spectrum from the eventual AGN.

  • do_Radio (bool) – (default = False) build the components modelling Radio emission if necessary. If the chosen SSP library does not already include nebular and synchrotron emission, do_Radio = True will build objects of type NFF (nebular free-free) and SNSYN (super-nova synchrotron). The eventual parameters of both components can be passed through arguments nff and syn, respectively. Note that if do_Radio = False both nff and syn will be ignored.

  • do_AGN (bool) – (default = False) build AGN component. The AGN spectrum is modelled through templates from Fritz et al. (2006) in the UV to IR bands and with an attenuated power-law in the X-ray band (only if also argument do_Xray is set to True. Tunable parameters can be set by passing them to the argument agn.

  • do_IGM (bool) – (default = False) build attenuation due to inter-galactic hydrogen with model from Inoue et al., (2014).

  • sfh (dict) – (default = None) arguments passed to the SFH object builder, if None is passed uses the default parameter of class galapy.StarFormationHistory.SFH

  • csp (dict) – (default = None) arguments passed to the CSP object builder, if None is passed uses the default parameter of class galapy.CompositeStellarPopulation.CSP

  • ism (dict) – (default = None) arguments passed to the ISM object builder, if None is passed uses the default parameter of class galapy.InterStellarMedium.ISM

  • agn (dict) – (default = None) arguments passed to the AGN object builder, if None is passed uses the default parameter of class galapy.ActiveGalacticNucleus.AGN. If do_AGN=False this argument is ignored.

  • nff (dict) – (default = None) arguments passed to the NFF object builder, if None is passed uses the default parameter of class galapy.NebularFreeFree.NFF. If do_Radio=False this argument is ignored.

  • syn (dict) – (default = None) arguments passed to the SNSYN object builder, if None is passed uses the default parameter of class galapy.Synchrotron.SNSYN. If do_Radio=False this argument is ignored.

wl(obs=False)

Returns the wavelength grid with the mask applied.

Parameters:

obs (bool) – (Optional, default = False) if set to True returns the observer’s frame wavelength grid, otherwise the rest-frame grid is returned.

Returns:

  • 1d-array

  • wavelength grid in observer’s frame (obs = True) or in

  • rest frame (obs = False)

get_wavelength_grid(lstep)

Reduces the granularity of the wavelength grid [optimization]

Parameters:

lstep (scalar int or boolean mask) – If the input is scalar, select one any lstep wavelengths. If the input is a boolean mask, only select the indexes of the array masked with lstep.

Returns:

A list of indices of the wavelength grid

Return type:

integer array

set_parameters(age=None, redshift=None, sfh=None, ism=None, agn=None, nff=None, syn=None)

Preferred method for changing the value of free-parameters.

Parameters:
  • age (float) – age of the galaxy

  • redshift (float) – redshift of the galaxy

  • sfh (dict) – (default = None) arguments passed to SFH.set_parameters method, if None is passed maintains the current parameterisation

  • csp (dict) – (default = None) arguments passed to the CSP.set_parameters method, if None is passed maintains the current parameterisation

  • ism (dict) – (default = None) arguments passed to the ISM.set_parameters method, if None is passed maintains the current parameterisation

  • agn (dict) – (default = None) arguments passed to the AGN.set_parameters method, if None is passed maintains the current parameterisation If the GXY object has been built without AGN component this argument is ignored.

  • nff (dict) – (default = None) arguments passed to the NFF.set_parameters method, if None is passed maintains the current parameterisation If the GXY object has been built without AGN component this argument is ignored.

  • syn (dict) – (default = None) arguments passed to the SNSYN.set_parameters method, if None is passed maintains the current parameterisation If the GXY object has been built without AGN component this argument is ignored.

Note

The method is built to optimise the number of computations performed. Do not pass arguments that would not change the current parameterisation. e.g.1 passing sfh = {} is less performant than sticking to the default sfh = None. e.g.2 passing the same value at each call is a waste of computational time: redshift = 2.0 passed at each call will considerably slow down execution. Bottom line: pass an argument only when necessary.

Lstellar()

Unattenuated stellar emission in solar luminosities. Approximates the integral

\[L_\lambda^\text{CSP}(\tau') = \int_0^{\tau'}\text{d}\tau L_\lambda^\text{SSP}\bigl[\tau, Z_\ast(\tau'-\tau)\bigr]\psi(\tau'-\tau)\]

where \(\tau'\) is the age of the galaxy, \(L_\lambda^\text{SSP}[\tau, Z\ast]\) is the luminosity of the Simple Stellar Population at given time \(\tau\) and at given stellar metallicity \(Z_\ast\), \(\psi(\tau)\) is the Star Formation History at time \(\tau\)

Returns:

  • 1-d array

  • stellar emission on the default wavelength grid

get_emission(store_attenuation=False, **kwargs)

Computes the overall emission coming from a galaxy with given parameterisation. The resulting shape of the SED depends on the components the galaxy object has been built with. This authomatically deals with the interplay among the different active components.

Parameters:
  • store_attenuation (bool) – (Optional, default = False) if set to True, stores the total, wavelength dependent, attenuation due to ISM in an internal variable (Aavg)

  • **kwargs (dictionary, optional) – arguments passed to function set_parameters()

Returns:

the emission on the selected wavelength grid in units of solar luminosities (\([L_\odot]\))

Return type:

1d-array

Note

Even though only the overall emission is returned, the contribution of each component is stored in the internal dictionary attribute components

get_avgAtt()

Returns the average attenuation in absolute magnitudes

get_SED()

Returns the flux at given distance in units of milli-Jansky [mJy].

\[S(\lambda_O) = \lambda_R^2 \cdot \dfrac{L_\text{tot}(\lambda_R) (1 + z)}{4\;\pi\;c\;D_L^2(z)} \cdot e^{-\tau_\text{IGM}(z)}\]

where \(\lambda_R\) and \(\lambda_O\) are the rest-frame and observer’s frame wavelength, respectively, \(L_\text{tot}\) is the total luminosity (as returned by function get_emission()), \(z\) is the redshift, \(c\) the light-speed, \(D_L(z)\) is the luminosity distance at observed redshift and, if it is included in the model, \(e^{-\tau_\text{IGM}(z)}\) is the IGM transmission at observed redshift.

Returns:

  • 1d-array

  • SED flux in milli-Jansky

components_to_flux()

Utility function converting emissions in the internal components dictionary to fluxes.

Returns:

  • dict

  • copy of the internal dictionary components with the emission conveted

  • to fluxes in milli-Jansky.

class galapy.Galaxy.PhotoGXY(*args, pms=None, **kwargs)

Galaxy class including photometric system. This is a class derived from galapy.Galaxy which implements authomatic computation of fluxes convolved with bandpass transmission filters.

Parameters:
  • pms (PMS instance, optional) – an instance of type galapy.PhotometricSystem.PMS, if not passed, it should be built using function build_photometric_system

  • *args (tuple, optional) – arguments to be passed to the constructor of the base class GXY

  • **kwargs (dictionary, optional) – keyword arguments to be passed to the constructor of the base class GXY

See also

galapy.PhotometricSystem.PMS

class implementing the photometric system

build_photometric_system

method for building a photometric system internally

GXY

base class

build_photometric_system(*args, **kwargs)

Forwards args and kwargs to the constructor of the photometric system.

Parameters:
  • *args (sequence) – positional arguments of the PMS contructor

  • **kwargs (dictionary) – keyword arguments of the PMS constructor

See also

galapy.PhotometricSystem.PMS

constructor of the photometric system class

photoSED()

Computes and returns the photometric bandpass fluxes in milliJansky.

class galapy.Galaxy.SpectralGXY(age=None, redshift=None, cosmo='Planck18', lstep=None, do_Xray=False, do_Radio=False, do_AGN=False, do_IGM=False, sfh=None, csp=None, ism=None, agn=None, nff=None, syn=None)

Not implemented yet