Wrappers

These are wrappers to simplify some of the more complicated (and even some of the simpler) functions in PySpecKit

Cube Fitting

Complicated code for fitting of a whole data cube, pixel-by-pixel

pyspeckit.wrappers.cube_fit.cube_fit(cubefilename, outfilename, errfilename=None, scale_keyword=None, vheight=False, verbose=False, signal_cut=3, verbose_level=2, overwrite=True, **kwargs)[source] [github] [bitbucket]

Light-weight wrapper for cube fitting

Takes a cube and error map (error will be computed naively if not given) and computes moments then fits for each spectrum in the cube. It then saves the fitted parameters to a reasonably descriptive output file whose header will look like

PLANE1  = 'amplitude'
PLANE2  = 'velocity'
PLANE3  = 'sigma'
PLANE4  = 'err_amplitude'
PLANE5  = 'err_velocity'
PLANE6  = 'err_sigma'
PLANE7  = 'integral'
PLANE8  = 'integral_error'
CDELT3  = 1
CTYPE3  = 'FITPAR'
CRVAL3  = 0
CRPIX3  = 1
Parameters

errfilename: [ None | string name of .fits file ] :

A two-dimensional error map to use for computing signal-to-noise cuts

scale_keyword: [ None | Char ] :

Keyword to pass to the data cube loader - multiplies cube by the number indexed by this header kwarg if it exists. e.g., if your cube is in T_A units and you want T_A*

vheight: [ bool ] :

Is there a background to be fit? Used in moment computation

verbose: [ bool ] :

verbose_level: [ int ] :

How loud will the fitting procedure be? Passed to momenteach and fiteach

signal_cut: [ float ] :

Signal-to-Noise ratio minimum. Spectra with a peak below this S/N ratio will not be fit and will be left blank in the output fit parameter cube

overwrite: [ bool ] :

Overwrite parameter .fits cube if it exists?

`kwargs` are passed to :class:`pyspeckit.Spectrum.specfit` :

pyspeckit.wrappers.fit_gaussians_to_simple_spectra.fit_gaussians_to_simple_spectra(filename, unit='km/s', doplot=True, baseline=True, plotresiduals=False, figuresavename=None, croprange=None, savename=None, **kwargs)[source] [github] [bitbucket]

As stated in the name title, will fit Gaussians to simple spectra!

kwargs will be passed to specfit

figuresavename [ None | string ]

After fitting, save the figure to this filename if specified

croprange [ list of 2 floats ]

Crop the spectrum to (min,max) in the specified units

savename [ None | string ]

After fitting, save the spectrum to this filename

Note that this wrapper can be used from the command line:

python fit_gaussians_to_simple_spectra.py spectrum.fits

NH3 fitter wrapper

Wrapper to fit ammonia spectra. Generates a reasonable guess at the position and velocity using a gaussian fit

Example use:

import pyspeckit
sp11 = pyspeckit.Spectrum('spec.nh3_11.dat', errorcol=999)
sp22 = pyspeckit.Spectrum('spec.nh3_22.dat', errorcol=999)
sp33 = pyspeckit.Spectrum('spec.nh3_33.dat', errorcol=999)
sp11.xarr.refX = pyspeckit.spectrum.models.ammonia.freq_dict['oneone']
sp22.xarr.refX = pyspeckit.spectrum.models.ammonia.freq_dict['twotwo']
sp33.xarr.refX = pyspeckit.spectrum.models.ammonia.freq_dict['threethree']
input_dict={'oneone':sp11, 'twotwo':sp22, 'threethree':sp33}
spf = pyspeckit.wrappers.fitnh3.fitnh3tkin(input_dict)

Note that if you want to use the plotter wrapper with cubes, you need to do something like the following, where the plot_special method of the stacked cubes object is set to the plotter_override function defined in the fitnh3_wrapper code:

cubes.plot_special = pyspeckit.wrappers.fitnh3.plotter_override
cubes.plot_special_kwargs = {'fignum':3, 'vrange':[55,135]}
cubes.plot_spectrum(160,99)
pyspeckit.wrappers.fitnh3.BigSpectrum_to_NH3dict(sp, vrange=None)[source] [github] [bitbucket]

