pytmosph3r.observations.lightcurve

Module Contents

class Lightcurve(rotate=None, n_transmittances=None, n_phases=None, n_egress=None, phases=None, times=None, start_egress=None, end_egress=None, wns=None, wls=None, atmosphere_only=False, store_transmittances=None, transmittance_surfaces=None, rays=None, **kwargs)[source]

Bases: pytmosph3r.observations.transmission.Transmission, pytmosph3r.util.geometry.CircleIntersection, exo_k.util.spectral_object.Spectral_object

The lightcurve model calculates the transit depth \((R_P/R_S)^2\) for all phases of the transit via compute(). It is possible to take into account the tidally-locked rotation of the planet.

Note

rotate will be considerably slower since we have to recompute the transmittance for all phases.

partial_transit() calculates the intersection of the transmittance maps (or rays opacities) with the star.

The parameters of the lightcurve model are below.

Parameters:
  • rotate (bool) – Activates the planet rotation (tidally-locked) during transit (requires recomputing of the whole transmittance at each phase, so this option is deactivated by default). If False, we simply use the default transmittance map for all phases (decided by the longitude of the Observer).

  • n_phases (int) – Total number of phases to calculate (not used if phases is given). Defaults to 30.

  • n_egress (int) – Number of phases during partial transit (ingress+egress). ‘auto’ will set it to 2/3 of n_phases.

  • phases (array, float) – List of phases (in \(rad\) or astropy units). Overwrites n_phases.

  • times (array, float) – List of timesteps (in \(s\) or astropy units). Needs the orbital period to be set (see Orbit). Overwrites phases and n_phases.

  • n_transmittances (int) – Number of phases for which to effectively compute the transmittance (the rest will be interpolated). WARNING: a larger number will considerably be slower (same reasons as for rotate).

  • wns (array, float) – NOT USED FOR NOW. Use ~pytmosph3r.opacity.Opacity.wn_range instead. List of wavenumbers for which to calculate the lightcurve. Wavenumbers (in \(cm^{-1}\)) for which to compute lightcurves. Defaults to the full list of wavenumbers of the input data.

  • wls (array, float) – NOT USED FOR NOW. Wavelengths (in \(\mu m\)) for which to compute lightcurves.

  • atmosphere_only (bool) – Outputs the effect of the atmosphere only (without the planet core radius), useful mainly for debugging and checking the signal of the atmosphere.

property n_in_star

Number of phases during full transit (equal to n_phases - n_egress).

property times
property end_transit

Returns phase corresponding to the end of the transit (simulation is completely out of the star.

store_transmittance = True

Store transmittance in Transmission().

store_transmittances

Store the list of rays_opacities in (in particular when rotate is used).

flux

Main output of the module. (Rp/Rs)**2 for all phases.

build(model)[source]

Builds an atmospheric grid based on altitude coordinates. See AltitudeAtmosphere.

default_values()[source]

Set default values. Called by compute().

star_rays_opacity(phase, rays_opacity, Rs=None)[source]

Returns the rays opacities to the star flux (using star radius), and the distance of each ray to the edge of the star ( negative values are ‘inside’ the star).

Parameters:
  • phase (float) – Phase to calculate star opacity for.

  • opacity (array) – Original rays opacity at this phase.

  • Rs (optional) – Star radius

Returns:

opacity, distance to star edge (both (N_r x N_a))

Return type:

tuple

partial_transit(phase, rays_opacity)[source]

Computes the value of the spectrum when the planet is partially (or completely) in front of the star, using intersections().

Parameters:

phase (float) – phase

Returns:

rays_opacity (N_r x N_a)

Return type:

array

get_spectral_chunks(wns)[source]

Get spectral chunks corresponding to wns in opacity database.

compute_rays_opacities()[source]

Compute rays opacities (transmittance) at n_transmittances stages of the transit. These will server as a basis for interpolation in compute().

compute(model)[source]

Compute the lightcurve. Transmittances are calculated by compute_rays_opacities() on n_transmittances points and interpolated over n_phases. The surface area of the star covered by the transmittance map is calculated by intersections().

inputs()[source]
outputs()[source]