pyFIT3D.modelling package

pyFIT3D.modelling.dust module

pyFIT3D.modelling.dust.Calzetti_extlaw(l, R_V=4.05)

Calculates the reddening law by Calzetti et al. (1994). Original comments in the .for file:

q = A_lambda/A_V for Calzetti et al reddening law (formula from hyperz-manual). l = lambda, in Angstrons x = 1/lambda in 1/microns

Parameters
  • l (array like) – Wavelength in Angstroms.

  • R_V (float, optional) – Selective extinction parameter (roughly “slope”). Default value is 4.05.

Returns

Extinction A_lamda/A_V. Array of same length as l.

Return type

array like

pyFIT3D.modelling.dust.Cardelli_extlaw(l, R_V=3.1)

Calculates q(lambda) to the Cardelli, Clayton & Mathis 1989 reddening law.

Parameters
  • l (array like) – Wavelenght vector (Angstroms)

  • R_V (float, optional) – Selective extinction parameter (roughly “slope”). Default is R_V = 3.1.

Returns

q(lambda) of the reddedning law.

Return type

array like

pyFIT3D.modelling.dust.Charlot_extlaw(l, mu=0.3)

Returns two-component dust model by Charlot and Fall 2000.

Parameters
  • l (array like) – The lambda wavelength array (in Angstroms).

  • mu (float, optional) – Fraction of extinction contributed by the ambient ISM. Default value is 0.3.

Returns

  • array like – au_lambda for t <= 10^7 yr.

  • array like – au_lambda for t > 10^7 yr.

pyFIT3D.modelling.dust.calc_extlaw(l, R_V, extlaw=None)

A wrapper to select which extintion law to use if needed by user.

Parameters
  • l (array like) – Wavelength in Angstroms.

  • R_V (float) – Selective extinction parameter (roughly “slope”).

  • extlaw (str {'CCM', 'CAL'}, optional) – Which extinction function to use. CCM will call Cardelli_extlaw. CAL will call Calzetti_extlaw. Default value is CCM.

Returns

q(lambda) of the extlaw extinction law.

Return type

array like

pyFIT3D.modelling.dust.spec_apply_dust(wave, flux, AV, R_V=None, extlaw=None)

Apply the dust extinction law factor to a spectrum flux at wavelenghts wave.

flux_intrinsic = flux * 10^(-0.4 * AV * q_l)

Parameters
  • wave (array like) –

    Wavelenghts at restframe.

    (1 + redshift) = wave / wave_at_restframe

  • flux (array like) – Fluxes to be dust-corrected.

  • AV (float or array like) – Dust extinction in mag. If AV is an array, will create an (N_AV, N_wave) array of dust spectra.

  • R_V (float or None) – Selective extinction parameter (roughly “slope”). Default value 3.1.

  • extlaw (str {'CCM', 'CAL'} or None) – Which extinction function to use. CCM will call Cardelli_extlaw. CAL will call Calzetti_extlaw. Default value is CCM.

Returns

The flux correct by dust extinction.

Return type

array like

pyFIT3D.modelling.gas module

class pyFIT3D.modelling.gas.EmissionLines(wavelength, flux, sigma_flux, config)

Bases: object

Class created to fit combined models for a system of emission lines in a spectrum.

If the user wants to add other models should create a method like self._elines_model() or self._poly1d_model().

config
best_config
_config_list
wavelength
flux
sigma_flux
update_config()
create_system_model()
_create_single_model()
create_system_model(normed=False, ignore_poly1d=False, parameters=None)

Return the emission system model evaluated at the given wavelengths at the given parameters.

This function evaluate the all the models found in the config file, i.e., the models defined for the emission line system, taking into account the links and fixed parameters.

Parameters
  • normed (bool) – Sets the peak of eline models to 1.

  • ignore_poly1d (bool) – Ignore any poly1d model.

  • parameters (array like) – The list of parameters to define the current model.

Returns

The model evaluated at the given wavelengths using the given parameters.

Return type

array like

redefine_max_flux(flux_max=None, redefine_min=False, max_factor=1.2, min_factor=0.012)

Redefines the max value of the flux in eline models.

