The QSO-Extension Module

The Sculptor Quasar Extension

This module defines models, masks and analysis routines specific for the analysis of type-I QSOs (quasars, type-I AGN).

At first we define the basic model functions and their setups. The setup functions initialize LMFIT models and parameters using the basic model functions defined here.

Complex models can be constructed by combining multiple of the basic model functions. For example, we define a setup function to initialize a model for the Hbeta and [OIII] lines consistent of six Gaussian emission lines.

sculptor_extensions.qso.CIII_complex_model_func(x, z, cen, cen_alIII, cen_siIII, flux, fwhm_km_s, shift_km_s, flux_alIII, fwhm_km_s_alIII, shift_km_s_alIII, flux_siIII, fwhm_km_s_siIII, shift_km_s_siIII)

Model function for the CIII] emission line complex, consisting of the Gaussian line models with a combined redshift parameter.

The width of the line is set by the FWHM in km/s.

The Gaussians are not normalized.

Parameters
  • x (np.ndarray) – Dispersion of the template model

  • z (float) – Redshift for CIII], AlIII, SiIII]

  • cen (float) – CIII] central wavelength

  • cen_alIII (float) – AlIII central wavelength

  • cen_siIII (float) – SiIII] central wavelength

  • amp (float) – Amplitude of the CIII] line

  • fwhm_km_s (float) – Full Width at Half Maximum (FWHM) of CIII] in km/s

  • shift_km_s (float) – Doppler velocity shift of the central wavelength

  • amp_alIII – Amplitude of the AlIII line

  • fwhm_km_s_alIII (float) – Full Width at Half Maximum (FWHM) of AlIII in km/s

  • shift_km_s_alIII – Doppler velocity shift of the central wavelength

  • amp_siIII – Amplitude of the SiIII] line

  • fwhm_km_s_siIII (float) – Full Width at Half Maximum (FWHM) of SiIII] in km/s

  • shift_km_s_siIII – Doppler velocity shift of the central wavelength

Returns

CIII] complex model

sculptor_extensions.qso.add_redshift_param(redshift, params, prefix)

Add the redshift parameter to the LMFIT parameters.

Parameters
  • redshift (float) – Redshift

  • params (lmfit.Parameters) – Model parameters

  • prefix (string) – Model prefix

Returns

sculptor_extensions.qso.balmer_continuum_model(x, z, amp_be, Te, tau_be, lambda_be)

Model of the Balmer continuum (Dietrich 2003)

This model is defined for a spectral dispersion axis in Angstroem.

This functions implements the Balmer continuum model presented in Dietrich 2003. The model follows a blackbody below the Balmer edge ( 3646A) and is zero above.

The amplitude parameter amp_be defines the amplitude at the balmer edge. The strength of the Balmer continuum can be estimated from the fluxden density at 3675A after subtraction of the power-law continuum component for reference see Grandi et al.(1982), Wills et al.(1985) or Verner et al.( 1999).

At wavelengths of 3646A higher order Balmer lines are merging. This has not been included in this model and thus it will produce a sharp break at the Balmer edge.

Parameters
  • x (np.ndarray) – Dispersion of the Balmer continuum

  • z (float) – Redshift

  • amp_be (float) – Amplitude of the Balmer continuum at the Balmer edge

  • Te (float) – Electron temperature

  • tau_be (float) – Optical depth at the Balmer edge

  • lambda_be (float) – Wavelength of the Balmer edge

Returns

Balmer continuum model

Return type

np.ndarray

sculptor_extensions.qso.calc_bolometric_luminosity(cont_lwav, cont_wav, reference='Shen2011')

Calculate the bolometric luminosity from the monochromatic continuum luminosity (erg/s/A) using bolometric correction factors from the literature.

The following bolometric corrections are available cont_wav = 1350, reference = Shen2011 cont_wav = 3000, reference = Shen2011 cont_wav = 5100, reference = Shen2011

