Formaldehyde Models¶
The Formaldehyde model is based on the Hyperfine Model. There is also a mm-line fitting module that uses preconstructed model grids to fit temperature from multiple transitions.
There is a formaldehyde fitter wrapper; see Wrappers.
This is a formaldehyde 1_11-1_10 / 2_12-2_11 fitter. It includes hyperfine components of the formaldehyde lines and has both LTE and RADEX LVG based models
Module API¶
- pyspeckit.spectrum.models.formaldehyde.formaldehyde(xarr, amp=1.0, xoff_v=0.0, width=1.0, return_hyperfine_components=False, texscale=0.01, tau=0.01, **kwargs)[source] [github] [bitbucket]¶
Generate a model Formaldehyde spectrum based on simple gaussian parameters
the “amplitude” is an essentially arbitrary parameter; we therefore define it to be Tex given tau=0.01 when passing to the fitter The final spectrum is then rescaled to that value
- class pyspeckit.spectrum.models.formaldehyde.formaldehyde_model(modelfunc, npars, shortvarnames=('A', '\\Delta x', '\\sigma'), fitunit=None, centroid_par=None, fwhm_func=None, fwhm_pars=None, integral_func=None, use_lmfit=False, guess_types=('amplitude', 'center', 'width'), **kwargs)[source] [github] [bitbucket]¶
Spectral Model Initialization
Create a Spectral Model class for data fitting
- Parameters
modelfunc : function
the model function to be fitted. Should take an X-axis (spectroscopic axis) as an input followed by input parameters. Returns an array with the same shape as the input X-axis
npars : int
number of parameters required by the model
use_lmfit: bool :
Use lmfit instead of mpfit to do the fitting
parnames : list (optional)
a list or tuple of the parameter names
parvalues : list (optional)
the initial guesses for the input parameters (defaults to ZEROS)
parlimits : list (optional)
the upper/lower limits for each variable (defaults to ZEROS)
parfixed : list (optional)
Can declare any variables to be fixed (defaults to ZEROS)
parerror : list (optional)
technically an output parameter. Specifying it here will have no effect. (defaults to ZEROS)
partied : list (optional)
not the past tense of party. Can declare, via text, that some parameters are tied to each other. Defaults to zeros like the others, but it’s not clear if that’s a sensible default
fitunit : str (optional)
convert X-axis to these units before passing to model
parsteps : list (optional)
minimum step size for each paremeter (defaults to ZEROS)
npeaks : list (optional)
default number of peaks to assume when fitting (can be overridden)
shortvarnames : list (optional)
TeX names of the variables to use when annotating
amplitude_types : tuple
A tuple listing the types of the different parameters when guessing. The valid values are ‘amplitude’, ‘width’, and ‘center’. These are handled by parse_3par_guesses, which translate these into input guess lists for the fitter. For a “standard” 3-parameter Gaussian fitter, nothing changes, but for other models that have more than 3 parameters, some translation is needed.
- Returns
A tuple containing (model best-fit parameters, the model, parameter :
errors, chi^2 value) :
- formaldehyde_integral(modelpars, linename='oneone')[source] [github] [bitbucket]¶
Return the integral of the individual components (ignoring height)
- pyspeckit.spectrum.models.formaldehyde.formaldehyde_pyradex(xarr, density=4, column=13, temperature=20, xoff_v=0.0, opr=1.0, width=1.0, tbackground=2.73, grid_vwidth=1.0, debug=False, verbose=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
grid_vwidth is the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230)
- pyspeckit.spectrum.models.formaldehyde.formaldehyde_radex(xarr, density=4, column=13, xoff_v=0.0, width=1.0, grid_vwidth=1.0, grid_vwidth_scale=False, texgrid=None, taugrid=None, hdr=None, path_to_texgrid='', path_to_taugrid='', temperature_gridnumber=3, debug=False, verbose=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
grid_vwidth is the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230) grid_vwidth_scale is True or False: False for LVG, True for Sphere
- pyspeckit.spectrum.models.formaldehyde.formaldehyde_radex_orthopara_temp(xarr, density=4, column=13, orthopara=1.0, temperature=15.0, xoff_v=0.0, width=1.0, Tbackground1=2.73, Tbackground2=2.73, grid_vwidth=1.0, grid_vwidth_scale=False, texgrid=None, taugrid=None, hdr=None, path_to_texgrid='', path_to_taugrid='', debug=False, verbose=False, getpars=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
grid_vwidth is the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230) grid_vwidth_scale is True or False: False for LVG, True for Sphere
- pyspeckit.spectrum.models.formaldehyde.formaldehyde_radex_tau(xarr, density=4, column=13, xoff_v=0.0, width=1.0, grid_vwidth=1.0, grid_vwidth_scale=False, taugrid=None, hdr=None, path_to_taugrid='', temperature_gridnumber=3, debug=False, verbose=False, return_hyperfine_components=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
uses hyperfine components
assumes tau varies but tex does not!
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
grid_vwidth is the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230) grid_vwidth_scale is True or False: False for LVG, True for Sphere
This is a formaldehyde 3_03-2_02 / 3_22-221 and 3_03-2_02/3_21-2_20 fitter. It is based entirely on RADEX models.
Module API¶
- pyspeckit.spectrum.models.formaldehyde_mm.build_despotic_grids(gridfile='ph2co_grid_despotic.fits', ph2coAbund=1e-08, nDens=21, logDensLower=2.0, logDensUpper=6.0, nCol=21, logColLower=11.0, logColUpper=15.0, nTemp=51, Tlower=10.0, Tupper=300.0, nDv=5, DvLower=1.0, DvUpper=5.0)[source] [github] [bitbucket]¶
Generates grids of p-H2CO line intensities using Despotic. Outputs a astropy Table.
- Parameters
gridfile : string
Name of grid file to output.
ph2coAbund : float
Fractional abundance of p-H2CO
nDens : int
Number of grid points in the volume density
logDensLower : float
log of volume density at lower bound of grid (log(n/cm**-3))
logDensUpper : float
log of volume density at upper bound of grid (log(n/cm**-3))
nCol : int
Number of grid points in the column density
logColLower : float
log of column density of p-H2CO at lower bound of grid (log(N/cm**-2))
logColUpper : float
log of column density of p-H2CO at upper bound of grid (log(N/cm**-2))
nTemp : int
Number of grid points in the temperature grid
Tower : float
temperature at lower bound of grid (K)
Tupper : float
temperature at upper bound of grid (K)
nDv : int
Number of grid points in the line width
DvLower : float
line width (non-thermal) at lower bound of grid (km/s)
DvUpper : float
line width (non-thermal) at upper bound of grid (km/s)
- pyspeckit.spectrum.models.formaldehyde_mm.formaldehyde_mm(xarr, amp=1.0, xoff_v=0.0, width=1.0, return_components=False)[source] [github] [bitbucket]¶
Generate a model Formaldehyde spectrum based on simple gaussian parameters
the “amplitude” is an essentially arbitrary parameter; we therefore define it to be Tex given tau=0.01 when passing to the fitter The final spectrum is then rescaled to that value
The components are independent, but with offsets set by frequency… in principle.
- pyspeckit.spectrum.models.formaldehyde_mm.formaldehyde_mm_despotic(xarr, temperature=25, column=13, density=4, xoff_v=0.0, width=1.0, grid_vwidth=1.0, h2co_303_202=None, h2co_322_221=None, h2co_321_220=None, debug=False, verbose=False, **kwargs)[source] [github] [bitbucket]¶
Fitter to p-H2CO using despotic grids. Requires building grids and passing in functions for interpolating the h2co transition optical depth and excitation temperatures.
- pyspeckit.spectrum.models.formaldehyde_mm.formaldehyde_mm_despotic_functions(gridtable)[source] [github] [bitbucket]¶
This builds interpolation functions for use in fitting.
- Parameters
gridtable : str
Name of grid in astropy table
- Returns
h2co_303_202, h2co_322_221, h2co_321_220 : function
Functions that return the excitation temperature and optical depth given input density, temperature, column density and line width.
- class pyspeckit.spectrum.models.formaldehyde_mm.formaldehyde_mm_model(modelfunc, npars, shortvarnames=('A', '\\Delta x', '\\sigma'), fitunit=None, centroid_par=None, fwhm_func=None, fwhm_pars=None, integral_func=None, use_lmfit=False, guess_types=('amplitude', 'center', 'width'), **kwargs)[source] [github] [bitbucket]¶
Spectral Model Initialization
Create a Spectral Model class for data fitting
- Parameters
modelfunc : function
the model function to be fitted. Should take an X-axis (spectroscopic axis) as an input followed by input parameters. Returns an array with the same shape as the input X-axis
npars : int
number of parameters required by the model
use_lmfit: bool :
Use lmfit instead of mpfit to do the fitting
parnames : list (optional)
a list or tuple of the parameter names
parvalues : list (optional)
the initial guesses for the input parameters (defaults to ZEROS)
parlimits : list (optional)
the upper/lower limits for each variable (defaults to ZEROS)
parfixed : list (optional)
Can declare any variables to be fixed (defaults to ZEROS)
parerror : list (optional)
technically an output parameter. Specifying it here will have no effect. (defaults to ZEROS)
partied : list (optional)
not the past tense of party. Can declare, via text, that some parameters are tied to each other. Defaults to zeros like the others, but it’s not clear if that’s a sensible default
fitunit : str (optional)
convert X-axis to these units before passing to model
parsteps : list (optional)
minimum step size for each paremeter (defaults to ZEROS)
npeaks : list (optional)
default number of peaks to assume when fitting (can be overridden)
shortvarnames : list (optional)
TeX names of the variables to use when annotating
amplitude_types : tuple
A tuple listing the types of the different parameters when guessing. The valid values are ‘amplitude’, ‘width’, and ‘center’. These are handled by parse_3par_guesses, which translate these into input guess lists for the fitter. For a “standard” 3-parameter Gaussian fitter, nothing changes, but for other models that have more than 3 parameters, some translation is needed.
- Returns
A tuple containing (model best-fit parameters, the model, parameter :
errors, chi^2 value) :
- pyspeckit.spectrum.models.formaldehyde_mm.formaldehyde_mm_radex(xarr, temperature=25, column=13, density=4, xoff_v=0.0, width=1.0, grid_vwidth=1.0, texgrid=None, taugrid=None, hdr=None, path_to_texgrid='', path_to_taugrid='', debug=False, verbose=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
- Parameters
grid_vwidth : float
the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230)
density : float
Density!
Erik R’s “fork” of the mm fitter:
This is a formaldehyde 3_03-2_02 / 3_22-221 and 3_03-2_02/3_21-2_20 fitter. It is based entirely on RADEX models.
This is the EWR fork of the fitter in pyspeckit.
Module API¶
- pyspeckit.spectrum.models.h2co_mm.formaldehyde_mm(xarr, amp=1.0, xoff_v=0.0, width=1.0, return_components=False)[source] [github] [bitbucket]¶
Generate a model Formaldehyde spectrum based on simple gaussian parameters
the “amplitude” is an essentially arbitrary parameter; we therefore define it to be Tex given tau=0.01 when passing to the fitter The final spectrum is then rescaled to that value
The components are independent, but with offsets set by frequency… in principle.
- class pyspeckit.spectrum.models.h2co_mm.formaldehyde_mm_model(modelfunc, npars, shortvarnames=('A', '\\Delta x', '\\sigma'), fitunit=None, centroid_par=None, fwhm_func=None, fwhm_pars=None, integral_func=None, use_lmfit=False, guess_types=('amplitude', 'center', 'width'), **kwargs)[source] [github] [bitbucket]¶
Spectral Model Initialization
Create a Spectral Model class for data fitting
- Parameters
modelfunc : function
the model function to be fitted. Should take an X-axis (spectroscopic axis) as an input followed by input parameters. Returns an array with the same shape as the input X-axis
npars : int
number of parameters required by the model
use_lmfit: bool :
Use lmfit instead of mpfit to do the fitting
parnames : list (optional)
a list or tuple of the parameter names
parvalues : list (optional)
the initial guesses for the input parameters (defaults to ZEROS)
parlimits : list (optional)
the upper/lower limits for each variable (defaults to ZEROS)
parfixed : list (optional)
Can declare any variables to be fixed (defaults to ZEROS)
parerror : list (optional)
technically an output parameter. Specifying it here will have no effect. (defaults to ZEROS)
partied : list (optional)
not the past tense of party. Can declare, via text, that some parameters are tied to each other. Defaults to zeros like the others, but it’s not clear if that’s a sensible default
fitunit : str (optional)
convert X-axis to these units before passing to model
parsteps : list (optional)
minimum step size for each paremeter (defaults to ZEROS)
npeaks : list (optional)
default number of peaks to assume when fitting (can be overridden)
shortvarnames : list (optional)
TeX names of the variables to use when annotating
amplitude_types : tuple
A tuple listing the types of the different parameters when guessing. The valid values are ‘amplitude’, ‘width’, and ‘center’. These are handled by parse_3par_guesses, which translate these into input guess lists for the fitter. For a “standard” 3-parameter Gaussian fitter, nothing changes, but for other models that have more than 3 parameters, some translation is needed.
- Returns
A tuple containing (model best-fit parameters, the model, parameter :
errors, chi^2 value) :
- pyspeckit.spectrum.models.h2co_mm.h2co_mm_radex(xarr, Temperature=25, logColumn=13, logDensity=4, xoff_v=0.0, width=1.0, grid_vwidth=1.0, gridbundle=None, debug=False, verbose=False, **kwargs)[source] [github] [bitbucket]¶
Use a grid of RADEX-computed models to make a model line spectrum
The RADEX models have to be available somewhere. OR they can be passed as arrays. If as arrays, the form should be: texgrid = ((minfreq1,maxfreq1,texgrid1),(minfreq2,maxfreq2,texgrid2))
xarr must be a SpectroscopicAxis instance xoff_v, width are both in km/s
- Parameters
grid_vwidth : float
the velocity assumed when computing the grid in km/s this is important because tau = modeltau / width (see, e.g., Draine 2011 textbook pgs 219-230)
density : float
Density!