Parameters
  • flux_max (double, optional) – Max flux value. Defaults to self.flux.max().

  • redefine_min (bool, optional) – If True redefines the minimum treshold of the parameters range. Defaults to False.

  • max_factor (double, optional) – Defines the factor to multiply the maximum value of the integrated flux. It defaults to 1.2.

  • min_factor (double, optional) – Defines the factor to multiply the minimum value of the integrated flux. It defaults to 0.012.

update_chi_square(residuals, inplace=False)

Return the reduced chi square computed from the given residuals

The residuals are assumed to be weighted by the standard deviation of the observations, so that:

chi**2 = sum( residuals**2 ) / (n_residuals - 1 - n_free_param)
Parameters
  • residuals (array like) – A vector of residuals between an assumed model and a set of observations, weighted by the standard deviation on those observations.

  • inplace (boolean) – Whether to update self.latest_chi_sq and self.chi_sq_chain inplace or not. Defaults to False.

Returns

The reduced chi square score

Return type

float

update_config(new_guess=None, update_ranges=False, frac_range=0.2)

Update the original config parameters guess to the given fitted parameters

This method will update the guess parameters so that the config object can be used as a guess of a new fitting procedure.

Parameters
  • new_guess (array like) – A list of the fitted parameters to be used as new guess. If None (default), the latest fit will be used as guess.

  • update_ranges (array like, optional) – Updates the range of the generated config. It will use the frac_range option to set the parameter range.

  • frac_range (float, optional) –

    Set the half range of the parameter interval in the config. I.e.:

    min = (1 - frac_range)*guess
    max = (1 + frac_range)*guess
    

Returns

An updated version of the initial configuration object using the latest fit as the guess.

Return type

ConfigEmissionModel class

update_fitting(parameters, placeholder=None, inplace=False)

Return a copy of the fitting parameters in each model updated

This method takes into account fixed parameters and links, so that only parameters meant to be fitted are updated and linked accordingly. If inplace is True self.fitting_chain will be updated with the given parameters inplace.

Parameters
  • parameters (array like) – A list of the parameter arrays to update the current values in the parameters to be fitted.

  • placeholder (array like or None) – A list of parameters like parameters to be used as the initial values to the updated parameters instead self.latest_fit. Default value is None.

  • inplace (boolean) – Whether to update self.latest_fit and self.fitting_chain inplace or not. Defaults to False.

Returns

A list of the full parameter space (including fixed parameters) updated taking into account links.

Return type

array like

class pyFIT3D.modelling.gas.EmissionLinesLM(wavelength, flux, sigma_flux, config, n_MC, n_loops, scale_ini, plot)

Bases: pyFIT3D.modelling.gas.EmissionLines

FITfunc(x, config, wavelength, flux, sigma_flux)

Return the models evaluated at the current state of the parameters

This function deals with parameter linkings & fixed parameters

Parameters
  • x (array like) – The list of parameters to fit

  • all_parameters (array like) – The list of all parameters (fixed and not fixed)

  • is_fixed (array like) –

    A boolean mask for the fixed parameters so that:

    x = map(lambda ar: ar[is_fixed])
    

  • wavelength (array like) – The wavelength array in which the model will be evaluated

  • flux (array like) – The observed fluxes

  • sigma_flux (array like) – The standard deviation on the observed fluxes

Returns

The array of residuals resulting from evaluating the model in at the given parameters and the observed independent variable, i.e., the wavelength

Return type

array like

fit(check_stats=True)
output(filename_output, filename_spectra=None, filename_config=None, append_output=False, spec_id=None)

Outputs the final data from the LM method.

Parameters
  • filename_spectra (str) – The output filename for the input spectra, output model and output residuals.

  • filename_output (str or None) – The output filename for the final result.

  • filename_config (str or None) – The output filename for the final configuration of the emission lines system. This file will be formatted like the input config file.

  • append_output (bool) – When True, open filename_output for append. Otherwise, rewrite the file. It defaults to rewrite the file.

  • spec_id (int or tuple) – The spec coordinate id. For cubes is a tuple with (ix, iy).

output_to_screen(spec_id=None)
class pyFIT3D.modelling.gas.EmissionLinesRND(wavelength, flux, sigma_flux, config, n_MC, n_loops, scale_ini, plot, fine_search=None)

Bases: pyFIT3D.modelling.gas.EmissionLines

Fits a system of emission lines (EL) with a pseudo-random search of the input parameters.