A rather complicated way to make the spdicts above given a spectrum…

pyspeckit.wrappers.fitnh3.fitnh3(spectrum, vrange=[-100, 100], vrangeunit='km/s', quiet=False, Tex=20, trot=15, column=1000000000000000.0, fortho=1.0, tau=None, Tkin=None, fittype='ammonia', spec_convert_kwargs={})[source] [github] [bitbucket]
pyspeckit.wrappers.fitnh3.fitnh3tkin(input_dict, dobaseline=True, baselinekwargs={}, crop=False, cropunit=None, guessline='twotwo', tex=15, trot=20, column=15.0, fortho=0.66, tau=None, thin=False, quiet=False, doplot=True, fignum=1, guessfignum=2, smooth=False, scale_keyword=None, rebase=False, tkin=None, npeaks=1, guesses=None, fittype='ammonia', guess_error=True, plotter_wrapper_kwargs={}, **kwargs)[source] [github] [bitbucket]

Given a dictionary of filenames and lines, fit them together e.g. {‘oneone’:’G000.000+00.000_nh3_11.fits’}

Parameters

input_dict : dict

A dictionary in which the keys are the ammonia line names (e.g., ‘oneone’, ‘twotwo’, etc) and the values are either Spectrum objects or filenames of spectra

dobaseline : bool

Fit and subtract a baseline prior to fitting the model? Keyword arguments to pyspeckit.spectrum.Spectrum.baseline are specified in baselinekwargs.

baselinekwargs : dict

The keyword arguments for the baseline

crop : bool or tuple

A range of values to crop the spectrum to. The units are specified by cropunit; the default None will use pixels. If False, no cropping will be performed.

cropunit : None or astropy unit

The unit for the crop parameter

guess_error : bool

Use the guess line to estimate the error in all spectra?

plotter_wrapper_kwargs : dict

Keyword arguments to pass to the plotter

fittype: ‘ammonia’ or ‘cold_ammonia’ :

The fitter model to use. This is overridden if tau is specified, in which case one of the ammonia_tau models is used (see source code)

pyspeckit.wrappers.fitnh3.make_axdict(splist, spdict)[source] [github] [bitbucket]
pyspeckit.wrappers.fitnh3.plot_nh3(spdict, spectra, fignum=1, show_components=False, residfignum=None, show_hyperfine_components=True, annotate=True, axdict=None, figure=None, **plotkwargs)[source] [github] [bitbucket]

Plot the results from a multi-nh3 fit

spdict needs to be dictionary with form:

‘oneone’: spectrum, ‘twotwo’: spectrum, etc.

pyspeckit.wrappers.fitnh3.plotter_override(sp, vrange=None, **kwargs)[source] [github] [bitbucket]

Do plot_nh3 with syntax similar to plotter()

N2H+ fitter wrapper

Wrapper to fit N2H+ using RADEX models. This is meant to be used from the command line, e.g.

python n2hp_wrapper.py file.fits

and therefore has no independently defined functions.

pyspeckit.wrappers.n2hp_wrapper.make_n2hp_fitter(path_to_radex='/Users/adam/work/n2hp/', fileprefix='1-2_T=5to55_lvg')[source] [github] [bitbucket]

Create a n2hp fitter using RADEX data cubes. The following files must exist:

path_to_radex+fileprefix+'_tex1.fits'
path_to_radex+fileprefix+'_tau1.fits'
path_to_radex+fileprefix+'_tex2.fits'
path_to_radex+fileprefix+'_tau2.fits'

e.g. /Users/adam/work/n2hp/1-2_T=5to55_lvg_tau1.fits

N2H+ extras

In place of the actual contents of N2H+ fitter, here are the modules used to make the wrapper

class pyspeckit.spectrum.models.model.SpectralModel(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]

A wrapper class for a spectra model. Includes internal functions to generate multi-component models, annotations, integrals, and individual components. The declaration can be complex, since you should name individual variables, set limits on them, set the units the fit will be performed in, and set the annotations to be used. Check out some of the hyperfine codes (hcn, n2hp) for examples.

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) :

n2hp.n2hp_radex(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) [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

[edit this page on github]

[edit this page on bitbucket]