Thin shell transit-model

Transit light-curves usually have one minimum and a “U” shape. This can be different, e.g., if optically thin chromospheric emission-lines are considered. In this case, most emission may come from the stellar limb resulting in a “W”-shaped transit-profile.

The model has been presented by Schlawin et al. 1.

The model class

class PyAstronomy.modelSuite.XTran.limBrightTrans.LimBrightTrans

Planetary transit light-curves for spherical shell model.

This class implements a model calculating the light curve of a planet transiting an optically thin spherical shell of negligible thickness (e.g., a stellar chromosphere).

The model provided by Schlawin et al. 2010 assumes that the thickness of the shell is much smaller than the size of the planet. The shell is optically thin and thus provides natural limb-brightening. The obscured part of the stellar surface is calculated based on computing the volume of the intersection of a sphere with a cylinder and then taking a partial derivative with respect to the radius of the sphere to find its surface area.

The code closely follows the IDL procedure located at http://www.astro.washington.edu/agol/.

Fit parameters:
  • p - Rp/Rs (ratio of planetary and stellar radius)

  • a - Semi-major axis of planetary orbit [stellar radii].

  • per - Orbital period [d]

  • T0 - Central transit time

  • i - Inclination of orbit [rad]

By default all parameters remain frozen.

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(time)

Calculate a light curve according to the analytical model

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(time)

Calculate a light curve according to the analytical model by Schlawin et al. 2010.

Parameters
timearray

An array of time points at which the light curve shall be calculated

Returns
Modelarray

The analytical light curve is stored in the property lightcurve.

Notes

Note

time = 0 -> Planet is exactly in the line of sight (phase = 0).

Example code - Calculate model light curve

# Import some unrelated modules
import matplotlib.pylab as plt
from numpy import arange
# ... and now the LimBrightTrans module
from PyAstronomy.modelSuite import LimBrightTrans


# Create LimBrightTransit object
lbt = LimBrightTrans()
# Set parameters
lbt["p"] = 0.08
lbt["a"] = 6.70
lbt["i"] = 87.84
lbt["T0"] = 4.0
lbt["per"] = 10.

# Choose some time axis and calculate model
time = arange(3., 5., 0.001)
y = lbt.evaluate(time)

# Let's see what happened...
plt.ylabel("Relative Flux")
plt.xlabel("Time")
plt.plot(time, y, 'r-')
plt.show()
1

Schlawin et al. 2010, “Exoplanetary Transits of Limb-brightened Lines: Tentative Si IV Absorption by HD 209458b”, 2010ApJ…722L..75S