wavelength

The observed wavelengths array.

Type

array like

flux

The observed spectrum array.

Type

array like

sigma_flux

The array with the standard deviation of the observed spectrum.

Type

array like

config

The class which configures the emission lines system.

Type

ConfigEmissionModel class

n_MC

The input number of Monte-Carlo realisations.

Type

int

n_loops

The maximum number of loops of the main search of the best EL parameters.

Type

int

scale_ini

Controls the size of the random search step.

Type

float

plot

Nice plot of the entire process.

Type

bool

Adds a final Monte-Carlo loop to self._MC_search with a fine search for the fitted parameters. The default value is False.

Type

bool

spec_id

The index coordinate(s) of the spectrum if looping from fits data. Default is None.

Type

int, tuple or None

fit :

The main function of the RND method. It makes a loop searching the fitting parameters narrowing the parameters range at each model that have a lower chi square.

output :

Prints the results to output files and to /dev/stdout.

See also

EmissionLines

fit(vel_fixed=False, sigma_fixed=False, randomize_flux=False, check_stats=True, redshift_flux_threshold=True, sigma_flux_threshold=False, oversize_chi=True)

The main function of the RND method. It makes a loop searching the fitting parameters narrowing the parameters range at each model that have a lower chi square.

output(filename_output, filename_spectra=None, filename_config=None, append_output=False, spec_id=None)

Outputs the final data from the rnd method.

Parameters
  • filename_output (str) – The output filename for the final result.

  • filename_spectra (str or None) – The output filename for the input spectra, output model and output residuals.

  • filename_config (str or None) – The output filename for the final configuration of the emission lines system. This file will be formatted like the input config file. If None does not print the config file.

  • append_output (bool) – When True, open filename_output for append. Otherwise, rewrite the file.

  • spec_id (int or tuple) – The spec coordinate id. For cubes is a tuple with (ix, iy).

output_to_screen(spec_id=None)
pyFIT3D.modelling.gas.output_config_final_fit(config, filename_output, chi_sq, parameters=None, e_parameters=None, append_output=False, spec_id=None)

Outputs the parameters and e_parameters formatted.

Parameters
  • config (ConfigEmissionModel class) – The ConfigEmissionModel instance.

  • filename_output (str) – The output filename for the final result.

  • parameters (array like, optional) – Parameters to be outputted. Defaults to config.guess.

  • e_parameters (array like, optional) – Sigma of the parameters to be outputted. Defaults to zeros_like(config.guess).

  • append_output (bool) – When True, open filename_output for append. Otherwise, rewrite the file.

  • spec_id (int or tuple) – The spec coordinate id. For cubes is a tuple with (ix, iy).

pyFIT3D.modelling.stellar module

class pyFIT3D.modelling.stellar.SSPModels(filename)

Bases: object

A helper class to deal with SSP models. It reads the models directly from FITS file.

wavelength

Wavelenghts at rest-frame for SSP models.

Type

array like

n_wave

Number of wavelengths.

Type

int

flux_models

Spectra of SSP models.

Type

array like

mass_to_light

Mass-to-light ratio of SSP models

Type

array like

normalization_wavelength

The normalization wavelength.

Type

int

parameters

The parameters of the SSP library. The keys are the parameters names. Example:

(...)
ssp = SSPModels(filename)
(...)
ages = ssp.parameters['AGE']
metallicities = ssp.parameters['MET']
Type

dict

get_par_average_from_coeffs()
get_models_extincted_obs_frame()
output_FITS()
get_model_from_coeffs(coeffs, wavelength, sigma, redshift, AV, sigma_inst, R_V=3.1, extlaw='CCM')

Shift and convolves SSP model fluxes (i.e. self.flux_models) to wavelengths wave_obs using sigma and sigma_inst. After this, applies dust extinction to the SSPs following the extinction law extlaw with AV attenuance. At the end, returns the SSP model spectra using coeffs.

