kinms_fitter

Submodules

Package Contents

Classes

kinms_fitter

Wrapper for easily kinematically modelling datacubes (typically from ALMA or VLA) using KinMS.

sb_profs

Creates flexible surface brightness profiles that can be combined together.

velocity_profs

warp_funcs

Functions

transformClouds(inclouds[, posAng, inc, cent, sbRad])

Transforms skySampler cloudlets for use in KinMS.

class kinms_fitter.kinms_fitter(filename, spatial_trim=None, spectral_trim=None, linefree_chans=[1, 5])

Wrapper for easily kinematically modelling datacubes (typically from ALMA or VLA) using KinMS. …

xc_guessfloat

Guess for Kinematic/morphological centre of the disc. RA, in degrees. [Default = Centre of supplied cube after clipping]

yc_guessfloat

Guess for Kinematic/morphological centre of the disc. Dec, in degrees. [Default = Centre of supplied cube after clipping]

vsys_guessfloat

Guess for Systemic velocity of the gas, in the radio velocity convention. km/s. [Default = Velocity of central channel of supplied cube after clipping]

totflux_guessfloat

Total flux guess. Units the same as those used in the cube (e.g. Jy km/s). [Default = total sum of the input cube] velDisp_guess : float Velocity dispersion guess. km/s. [Default = 8 km/s]

pa_guessfloat

Optional position angle guess, if no warp specified. Degrees. [Default = 0 deg]

inc_guessfloat

Optional Inclination angle guess, if no warp specified. Degrees. [Default = 45 deg]

inc_profilelist of warp_profile objects

Optional variation of the inclination with radius as described by warp_profile objects. [Default = None]

pa_profilelist of warp_profile objects

Optional variation of the inclination with radius as described by warp_profile objects. [Default = None]

sb_profilelist of sb_profs objects

Optional - Define the surface brightness profile with a list of sb_profs objects which will be evaluated in order. [Default = single exponential disc] expscale_guess : float Optional - defines the scalelength of the default single exponential disc used if sb_profile is not set. arcseconds. [Default = a fifth of the total extent of the cube] skySampClouds : ndarray Optional - Input cloud model made with the skySampler tool if you dont want to fit the gas distribution.

vel_profilelist of vel_profs objects

Define the circular velocity profile with a list of vel_profs objects. [Default = tilted rings]

radial_motionlist of radial_motion objects

Optional - Define the radial motion types with radial_motion objects (see base KinMS documentation). [Default: None]

nringsint

Optional - Number of tilted rings to use, spaced a beamwidth apart. [Default = fill the total extent of the cube]

vel_guessfloat

Optional - starting guess for the tilted ring asymptotic velocity. km/s [Default = 28percent of the velocity width coverved by the cube]

xcent_rangetwo element ndarray

Allowed range for fitting of the kinematic centre. RA in degrees [Default: minimum and maximum RA within the input cube]

ycent_rangetwo element ndarray

Allowed range for fitting of the kinematic centre. Dec in degrees [Default: minimum and maximum Dec within the input cube] vsys_range : two element ndarray Allowed range for fitting of the systemic velocity. km/s. [Default: minimum and maximum velocity contained within the input cube]

totflux_rangetwo element ndarray

Allowed range for fitting of the total flux. Units the same as input cube. [Default: Min=0 and Max=3x the total flux in the input cube]

velDisp_rangetwo element ndarray

Allowed range for fitting of the velocity dispersion. km/s. [Default: Min=0 and Max=50 km/s]

pa_rangetwo element ndarray

Optional - Allowed range for fitting of the position angle. Degrees. [Default: Min=0 and Max=360 degrees]

inc_rangetwo element ndarray

Optional - Allowed range for fitting of the inclination angle. Degrees. [Default: Min=0 and Max=89 degrees]

expscale_rangetwo element ndarray

Optional - Allowed range for fitting of the default exponential scale radius. arcseconds. [Default: Min=0 and Max=the size of your input cube]