The Shen et al. 2011 (ApJ, 194, 45) bolometric corrections are based on the composite spectral energy distribution (SED) in Richards et al. 2006 ( ApJ, 166,470).

Parameters
  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

  • reference (string) – A reference string to select from the available bolometric corrections.

Returns

Returns a tuple of the bolometric luminosity in erg/s and a reference string indicating the publication and continuum wavelength of the bolometric correction.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.calc_eddington_luminosity(bh_mass)

Calculate the Eddington luminosity for a given black hole mass.

Parameters

bh_mass (u.Quantity) – Black hole mass as an astropy Quantity

Returns

Returns the Eddington luminosity in cgs units (erg/s).

Return type

u.Quantity

sculptor_extensions.qso.calc_eddington_ratio(lbol, bh_mass)

Calculate the Eddington ratio for a provided bolometric luminosity and black hole mass.

Parameters
  • lbol (u.Quantity) – Bolometric luminosity

  • bh_mass (u.Quantity) – Black hole mass

Returns

sculptor_extensions.qso.calc_velocity_shifts(z, z_sys, relativistic=True)

Calculate the velocity difference of a feature redshift with respect to the systemic redshift.

This function is currently simply a wrapper around the linetools functions calculating the velocity difference.

Parameters
  • z (float) – The redshift of the spectroscopic feature (e.g., absorption or emission line).

  • z_sys (float) – The systemic redshift

  • relativistic – Boolean indicating whether the doppler velocity is calculated assuming relativistic velocities.

Type

bool

Returns

Returns the velocity difference in km/s.

Return type

u.Quantity

sculptor_extensions.qso.correct_CIV_fwhm_for_blueshift(civ_fwhm, blueshift)

Correct the CIV FWHM for the CIV blueshifts using the relation determined in Coatman et al. (2017, MNRAS, 465, 2120).

The correction follows Eq. 4 of Coatman et al. (2017).

Parameters
  • civ_fwhm (astropy.units.Quantity) – FWHM of the CIV line in km/s

  • blueshift (astropy.units.Quantity) – Blueshift of the CIV line in km/s

Returns

Corrected FWHM in km/s

Return type

astropy.units.Quantity

sculptor_extensions.qso.line_model_gaussian(x, z, flux, cen, fwhm_km_s)

Gaussian line model

The central wavelength of the Gaussian line model is determined by the central wavelength cen and the redshift, z. These parameters are degenerate in a line fit and it is adviseable to fix one of them (to predetermined values e.g., the redshift or the central wavelength).

The width of the line is set by the FWHM in km/s.

The Gaussian is normalized.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • z (float) – Redshift

  • flux (float) – Amplitude of the Gaussian

  • cen (float) – Central wavelength

  • fwhm_km_s (float) – FWHM of the Gaussian in km/s

Returns

Gaussian line model

Return type

np.ndarray

sculptor_extensions.qso.line_model_gaussian_nii_doublet(x, z, flux_a, flux_b, fwhm_km_s_a, fwhm_km_s_b)

Doublet line model for the [NII] lines at 6549.85 A and 6585.28 A.

This model ties the redshift of the forbidden transitions of [NII] at 6549.85 A and 6585.28 A together. FWHM and fluxes are free parameters for each Gaussian line model.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • z (float) – Redshift

  • flux_a (float) – Amplitude of the 1st Gaussian component

  • flux_b (float) – Amplitude of the 2nd Gaussian component

  • fwhm_km_s_a (float) – FWHM of the 1st Gaussian component in km/s

  • fwhm_km_s_b (float) – FWHM of the 2nd Gaussian component in km/s

Returns

Gaussian doublet line model

Return type

np.ndarray

sculptor_extensions.qso.line_model_gaussian_oiii_doublet(x, z, flux, fwhm_km_s, fluxratio)

Doublet line model for the [OIII] lines at 4960.30 A and 5008.24 A.

