:py:mod:`pytmosph3r.observations.lightcurve` ============================================ .. py:module:: pytmosph3r.observations.lightcurve Module Contents --------------- .. py: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) Bases: :py:obj:`pytmosph3r.observations.transmission.Transmission`, :py:obj:`pytmosph3r.util.geometry.CircleIntersection`, :py:obj:`exo_k.util.spectral_object.Spectral_object` The lightcurve model calculates the transit depth :math:`(R_P/R_S)^2` for all phases of the transit via :func:`compute`. It is possible to take into account the tidally-locked rotation of the planet. .. note:: :attr:`rotate` will be considerably slower since we have to recompute the transmittance for all phases. :func:`partial_transit` calculates the intersection of the transmittance maps (or rays opacities) with the star. The parameters of the lightcurve model are below. :param rotate: 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). :type rotate: bool :param n_phases: Total number of phases to calculate (not used if :attr:`phases` is given). Defaults to 30. :type n_phases: int :param n_egress: Number of phases during partial transit (ingress+egress). 'auto' will set it to 2/3 of :attr:`n_phases`. :type n_egress: int :param phases: List of phases (in :math:`rad` or astropy units). Overwrites :attr:`n_phases`. :type phases: array, float :param times: List of timesteps (in :math:`s` or astropy units). Needs the orbital period to be set (see :class:`~pytmosh3r.orbit.Orbit`). Overwrites :attr:`phases` and :attr:`n_phases`. :type times: array, float :param n_transmittances: 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 :attr:`rotate`). :type n_transmittances: int :param wns: NOT USED FOR NOW. Use `~pytmosph3r.opacity.Opacity.wn_range` instead. List of wavenumbers for which to calculate the lightcurve. Wavenumbers (in :math:`cm^{-1}`) for which to compute lightcurves. Defaults to the full list of wavenumbers of the input data. :type wns: array, float :param wls: NOT USED FOR NOW. Wavelengths (in :math:`\mu m`) for which to compute lightcurves. :type wls: array, float :param atmosphere_only: Outputs the effect of the atmosphere only (without the planet core radius), useful mainly for debugging and checking the signal of the atmosphere. :type atmosphere_only: bool .. py:property:: n_in_star Number of phases during full transit (equal to :attr:`n_phases` - :attr:`n_egress`). .. py:property:: times .. py:property:: end_transit Returns phase corresponding to the end of the transit (simulation is completely out of the star. .. py:attribute:: store_transmittance :value: True Store transmittance in Transmission(). .. py:attribute:: store_transmittances Store the list of rays_opacities in (in particular when :attr:`rotate` is used). .. py:attribute:: flux Main output of the module. (Rp/Rs)**2 for all :attr:`phases`. .. py:method:: build(model) Builds an atmospheric grid based on altitude coordinates. See :class:`~pytmosph3r.atmosphere.AltitudeAtmosphere`. .. py:method:: default_values() Set default values. Called by :func:`~pytmosph3r.transmission.Transmission.compute()`. .. py:method:: star_rays_opacity(phase, rays_opacity, Rs=None) 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). :param phase: Phase to calculate star opacity for. :type phase: float :param opacity: Original rays opacity at this phase. :type opacity: array :param Rs: Star radius :type Rs: optional :returns: opacity, distance to star edge (both (N_r x N_a)) :rtype: tuple .. py:method:: partial_transit(phase, rays_opacity) Computes the value of the spectrum when the planet is partially (or completely) in front of the star, using :func:`~pytmosph3r.util.geometry.CircleIntersection.intersections`. :param phase: phase :type phase: float :returns: rays_opacity (N_r x N_a) :rtype: array .. py:method:: get_spectral_chunks(wns) Get spectral chunks corresponding to `wns` in opacity database. .. py:method:: compute_rays_opacities() Compute rays opacities (transmittance) at :attr:`n_transmittances` stages of the transit. These will server as a basis for interpolation in :func:`compute`. .. py:method:: compute(model) Compute the lightcurve. Transmittances are calculated by :func:`compute_rays_opacities` on :attr:`n_transmittances` points and interpolated over :attr:`n_phases`. The surface area of the star covered by the transmittance map is calculated by :func:`~pytmosph3r.util.geometry.CircleIntersection.intersections`. .. py:method:: inputs() .. py:method:: outputs()