Parameters
  • coeffs (array like) – Coefficients of each SSP model.

  • wavelength (array like) – Wavelenghts at observed frame.

  • sigma (float) – Velocity dispersion (i.e. sigma) at observed frame.

  • redshift (float) – Redshift of the Observed frame.

  • AV (float or array like) – Dust extinction in mag. TODO: If AV is an array, will create an (n_AV, n_wave) array of dust spectra.

  • sigma_inst (float or None) – Sigma instrumental.

  • R_V (float, optional) – Selective extinction parameter (roughly “slope”). Default value 3.1.

  • extlaw (str {'CCM', 'CAL'}, optional) – Which extinction function to use. CCM will call Cardelli_extlaw. CAL will call Calzetti_extlaw. Default value is CCM.

Returns

SSP model spectrum created by coeffs.

Return type

array like

get_models_extincted_obs_frame(wavelength, sigma, redshift, AV, sigma_inst, R_V=3.1, extlaw='CCM', coeffs=None)

Shift and convolves SSP model fluxes (i.e. self.flux_models) to wavelength using the observed kinetic parameters redshift and sigma (velocity dispersion). The convolution considers the instrumental dispersion of the observed data if sigma_inst is not None. After the shift and convolution of the models, it applies dust extinction to the SSPs following the extinction law extlaw with AV attenuance. It uses the vector population coeffs in order to accelerate the process only performing the calculation over the models where coeffs != 0.

Parameters
  • wavelength (array like) – Wavelenghts at observed frame.

  • sigma (float) – Velocity dispersion (i.e. sigma) at observed frame.

  • redshift (float) – Redshift of the Observed frame.

  • AV (float or array like) – Dust extinction in mag. TODO: If AV is an array, will create an (n_AV, n_wave) array of dust spectra.

  • sigma_inst (float or None) – Sigma instrumental.

  • R_V (float, optional) – Selective extinction parameter (roughly “slope”). Default value 3.1.

  • extlaw (str {'CCM', 'CAL'}, optional) – Which extinction function to use. CCM will call Cardelli_extlaw. CAL will call Calzetti_extlaw. Default value is CCM.

  • coeffs (array like) – Coefficients of each SSP model.

Returns

SSP spectra shift and convolved at the observed frame extincted by dust. Has the same dimension as (SSPModels.n_models, wavelenght.size).

Return type

array like

get_par_average_from_coeffs(coeffs, par_name, log_avg=False)

Return the value from the age and metallicity weighted by light and mass from the age and metallicity of each model weighted by the coeffs.

Parameters
  • coeffs (array like) – Coefficients of each SSP model.

  • par_name (str) – Name of parameter to be averaged.

Returns

  • float – Parameter average.

  • float – Parameter average weighted by mass.

output_FITS(label='new', path='.')

Outputs the SSP library FITS file.

Parameters
  • label (str, optional) –

    Labels the output filename as:

    output_filename = f'SSPModels-{label}.fits.gz'
    

    Defaults to ‘new’.

  • path (str, optional) – The output path. Defaults to the local running directory.

Returns

The output HDU list class.

Return type

astropy.io.fits.HDUList

class pyFIT3D.modelling.stellar.SSPModels_AVgrid(filename, AV_min=0, AV_max=2, n=5)

Bases: pyFIT3D.modelling.stellar.SSPModels

class pyFIT3D.modelling.stellar.StPopSynt(config, wavelength, flux, eflux, ssp_file, out_file, ssp_nl_fit_file=None, sigma_inst=None, min=None, max=None, w_min=None, w_max=None, nl_w_min=None, nl_w_max=None, elines_mask_file=None, mask_list=None, R_V=3.1, extlaw='CCM', spec_id=None, guided_errors=None, plot=None, verbose=False, ratio_master=None, fit_gas=True)

Bases: object

This class is created to fit an observed spectrum using a decomposition of the underlying stellar population through a set of simple stellar population (SSP) models.

The process starts with the search for the redshift, sigma and dust extinction (AV) of the observed spectrum. After this kinemactic evaluation, the gas spectrum is derived using a set of emission line models preseted in config and fitted by pyFIT3D.common.gas_tools.fit_elines_main.

Then, the gas spectrum is subtracted from the original spectrum and the stellar population synthesis is performed, with the calculation of the coefficients of a set of SSP models. The search for the best coefficients is made through a Monte-Carlo perturbation loop of the observed spectrum.

This class uses :class: SSPModels to manage the simple stellar population models and pyFIT3D.common.gas_tools.fit_elines_main for the emission line fitting.

An example of utilization of :class: StPopSynt can be found in: pyFIT3D.common.auto_ssp_tools.auto_ssp_elines_single_main