This model ties the redshift, the FWHM and the fluxes (ratio 1:3) of the forbidden transitions of [OIII] at 4960.30 A and 5008.24 A together.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • z (float) – Redshift

  • flux (float) – Amplitude of the Gaussian

  • cen (float) – Central wavelength

  • fwhm_km_s (float) – FWHM of the Gaussian in km/s

Returns

Gaussian doublet line model

Return type

np.ndarray

sculptor_extensions.qso.line_model_gaussian_sii_doublet(x, z, flux_a, flux_b, fwhm_km_s_a, fwhm_km_s_b)

Doublet line model for the [SII] lines at 6718.29 A and 6732.67 A.

This model ties the redshift, of the forbidden transitions of [SII] at 6718.29 A and 6732.67 A together. FWHM and fluxes are free parameters for each Gaussian line model.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • z (float) – Redshift

  • flux_a (float) – Amplitude of the 1st Gaussian component

  • flux_b (float) – Amplitude of the 2nd Gaussian component

  • fwhm_km_s_a (float) – FWHM of the 1st Gaussian component in km/s

  • fwhm_km_s_b (float) – FWHM of the 2nd Gaussian component in km/s

Returns

Gaussian doublet line model

Return type

np.ndarray

sculptor_extensions.qso.power_law_at_2500(x, amp, slope, z)

Power law model anchored at 2500 AA

This model is defined for a spectral dispersion axis in Angstroem.

Parameters
  • x (np.ndarray) – Dispersion of the power law

  • amp (float) – Amplitude of the power law (at 2500 A)

  • slope (float) – Slope of the power law

  • z (float) – Redshift

Returns

Power law model

Return type

np.ndarray

sculptor_extensions.qso.power_law_at_2500_plus_bc(x, amp, slope, z, amp_be, Te, tau_be, lambda_be)

QSO continuum model consisting of a power law anchored at 2500 A and a balmer continuum contribution.

This model is defined for a spectral dispersion axis in Angstroem.

The amplitude of the Balmer continuum is set independently of the power law component by the amplitude of the balmer continuum at the balmer edge amp_be.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • amp (float) – Amplitude of the power law (at 2500 A)

  • slope (float) – Slope of the power law

  • z (float) – Redshift

  • amp_be (float) – Amplitude of the Balmer continuum at the Balmer edge

  • Te (float) – Electron temperature

  • tau_be (float) – Optical depth at the Balmer edge

  • lambda_be (float) – Wavelength of the Balmer edge

Returns

QSO continuum model

Return type

np.ndarray

sculptor_extensions.qso.power_law_at_2500_plus_fractional_bc(x, amp, slope, z, f, Te, tau_be, lambda_be)

QSO continuum model consisting of a power law anchored at 2500 A and a balmer continuum contribution.

This model is defined for a spectral dispersion axis in Angstroem.

The amplitude of the Balmer continuum is set to be a fraction of the power law component at the Balmer edge (3646A) using the variabe f.

Parameters
  • x (np.ndarray) – Dispersion of the continuum model

  • amp (float) – Amplitude of the power law (at 2500 A)

  • slope (float) – Slope of the power law

  • z (float) – Redshift

  • f (float) – Amplitude of the Balmer continuum as a fraction of the power law component

  • Te (float) – Electron temperature

  • tau_be (float) – Optical depth at the Balmer edge

  • lambda_be (float) – Wavelength of the Balmer edge

Returns

QSO continuum model

Return type

np.ndarray

sculptor_extensions.qso.se_bhmass_civ_c17_fwhm(civ_fwhm, cont_lwav, cont_wav)

Calculate the single-epoch virial BH mass based on the CIV FWHM and monochromatic continuum luminosity at 1350A.

This relationship follows Eq.6 of from Coatman et al. (2017, MNRAS, 465, 2120)

