Rotational broadening profile

class PyAstronomy.modelSuite.RotBroadProfile

Implements rotational broadening with linear limb-darkening.

Fit Parameters:
  • xmax - Maximal extent of the profile

  • eps - Linear limb-darkening coefficient

  • A - Area under the profile (negative for absorption)

  • off - An offset applied to the profile

  • lin - Gradient of a linear term to adjust the ‘continuum’

  • mu - Center of the profile (same units as xmax)

  • gsig - The standard deviation of a Gaussian with which the

    rotational profile is convoluted, e.g., to model instrumental resolution.

The profile is given by:

\[G(x) = A \left(c_1\sqrt(1-(x/x_{max})^2) + c_2(1-(x/x_{max})^2)\right) + off + lin \cdot x\]

with the constants given by:

\[c_1 = \frac{2(1-\epsilon)}{\pi x_{max} (1-\epsilon/3)} \;\;\;\; c_2 = \frac{\epsilon}{2 x_{max} (1-\epsilon/3)}\]

Here, x can either denote a velocity or a wavelength shift. Thus, xmax can be given in km/s or Angstrom depending on the input; see, e.g., “Stellar Photospheres” by D.F. Gray for a derivation.

Note that the profile is normalized, i.e., the area amounts to one. This may, however, depend on the sampling. If the profile is undersampled, errors can become significant.

Methods

MCMCautoParameters(ranges[, picky, ...])

Convenience function to generate parameters for MCMC fit.

addConditionalRestriction(*args)

Define a conditional restriction.

assignValue(specval)

Assign new values to variables.

assignValues(specval)

Assign new values to variables.

autoFitMCMC(x, y, ranges[, picky, stepsize, ...])

Convenience function to using auto-generated sampling parameters in MCMC.

availableParameters()

Provides a list of existing parameters.

delRestriction(parName)

Delete restriction

description([parenthesis])

Returns a description of the model based on the names of the individual components.

errorConfInterval(par[, dstat, statTol, ...])

Calculate confidence interval for a parameter.

evaluate(v)

Calculates the rotational broadening profile according to current parameter values.

fit(x, y[, yerr, X0, minAlgo, mAA, ...])

Carries out a fit.

fitEMCEE([x, y, yerr, nwalker, priors, ...])

MCMC sampling using emcee package.

fitMCMC(x, y, X0, Lims, Steps[, yerr, ...])

Carry out MCMC fit/error estimation.

freeParamNames()

Get the names of the free parameters.

freeParameters()

Get names and values of free parameters.

freeze(specifiers)

Consider variables free to float.

frozenParameters()

Get names and values of frozen parameters.

getRelationsOf(specifier)

Return relations of a variable.

getRestrictions()

Get all restrictions.

hasVariable(specifier)

Determine whether the variable exists.

numberOfFreeParams()

Get number of free parameters.

parameterSummary([toScreen, prefix, sorting])

Writes a summary of the parameters in text form.

parameters()

Obtain parameter names and values.

relate(dependentVar, independentVars[, func])

Define a relation.

removeConditionalRestriction(*args)

Remove an existing conditional constraint.

renameVariable(oldName, newName)

Change name of variable.

restoreState(resource)

Restores parameter values from file or dictionary.

saveState(*args, **kwargs)

Save the state of the fitting object.

setObjectiveFunction([miniFunc])

Define the objective function.

setPenaltyFactor(penalFac)

Change the penalty factor.

setRestriction(restricts)

Define restrictions.

setRootName(root[, rename])

Define the root name of the model.

showConditionalRestrictions(**kwargs)

Show conditional restrictions.

steppar(pars, ranges[, extractFctVal, quiet])

Allows to step a parameter through a specified range.

thaw(specifiers)

Consider variables fixed.

untie(parName[, forceFree])

Remove all relations of parameter parName, i.e., the parameter is not dependend on other parameters.

updateModel()

Recalculate the model using current settings.

evaluate(v)

Calculates the rotational broadening profile according to current parameter values.

Parameters: x : array

Wavelength or velocity.

Example of usage

import numpy as np
import matplotlib.pylab as plt
from PyAstronomy import modelSuite as ms

# Get an instance of the model ...
x = ms.RotBroadProfile()
# ... and define some starting value
x["xmax"] = 60.0
x["A"] = 1.0
x["eps"] = 0.8
x["off"] = 0.0

# Define a radial velocity axis
vv = np.linspace(-90., 90., 200)

# Construct some "data" and ...
data = x.evaluate(vv)
# ... add noise
data += np.random.normal(0.0, 1e-3, data.size)

# Fit the model using A, xmax, and eps as free
# parameters ...
x.thaw(["A", "xmax", "eps"])
x.fit(vv, data)
# ... and show the resulting parameter values.
x.parameterSummary()

# Plot the data and the model
plt.plot(vv, data, 'bp')
plt.plot(vv, x.model, 'r--')
plt.show()