config

The class which configures the whole SSP fit process.

Type

ConfigAutoSSP class

spectra

A dictionary with the observed spectral information (wavelengths, fluxes and the standard deviation of the fluxes, with the masks).

Type

dict

R_V

Selective extinction parameter (roughly “slope”). Default value 3.1.

Type

float, optional

extlaw

Which extinction function to use. CCM will call Cardelli_extlaw. CAL will call Calzetti_extlaw. Default value is CCM.

Type

str {‘CCM’, ‘CAL’}, optional

verbose

If True produces a nice text output.

Type

bools, optional

filename

Path to the SSP models fits file.

Type

str

filename_nl_fit

Path to the SSP models fits file used in the non-linear round of fit. The non-linear procedure search for the kinematics parameters (redshift and sigma) and the dust extinction (AV). If None, self.ssp_nl_fit is equal to self.ssp. Default value is None.

Type

str or None

n_loops

Counts the number of loops in the whole fit process.

Type

int

n_loops_nl_fit

Counts the number of loops in the non-linear fit process.

Type

int

n_loops_redshift

Counts the number of loops in the redshift search process.

Type

int

n_loops_sigma

Counts the number of loops in the sigma search process.

Type

int

n_loops_AV

Counts the number of loops in the A_V search process.

Type

int

best_coeffs_redshift

Best SSP Models coefficients of redshift fit.

best_chi_sq_redshift

Best Chi squared of redshift fit.

best_redshift

Best redshift.

e_redshift

Error in best redshift.

redshift_chain

Chain of inspected redshifts.

best_coeffs_sigma

Best SSP Models coefficients of sigma fit.

best_chi_sq_sigma

Best Chi squared of sigma fit.

best_sigma

Best sigma.

e_sigma

Error in best sigma.

sigma_chain

Chain of inspected sigmas.

best_coeffs_AV

Best SSP Models coefficients of AV fit.

best_chi_sq_AV

Best Chi squared of AV fit.

best_AV

Best AV.

e_AV

Error in best AV.

AV_chain

Chain of inspected AVs.

coeffs_ssp_chain

All coefficients of the ssp search process.

chi_sq_ssp_chain

All Chi squared of the ssp search process.

non_linear_fit :
gas_fit :
ssp_fit :
fit_WLS_invmat :
fit_WLS_invmat_MC :
output_coeffs_MC :
output_fits :
output :
fit_WLS_invmat :

Fits an observed spectrum a with a linear combination of SSP models using Weighted Least Squares (WLS) through matrix inversion. This process consider measured errors of the observed spectrum.

get_last_redshift :
get_last_chi_sq_redshift :
get_last_coeffs_redshift :
update_redshift_params :
plot_last :
get_last_sigma :
get_last_chi_sq_sigma :
get_last_coeffs_sigma :
update_sigma_params :
get_last_AV :
get_last_chi_sq_AV :
get_last_coeffs_AV :
update_AV_params :
get_last_chi_sq_ssp :
get_last_coeffs_ssp :
update_ssp_params :
EL_fit_plot_close()
EL_fit_plot_final_spec()
EL_fit_plot_init(model_min, n_lines_per_row=6)
EL_fit_plot_system(els, els_config, model_min)
calc_SN_norm_window(reference_wavelength=None, half_range=45)

Calculates the signal-to-noise inside the wavelength window defined by (reference_wavelength +/- half_range).

Parameters
  • reference_wavelength (float, optional) – Configures the central wavelength of the normalization window. If None uses the normalization wavelength of the stellar population models (SSPModels.wavenorm).

  • half_range (int, optional) – Sets the normalization window size in Angstroms. It defaults to 45 AA.

sel_SN_norm_window

An array of booleans selecting the wavelength range of the normalization window. The reference wavelength array is self.spectra[‘raw_wave’].

Type

array like

SN_norm_window

The signal-to-noise at the normalization window.

Type

float

correct_elines_mask(redshift, window_size=None, update_nl_range=False)
fit_AV(ssp=None)

Fits the dust extinction parameter (AV) of an observed spectra. This method updates the AV chain.

Parameters

ssp (SSPModels class or None) – The class with ssp models used during the fit process. Defaults to self.models_nl_fit.