The FWHM of the CIV line was corrected by the CIV blueshift using their Eq. 4. In this study the CIV line was modeled with sixth order Gauss-Hermite (GH) polynomials using the normalization of van der Marel & Franx (1993) and the functional forms of Cappelari et al. (2002). GH polynomials were chose because they are flexble enough to model the often very asymmetric CIV line profile.

The authors state that using commonly employed three Gaussian components, rather than the GH polynomials, resulted in only marginal differences in the line parameters.

Parameters
  • civ_fwhm (astropy.units.Quantity) – FWHM of the CIV line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity at 1350A in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

Returns

Returns a tuple of the BH mass estimate based on the CIV FWHM and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.se_bhmass_civ_vp06_fwhm(civ_fwhm, cont_lwav, cont_wav)

Calculate the single-epoch virial BH mass based on the CIV FWHM and monochromatic continuum luminosity at 1350A.

The monochromatic continuum luminosity at 1450A can be used without error or penalty in accuracy.

This relationship is taken from Vestergaard & Peterson 2006, ApJ 641, 689

The FWHM of the CIV line was measured with the methodology described in Peterson et al. 2004. The line width measurements to establish the CIV single-epoch relation are corrected for resolution effects as described in Peterson et al. 2004.

“The sample standard deviation of the weighted average zero point offset, which shows the intrinsic scatter in the saample is +-0.36 dex. This value is more representative of the uncertainty zer point than is the formal error.”

Parameters
  • civ_fwhm (astropy.units.Quantity) – FWHM of the CIV line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity at 1350A/1450A in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

Returns

Returns a tuple of the BH mass estimate based on the CIV FWHM and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.se_bhmass_civ_vp06_sigma(civ_sigma, cont_lwav, cont_wav)

Calculate the single-epoch virial BH mass based on the CIV line dispersion (sigma) and monochromatic continuum luminosity at 1350A.

The monochromatic continuum luminosity at 1450A can be used without error or penalty in accuracy.

This relationship is taken from Vestergaard & Peterson 2006, ApJ 641, 689

The FWHM of the CIV line was measured with the methodology described in Peterson et al. 2004. The line width measurements to establish the CIV single-epoch relation are corrected for resolution effects as described in Peterson et al. 2004.

Peterson et al. (2004) note a number of advantages to using sigma rathern than the FWHM as the line width measure.

“The sample standard deviation of the weighted average zero point offset, which shows the intrinsic scatter in the saample is +-0.33 dex. This value is more representative of the uncertainty zer point than is the formal error.”

Parameters
  • civ_sigma (astropy.units.Quantity) – Line dispersion (sigma) of the CIV line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity at 1350A/1450A in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

Returns

Returns a tuple of the BH mass estimate based on the CIV sigma and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.se_bhmass_hbeta_vp06(hbeta_fwhm, cont_lwav, cont_wav=<Quantity 5100. A>)

Calculate the single-epoch virial BH mass based on the Hbeta FWHM and monochromatic continuum luminosity at 5100A.

This relationship is taken from Vestergaard & Peterson 2006, ApJ 641, 689

Note that the Hbeta line width to establish this single-epoch virial estimator was established by using the FWHM of only the broad component. The line width was corrected for the spectral resolution as described in Peterson et al. 2004.

The relationship is based on line width measurements of quasars published in Boroson & Green 1992 and Marziani 2003.

“The sample standard deviation of the weighted average zero point offset, which shows the intrinsic scatter in the saample is +-0.43 dex. This value is more representative of the uncertainty zer point than is the formal error.”

Parameters
  • hbeta_fwhm (astropy.units.Quantity) – FWHM of the Hbeta line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity at 5100A in erg/s/A

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity (default = 5100A).

Returns

Returns a tuple of the BH mass estimate based on the Hbeta FWHM and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.se_bhmass_mgii_s11_fwhm(mgii_fwhm, cont_lwav, cont_wav)

Calculate the single-epoch virial BH mass based on the MgII FWHM and monochromatic continuum luminosity at 3000A.