vel_rangetwo element ndarray

Optional - Allowed range for fitting of the tilted ring velocities. km/s. [Default: Min=0 and Max=half the cube total velocity width]

pa_priorobject

Optional - supply an object defining the PA prior (see GAStimator documentation). [Default = boxcar between range limits]

xc_priorobject

Optional - supply an object defining the RA centre prior (see GAStimator documentation). [Default = boxcar between range limits]

yc_priorobject

Optional - supply an object defining the Dec centre prior (see GAStimator documentation). [Default = boxcar between range limits]

vsys_priorobject

Optional - supply an object defining the Vsys prior (see GAStimator documentation). [Default = boxcar between range limits]

inc_priorobject

Optional - supply an object defining the inc prior (see GAStimator documentation). [Default = boxcar between range limits]

totflux_priorobject

Optional - supply an object defining the totflux prior (see GAStimator documentation). [Default = boxcar between range limits]

velDisp_priorobject

Optional - supply an object defining the velocity dispersion prior (see GAStimator documentation). [Default = boxcar between range limits] initial_guesses : ndarray Override the initial guesses used for the fit. [Default = No override]

objname : string Name of the object you are fitting. [Default = read from FITS header OBJECT entry, or ‘object’] nprocesses : int Number of processors to use for MCMC fitting. [Default = fall back on GAStimator default choices]

nitersint

Number of iterations to use when fitting. [Default=3000]

output_initial_modelbool

Save the initially generated model as a FITS file. [Default=False]

pdfbool

Make PDF outputs. [Default = True]

silentbool

Fit silently. [Default = False]

show_cornerbool

Show corner plots if MCMC fitting. [Default = True]

nSampsint

Number of samples to use for KinMS model. [Default = 5e5] chi2_var_correct : bool Use chi-sqr variance correction. [Default = True] text_output : bool Save output text files. [Default = True]

save_all_acceptedbool

Save all the accepted guesses (allowing you to recreate corner plots etc). [Default = True]

interactivebool

Use interactive mode for plots. [Default = True]

show_plotsbool

Show the diagnotic plots. [Default = True]

output_cube_filerootstring

File root for output best fit datacube. [Default = False] lnlike_func : object Override the default likelihood function. [Default =. Dont override] tolerance : float Tolerance for simple fit. Smaller numbers are more stringent (longer runtime) [Default=0.1]

bunit : string Unit of datacube fluxes

cellsizefloat

Cube cellsize in arcsec timetaken : float Time taken for the fit in seconds.

chi2_correct_facfloat

chi2 correction factor applied.

mask_sum0

Number of detected pixels in the data mask.

labelsndarray of strings

Label for each fitted quantity

colorbar(mappable, ticks=None)

Add a colorbar to a given plot

param mappablematplotlib figure

Figure which to add the colorbar too

param ticksndarray

Ticks to show on the colourbar. [Default = None]

read_in_a_cube(path)

Reads in the datacube.

get_header_coord_arrays(hdr)

Get coordinate arrays from a FITS header

rms_estimate(cube, chanstart, chanend)

Estimate the RMS in the inner quarter of the datacube in a given channel range.

from_fits_history(hdr)

Get beam parameters if they happen to be in history keywords.

Stolen from radio_beam, with thanks!

read_primary_cube(cube)

Wrapper method for reading in datacube.

setup_params()

Setup the fit guesses, minima, maxima etc from the inputs.

clip_cube()

Clip the input datacube to size.

model_simple(param, fileName='')

Function to create a model from the given input parameters.

mcmc_fit(initial_guesses, labels, minimums, maximums, fixed, priors, precision)

Function to run the MCMC fit.

simple_chi2(theargs, info)

Likelihood function for the simple fit mode.

Returns the KinMS logo.

simple_fit(initial_guesses, labels, minimums, maximums, fixed)

Runs the simple fit.

plot(block=True, overcube=None, savepath=None, **kwargs)

Makes the plots.