fit_WLS_invmat(ssp, sigma_inst=None, sigma=None, redshift=None, AV=None, smooth_cont=False, sel_wavelengths=None)

A wrapper to _fit_WLS_invmat setting up the spectra, masks, and extinction law. If smooth_cont is True smooths the continuum modeled by sigma_mean, where:

sigma_mean^2 = sigma_inst^2 + (5000*(sigma/c))^2

Obs.: This method calls ssp.apply_dust_to_flux_models_obsframe() which will set the models spectra to the input AV ssp.flux_models_obsframe_dust.

Parameters
  • ssp (SSPModels class) – The class with ssp models used during the fit process.

  • sigma_inst (float or None) – Sigma instrumental. If the sigma_inst is None, is used self.sigma_inst.

  • redshift (float or None) – Redshift of the Observed frame. If the redshift is None, is used self.redshift.best[“value”].

  • sigma (float or None) – Velocity dispersion (i.e. sigma) at observed frame. If the sigma is not set None, is used self.sigma.best[“value”].

  • AV (float, array like or None) – Dust extinction in mag. If AV is an array, will create an (n_AV, n_wave) array of dust spectra. If None, is used self.AV.best[“value”]

  • smooth_cont (bool, optional) – If True smooths the continuum modeled by sigma_mean, where sigma_mean^2 = sigma_inst^2 + (5000*(sigma/c))^2

  • sel_wavelengths (array like or None) – The selection of valid wavelengths. If None uses self.spectra[‘sel_wl’] as the wavelength selection.

Returns

  • array like – Coefficients of the WLS fit.

  • float – Chi square of the fit.

fit_WLS_invmat_MC(ssp, n_MC=20, sigma_inst=None, sigma=None, redshift=None, AV=None, sel_wavelengths=None)

This method wraps a Monte-Carlo realisation of _fit_WLS_invmat() n_MC times used during the SSP fit of an observed spectra. This process consider measured errors of the observed spectrum. This method models the observed spectrum without emission lines, i.e. it uses self.spectra[‘raw_flux_no_gas’] as the input spectra).

Obs.: This method calls ssp.apply_dust_to_flux_models_obsframe() which will set the models spectra to the input AV ssp.flux_models_obsframe_dust.

Parameters
  • ssp (SSPModels class) – The class with ssp models used during the fit process.True

  • n_MC (int) – Number of Monte-Carlo loops

  • sigma_inst (float or None) – Sigma instrumental. If the sigma_inst is None, is used self.sigma_inst.

  • redshift (float or None) – Redshift of the Observed frame. If the redshift is None, is used self.redshift.best[“value”].

  • sigma (float or None) – Velocity dispersion (i.e. sigma) at observed frame. If the sigma is not set None, is used self.sigma.best[“value”].

  • AV (float, array like or None) – Dust extinction in mag. If AV is an array, will create an (n_AV, n_wave) array of dust spectra. If None, is used self.AV.best[“value”]

  • sel_wavelengths (array like or None) – The selection of valid wavelengths. If None uses self.spectra[‘sel_wl’] as the wavelength selection.

Returns

  • array like – Coefficients from the MC realisation.

  • array like – Chi squareds from the MC realisation.

  • array like – Models from the MC realisation.

fit_redshift(ssp=None, correct_wl_ranges=False)
fit_sigma(guided=False, ssp=None, calc_coeffs=True, medres_merit=False)
gas_fit(ratio=True, ratio_range=None, ratio_std=None, reference_wavelength=None, half_range=45, EL_fit_half_range_sysvel=None, EL_fit_correct_wl_range=2, EL_fit_guide_vel=True)

Prepares the observed spectra in order to fit systems of emission lines to the residual spectra.