This relationship is taken from Shen et al. 2011, ApJ, 194, 45

To model the FeII contribution beneath MgII line the authors use empirical FeII templates from Borosn & Grenn 1992, Vestergaard & Wilkes 2001, and Salviander 2007.

Salviander modified the Vestergaard & Wilkes (2001) template in the region of 2780–2830A centered on MgII, where the Vestergaard & Wilkes (2001) template is set to zero. For this region, they incorporate a theoretical FeII model ofSigut & Pradhan (2003) scaled to match the Vestergaard & Wilkes(2001) template at neighboring wavelengths.

As the subtraction of the underlying FeII continuum can have systematic effects on the measurement of the MgII FWHM and therefore the BH mass estimate it is advised to always employ the same continuum construction procedure as in the reference sample that established the single-epoch virial relationship.

Parameters
  • mgii_fwhm (astropy.units.Quantity) – FWHM of the MgII line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

Returns

Returns a tuple of the BH mass estimate based on the MgII FWHM and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.se_bhmass_mgii_vo09_fwhm(mgii_fwhm, cont_lwav, cont_wav)

Calculate the single-epoch virial BH mass based on the MgII FWHM and monochromatic continuum luminosity at 1350A, 2100A, 3000A, or 5100A.

This relationship is taken from Vestergaard & Osmer 2009, ApJ 641, 689

To determine the FWHM of the MgII line, the authors modeled the FeII emission beneath the MgII line with the Vestergaard & Wilkes 2001 and the Boroson & Green 1992 iron templates.

Most of the MgII lines were modeled with a single Gaussian component, in cases of high-quality spectra two Gaussian components were used. For the single-Gaussian components the authors adopted the measurements of the FWHM and uncertainties tabulated by Forster et al. (their Table 5). For multi-Gaussian components the FWHM of the full modeled profile was measured.

As the subtraction of the underlying FeII continuum can have systematic effects on the measurement of the MgII FWHM and therefore the BH mass estimate it is advised to always employ the same continuum construction procedure as in the reference sample that established the single-epoch virial relationship.

“The absolute 1 sigma scatter in the zero points is 0.55dex, which includes the factor ~2.9 uncertainties of the reverberation mapping masses to which these mass estimation relations are anchored (see Vestergaard & Peterson 2006 and Onken et al. 2004 for details)”

Parameters
  • mgii_fwhm (astropy.units.Quantity) – FWHM of the MgII line in km/s

  • cont_lwav (astropy.units.Quantity) – Monochromatic continuum luminosity in erg/s/A.

  • cont_wav (astropy.units.Quantity) – Wavelength of the monochromatic continuum luminosity in A.

Returns

Returns a tuple of the BH mass estimate based on the MgII FWHM and a reference string for the single-epoch scaling relationship.

Return type

astropy.units.Quantity, string

sculptor_extensions.qso.setup_SWIRE_Ell2_template(prefix, **kwargs)

Setup the SWIRE library Ell2 galaxy template model

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_SWIRE_NGC6090_template(prefix, **kwargs)

Setup the SWIRE library NGC6090 galaxy template model

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_doublet_line_model_nii(prefix, **kwargs)

Set up a double line model for the [NII] emission lines at 6549.85 A and 6585.28 A.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters

prefix – The name of the doublet line model. If prefix is None then a

pre-defined name will be assumed. :type prefix: string :param kwargs: :return: Return a list of LMFIT models and a list of LMFIT parameters :rtype: (list, list)

sculptor_extensions.qso.setup_doublet_line_model_oiii(prefix, **kwargs)

Set up a double line model for the [OIII] emission lines at 4960.30 A and 5008.24 A.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters

prefix – The name of the doublet line model. If prefix is None then a

pre-defined name will be assumed. :type prefix: string :param kwargs: :return: Return a list of LMFIT models and a list of LMFIT parameters :rtype: (list, list)

sculptor_extensions.qso.setup_doublet_line_model_sii(prefix, **kwargs)