write_text(bestvals, errup, errdown, units, fixed, runtime, mode, errors_warnings='None', fname='KinMS_fitter_output.txt')

Write output text files.

run(method='mcmc', justplot=False, savepath='./', **kwargs)

Run the fit.

Parameters
  • method (string - choice of 'mcmc', 'simple' or 'both'.) – Type of fit to run.

  • justplot (bool) – Just show the initial model plot without fitting. Helpful for setting initial guesses.

  • savepath (string) – Path to append to saved files.

  • **kwargs

    any additional keywords to pass to other methods.

class kinms_fitter.sb_profs

Creates flexible surface brightness profiles that can be combined together.

class expdisk(guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Creates an exponentially declining surface brightness profile.

guessesndarray of float

Initial guesses. If a single element, then the disc profile is normalised (peak=1). If two elements then elements should be [‘PeakFlux_exp’,’Rscale_exp’]. Rscale units of arcsec.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class gaussian(guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Creates an gaussian surface brightness profile.

guessesndarray of float

Initial guesses. If a two elements, then the gaussian profile is normalised (peak=1). If three elements then elements should be [‘PeakFlux_gauss’,’Mean_gauss’,’sigma_gauss’]. Units of mean/sigma are arcsec.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class cutoff(guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Creates an cutoff in a surface brightness profile between two radii.

guessesndarray of float

Initial guesses. Two elements [‘Start_cutoff’,’End_cutoff’]. Units of arcsec.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

eval(r, params)

Evaluates a list of surface brightness profiles.

modellistlist of objects

List of sb_prof objects, or objects that have the same methods/inputs.

r: ndarray of floats

Radius array, units of arcseconds.

params: ndarray of floats

Parameters to use in each model in the list.

outndarray of floats

Output combined surface brightness profile.

class kinms_fitter.velocity_profs

Circular velocity curves from physical and arbitary models that can be combined together.

class tilted_rings(bincentroids, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Arbitary velocity profile using the tilted ring formalism.

bincentroids: ndarray of float

Radii at which to constrain the profile. Linearly interpolated between these points.

guessesndarray of float

Initial guesses. Size of bincentroids, units of km/s.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class keplarian(distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Keplarian velocity profile for a point mass at r=0.

distance: ndarray of float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. One element log10 of the central mass (in Msun).

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class mge_vcirc(surf, sigma, qobs, distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Evaulate an MGE model of the potential, with or without a central point mass.

surfndarray of float

Luminosity of each gaussian component, units of Lsun.

sigmandarray of float

Width of each gaussian, units of arcsec.

qobsndarray of float

Axial ratio of each gaussian.

distance: ndarray of float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. One or two elements [M/L and optionally log10_CentralMass in Msun].

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class arctan(guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Arctangent velocity profile.

guessesndarray of float

Initial guesses. Vmax and Rturn in units of km/s and arcseconds.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class radial_barflow(guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Radial barflow from Spekkens+Sellwood.

guessesndarray of float

Initial guesses. Vtangential, Vradial, BarRadius and BarPA in units of km/s, km/s, arcsec, degrees.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class sersic(distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Velocity curve arising from a sersic mass distribution (spherically symmetric).

distance: float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. Total mass, effective radius and sersic index. Units of log10(Msun), arcsec and unitless.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

mass(r, theargs, norm=1)
__repr__()

Return repr(self).

__call__(x, theargs, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class bulge_disc(distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Velocity curve arising from a combination of a n=4 and n=1 sersic mass distributions (spherically symmetric).

distance: float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. Total mass, effective radius and sersic index. Units of log10(Msun), arcsec and unitless.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(x, theargs, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class nfw(distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None, hubbleconst=68.0)

Velocity curve arising from an NFW halo.

distance: float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. M200 and concentration. Units of log10(Msun) and unitless.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(xs, theargs, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class exponential_disc(distance, guesses, minimums, maximums, priors=None, precisions=None, fixed=None)

Velocity curve arising from a razor thin exponential disc. Based on eqn 8.74 in “Dynamics and Astrophysics of Galaxies” by Bovy.

distance: float

Distance to the object in Mpc.

guessesndarray of float

Initial guesses. Mdisk and Scaleradius. Units of log10(Msun) and arcsec.

minimumsndarray of float

Minimums for the given parameter.

maximumsndarray of float

Maximums for the given parameter.

priorsndarray of objects

Optional- Priors for the given parameter (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameter. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = False]

__repr__()

Return repr(self).

__call__(xs, theargs, **kwargs)

Returns the required profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

eval(r, params, inc=90)

Evaluates a list of velocity_profs objects, returning the total circular velocity profile.

modellistlist of objects

List of sb_prof objects, or objects that have the same methods/inputs.

r: ndarray of floats

Radius array, units of arcseconds.

params: ndarray of floats

Parameters to use in each model in the list.

inc: float

Inclination of system, in degrees.

kinms_fitter.transformClouds(inclouds, posAng=90.0, inc=0.0, cent=[0.0, 0.0], sbRad=None)

Transforms skySampler cloudlets for use in KinMS.

Calculate the galaxy co-ordinates of clouds from the sky plane. This MUST be used if any of the following conditions are true: inc != 0 posAng != 90 cent != [0,0]

This exists as a stand-alone routine since an MCMC fit to the galaxy will likely need to run this every step, and the sampleClouds routine is computationally expensive ============ Inputs:

clouds: np.ndarray The output of the sampleClouds array [x,y,I]

posAng: 0=<float<360 The position angle of the galaxy major axis, measured from the y-axis of the cube

inc: 0=<float<90 The inclination of the galaxy, with 90 being edge-on

cent: [float,float] The photometric centre of the galaxy relative to the centre of the cube

Outputs:

clouds: np.ndarray Positions and intensities of particles [x’,y’,I]

class kinms_fitter.warp_funcs

Creates position angle/inclination warps as a function of radius.

class linear(guesses, minimums, maximums, priors=None, precisions=None, fixed=None, labels_prefix='PA')

Creates a warp thats linear with radius (optionally: that also flattens after some radius).

guessesndarray of float

Initial guesses. Two or three elements [Gradient, Intercept, (optionally Cutoff_Radius)]. Units of deg/arcsec, deg, (arcsec).

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args)

Returns the required warp profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class flat(guesses, minimums, maximums, priors=None, precisions=None, fixed=None, labels=['PA'], units=['deg'])

No warp- inc/PA flat with radius.

guessesndarray of float

Initial guess. One element, units of degrees.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args)

Returns the required warp profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

class tilted_rings(bincentroids, guesses, minimums, maximums, priors=None, precisions=None, fixed=None, labels_prefix='PA')

Arbitary warp using the tilted ring formalism.

bincentroids: ndarray of float

Radii at which to constrain the profile. Linearly interpolated between these points.

guessesndarray of float

Initial guesses. Size of bincentroids, units of degrees.

minimumsndarray of float

Minimums for the given parameters.

maximumsndarray of float

Maximums for the given parameters.

priorsndarray of objects

Optional- Priors for the given parameters (see GAStimator priors).

precisionsndarray of float

Optional - Precision you want to reach for the given parameters. [Default = 10 percent of range]

fixed: ndarray of bool

Optional - Fix this parameter to the input value in guesses. [Default = All False]

__repr__()

Return repr(self).

__call__(x, args, **kwargs)

Returns the required warp profile.

xndarray of float

Input radial array in arcseconds

argsndarray of float

Input arguments to evalue the profile with

eval(r, params)

Evaluates a list of warp functions.

modellistlist of objects

List of warp_funcs objects, or objects that have the same methods/inputs.

r: ndarray of floats

Radius array, units of arcseconds.

params: ndarray of floats

Parameters to use in each model in the list.

outndarray of floats

Output combined warp profile.