Parameters
  • ratio (bool) –

  • ratio_range (float, optional) –

  • ratio_std (float, optional) –

  • reference_wavelength (float, optional) – The central wavelength of the normalization window at rest-frame. If None uses the normalization wavelength of the stellar population models (SSPModels.wavenorm).

  • half_range (int, optional) – Sets the normalization window size in Angstroms. It defaults to 45 AA.

  • EL_fit_half_range_sysvel (float, optional) –

    Defines the range of the velocity of the emission lines. Defaults to __selected_half_range_sysvel_auto_ssp__. .. code-block:: python

    systemic_velocity = light_velocity*redshift vel_range = [

    systemic_velocity - half_range_sysvel, systemic_velocity + half_range_sysvel,

    ]

  • EL_fit_correct_wl_range (int {0, 1, 2}, optional) –

    Determine the wavelength range redshift correction of the emission line system search. The accepted values are:

    0 : Do not correct wavelength range.
    1 : correct wavelength range using `self.best_redshift`
    2 : correct the range based on the system of emission lines.
        wl_range = [
            (eml_central_wave_min - (half_range_sysvel/10)),
            (eml_central_wave_max + (half_range_sysvel/10))
        ]
    

    Defaults to 0.

  • EL_fit_guide_vel (bool, optional) – Guide velocity guess of the fit around self.best_redshift. Defaults to True.

spectra['raw_flux_no_gas']

The raw observed spectrum without the model of the emission lines.

Type

array like

get_best_model_from_coeffs(ssp, coeffs)
get_last_chi_sq_ssp()
get_last_coeffs_ssp()
nl_fit_plot_close()
nl_fit_plot_init()
nl_fit_plot_last(sel, parameter)
non_linear_fit(guide_sigma=False, fit_sigma_rnd=True, sigma_rnd_medres_merit=False, correct_wl_ranges=False)

Do the non linear fit in order to find the kinematics parameters and the dust extinction. This procedure uses the set of SSP models, self.models_nl_fit. At the end will set the first entry to the ssp fit chain with the coefficients for the best model after the non-linear fit.

output(filename, write_header=True, block_plot=True)

Summaries the run in a csv file.

Parameters

filename (str) – Output filename.

output_coeffs_MC(filename, write_header=True)

Outputs the SSP coefficients table to the screen and to the output file filename.

Parameters

filename (str) – The output filename to the coefficients table.

output_coeffs_MC_to_screen()
output_fits(filename, seed=None)

Writes the FITS file with the output spectra (original, model, residual and joint).

Parameters

filename (str) – Output FITS filename.

output_gas_emission(filename, spec_id=None, append=True)
output_single_ssp(filename)
output_to_screen(block_plot=True)
plot_final_summary(cmap='Blues', vmin=0.1, vmax=50, sigma_unit='km/s', percent=True, Zsun=0.019, block_plot=False)
redshift_correct_masks(redshift, eline_half_range=None, correct_wl_ranges=False)
resume_results()
ssp_fit(n_MC=20)

Generates minimal ssp model through a Monte-Carlo search of the coefficients. I.e. go through fit_WLS_invmat() n_MC times (using fit.WLS_invmat_MC).

Parameters

n_MC (int) – Number of Monte-Carlo loops. Default value is 20.

See also

fit_WLS_invmat

ssp_init()
ssp_single_fit()
update_ssp_parameters(coeffs, chi_sq)

Save the SSP fit parameters chain of tries.

Parameters
  • coeffs (array like) – The last calculed coefficients for the SSP models.

  • chi_sq (float) – The chi squared of the fit.

class pyFIT3D.modelling.stellar.nl_fit_par(name, guess, delta, max_val, min_val, broad_fit_update_parameter_function=None, broad_fit_update_delta_function=None, narrow_fit_update_parameter_function=None, narrow_fit_update_delta_function=None, update_min_par_function=None, update_max_par_function=None, update_delta_function=None, verbose=0)

Bases: object

get_broad_chain()
get_broad_chain_key(key)
get_broad_hyperbolic_fit()
get_narrow_chain()
get_narrow_chain_key(key)
get_narrow_hyperbolic_fit()
loop(min_par, max_par, delta_par, update_parameter_function=None, update_delta_function=None)

Creates the array of tries of a parameter following an update_parameter_function and a update_delta_function.

Parameters
  • min_par (float) –

  • max_par (float) –

  • delta_par (float) –

  • update_parameter_function (function or None) –

  • update_delta_function (function or None) –

Returns

Returns an iterator that loops a parameter from min_par to max_par values following an update_parameter_function.

Return type

iterator

loop_broad()
loop_narrow()
report_best(level=1, print_coeffs=False)
report_broad_last(level=1, print_coeffs=False)
report_narrow_last(level=1, print_coeffs=False)
set_best(value, coeffs, chi_sq, model)
set_chain(value, coeffs, chi_sq, model)