Set up a double line model for the [SII] emission lines at 6718.29 A and 6732.67 A.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters

prefix – The name of the doublet line model. If prefix is None then a

pre-defined name will be assumed. :type prefix: string :param kwargs: :return: Return a list of LMFIT models and a list of LMFIT parameters :rtype: (list, list)

sculptor_extensions.qso.setup_galaxy_template_model(prefix, template_filename, templ_disp_unit, templ_fluxden_unit, fwhm=2500, redshift=0, amplitude=1, intr_fwhm=900, dispersion_limits=None)

Initialize a galaxy template model

Parameters
  • prefix (string) – Model prefix

  • template_filename (string) – Filename of the iron template

  • fwhm (float) – FWHM the template should be broadened to

  • redshift (float) – Redshift

  • amplitude (float) – Amplitude of the template model

  • intr_fwhm (float) – Intrinsic FWHM of the template

  • dispersion_limits ((float, float)) –

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_MgII_T06(prefix, **kwargs)

Setup the Tsuzuki 2006 iron template model around MgII (2200-3500A)

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_MgII_VW01(prefix, **kwargs)

Setup the Vestergaard & Wilkes 2001 iron template model around MgII (2200-3500A)

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_OPT_BG92(prefix, **kwargs)

Setup the Boroson & Green 1992 iron template model around Hbeta (3700-7480A)

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_T06(prefix, **kwargs)

Setup the Tsuzuki 2006 iron template model

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_UV_VW01(prefix, **kwargs)

Setup the Vestergaard & Wilkes 2001 iron template model around CIV (1200-2200A)

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_iron_template_model(prefix, template_filename, templ_disp_unit, templ_fluxden_unit, fwhm=2500, redshift=0, amplitude=1, intr_fwhm=900, dispersion_limits=None)

Initialize an iron template model

Parameters
  • prefix (string) – Model prefix

  • template_filename (string) – Filename of the iron template

  • fwhm (float) – FWHM the template should be broadened to

  • redshift (float) – Redshift

  • amplitude (float) – Amplitude of the template model

  • intr_fwhm (float) – Intrinsic FWHM of the template

  • dispersion_limits ((float, float)) –

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_line_model_CIII_complex(prefix, **kwargs)

Set up a 3 component Gaussian line model for the CIII], AlIII and SiIII] emission lines.

Note that a special model function exists as all three Gaussian line models share a common redshift parameter.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return model and model parameters

sculptor_extensions.qso.setup_line_model_CIV_2G(prefix, **kwargs)

Set up a 2 component Gaussian line model for the CIV emission line.

This model is defined for a spectral dispersion axis in Angstroem.

This setup models the broad CIV line emission as seen in type-I quasars and AGN. Due to the broad nature of the line CIV the CIV doublet is assumed to be unresolved.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_line_model_Halpha_2G(prefix, **kwargs)

Set up a 2 component Gaussian line model for the HAlpha emission line.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters

prefix – The name of the emission line model. If prefix is None then a

pre-defined name will be assumed. :type prefix: string :param kwargs: :return: Return a list of LMFIT models and a list of LMFIT parameters :rtype: (list, list)

sculptor_extensions.qso.setup_line_model_Hbeta_2G(prefix, **kwargs)

Set up a 2 component Gaussian line model for the HBeta emission line.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters

prefix – The name of the emission line model. If prefix is None then a

pre-defined name will be assumed. :type prefix: string :param kwargs: :return: Return a list of LMFIT models and a list of LMFIT parameters :rtype: (list, list)

sculptor_extensions.qso.setup_line_model_MgII_2G(prefix, **kwargs)

Set up a 2 component Gaussian line model for the MgII emission line.

This model is defined for a spectral dispersion axis in Angstroem.

This setup models the broad MgII line emission as seen in type-I quasars and AGN. Due to the broad nature of the line MgII the MgII doublet is assumed to be unresolved.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_line_model_SiIV_2G(prefix, **kwargs)

Set up a 2 component Gaussian line model for the SiIV emission line.

This setup models the broad SiIV line emission as seen in type-I quasars and AGN. Due to the broad nature of the line SiIV, the SiIV doublet is assumed to be unresolved and blended with the OIV emission line.

This model is defined for a spectral dispersion axis in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_line_model_gaussian(prefix, **kwargs)

Initialize the Gaussian line model.

Parameters
  • prefix (string) – Model prefix

  • kwargs – Keyword arguments

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_power_law_at_2500(prefix, **kwargs)

Initialize the power law model anchored at 2500 A.

Parameters
  • prefix (string) – Model prefix

  • kwargs – Keyword arguments

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_power_law_at_2500_plus_bc(prefix, **kwargs)

Initialize the quasar continuum model consistent of a power law anchored at 2500A and a balmer continuum contribution.

Parameters
  • prefix (string) – Model prefix

  • kwargs – Keyword arguments

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_power_law_at_2500_plus_fractional_bc(prefix, **kwargs)

Initialize the quasar continuum model consistent of a power law anchored at 2500A and a balmer continuum contribution.

The Balmer continuum amplitude at the Balmer edge is set to be a fraction of the power law component.

Parameters
  • prefix (string) – Model prefix

  • kwargs – Keyword arguments

Returns

LMFIT model and parameters

Return type

(lmfit.Model, lmfit.Parameters)

sculptor_extensions.qso.setup_split_iron_template_MgII_T06(prefix, **kwargs)

Setup the Tsuzuki 2006 iron template model subdivided into three separate models at 2200-2660, 2660-3000, 3000-3500.

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_split_iron_template_MgII_VW01(prefix, **kwargs)

Setup the Vestergaard & Wilkes 2001 iron template model subdivided into three separate models at 2200-2660, 2660-3000, 3000-3500.

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_split_iron_template_OPT_BG92(prefix, **kwargs)

Setup the Boroson & Green 1992 iron template model subdivided into three separate models at 3700-4700, 4700-5100, 5100-5600.

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_split_iron_template_UV_VW01(prefix, **kwargs)

Setup the Vestergaard & Wilkes 2001 iron template model subdivided into three separate models at 1200-1560, 1560-1875, 1875-2200.

The dispersion axis for this model is in Angstroem.

Parameters
  • prefix (string) – The input parameter exists for conformity with the Sculptor models, but will be ignored. The prefix is automatically set by the setup function.

  • kwargs

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.setup_subdivided_iron_template(templ_list, fwhm=2500, redshift=0, amplitude=1)

Setup iron template models from a predefined list of templates and dispersion ranges.

Parameters
  • templ_list – List of template names for which models will be set up

  • fwhm (float) – Goal FWHM of the template model

  • redshift (float) – Redshift of the template model

  • amplitude (float) – Amplitude of the template model

Type

list

Returns

Return a list of LMFIT models and a list of LMFIT parameters

Return type

(list, list)

sculptor_extensions.qso.template_model(x, amp, z, fwhm, intr_fwhm, templ_disp=None, templ_fluxden=None, templ_disp_unit_str=None, templ_fluxden_unit_str=None)

Template model

Parameters
  • x (np.ndarray) – Dispersion of the template model

  • amp (float) – Amplitude of the template model

  • z (float) – Redshift

  • fwhm (float) – FWHM the template should be broadened to

  • intr_fwhm (float) – Intrinsic FWHM of the template

  • templ_disp (np.ndarray) – Dispersion axis of the template. This must match the same dispersion unit as the model

  • templ_fluxden (templ_fluxden: np.ndarray) – Flux density of the template.

  • templ_disp_unit_str (str) – Dispersion unit of the template as a string in astropy cds format.

  • templ_fluxden_unit – Flux density unit of the template as a string in astropy cds format.

Returns

Template model as a Scipy interpolation