MulensModel package¶
Subpackages¶
Submodules¶
- MulensModel.b0b1utils module
- MulensModel.binarylens module
- MulensModel.binarylensimports module
- MulensModel.binarylenswithshear module
- MulensModel.causticsbinary module
- MulensModel.causticsbinarywithshear module
- MulensModel.causticspointwithshear module
- MulensModel.coordinates module
- MulensModel.elliputils module
- MulensModel.event module
- MulensModel.fitdata module
- MulensModel.horizons module
- MulensModel.limbdarkeningcoeffs module
- MulensModel.magnificationcurve module
- MulensModel.model module
- MulensModel.modelparameters module
- MulensModel.mulensdata module
- MulensModel.pointlens module
- MulensModel.pointlenswithshear module
- MulensModel.satelliteskycoord module
- MulensModel.trajectory module
- MulensModel.uniformcausticsampling module
- MulensModel.utils module
- MulensModel.version module
Module contents¶
- class MulensModel.BinaryLensPointSourceWM95Magnification(**kwargs)¶
Bases:
MulensModel.binarylens._BinaryLensPointSourceMagnification
Equations for calculating point-source–binary-lens magnification following the Witt & Mao 1995, ApJL, 447, L105 prescription assuming a POINT source.
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- trajectory:
- class MulensModel.BinaryLensPointSourceVBBLMagnification(**kwargs)¶
Bases:
MulensModel.binarylens._BinaryLensPointSourceMagnification
Equations for calculating point-source–binary-lens magnification using VBBL for point sources.
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- trajectory:
- class MulensModel.BinaryLensQuadrupoleMagnification(gamma=None, u_limb_darkening=None, **kwargs)¶
Bases:
MulensModel.binarylens.BinaryLensPointSourceMagnification
,MulensModel.binarylens._LimbDarkeningForMagnification
,MulensModel.binarylens._FiniteSource
Magnification in quadrupole approximation of the binary-lens/finite-source event - based on Gould 2008 ApJ 681, 1593.
The origin of the coordinate system is at the center of mass and both masses are on X axis with higher mass at negative X; this means that the higher mass is at (X, Y)=(-s*q/(1+q), 0) and the lower mass is at (s/(1+q), 0).
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- gamma: float
Linear limb-darkening coefficient in gamma convention.
- trajectory:
- class MulensModel.BinaryLensHexadecapoleMagnification(all_approximations=False, **kwargs)¶
Bases:
MulensModel.binarylens.BinaryLensQuadrupoleMagnification
Magnification in hexadecapole approximation of the binary-lens/finite-source event - based on Gould 2008 ApJ 681, 1593.
For coordinate system convention see
BinaryLensQuadrupoleMagnification
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- gamma: float
Linear limb-darkening coefficient in gamma convention.
- all_approximations: boolean, optional
If True, return hexadecapole, quadrupole, and point source approximations. Default is *False.
- trajectory:
- class MulensModel.BinaryLensVBBLMagnification(gamma=None, u_limb_darkening=None, accuracy=0.001, **kwargs)¶
Bases:
MulensModel.binarylens._BinaryLensPointSourceMagnification
,MulensModel.binarylens._LimbDarkeningForMagnification
,MulensModel.binarylens._FiniteSource
Binary lens finite source magnification calculated using VBBL library that implements advanced contour integration algorithm presented by Bozza 2010 MNRAS, 408, 2188. See also VBBL website by Valerio Bozza.
For coordinate system convention see
BinaryLensQuadrupoleMagnification
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- gamma: float, optional
Linear limb-darkening coefficient in gamma convention.
- u_limb_darkening: float, optional
Linear limb-darkening coefficient in u convention. Note that either gamma or u_limb_darkening can be set. If neither of them is provided then limb darkening is ignored.
- accuracy: float, optional
Requested accuracy of the result.
- trajectory:
- class MulensModel.BinaryLensAdaptiveContouringMagnification(gamma=None, u_limb_darkening=None, accuracy=0.1, ld_accuracy=0.001, **kwargs)¶
Bases:
MulensModel.binarylens._BinaryLensPointSourceMagnification
,MulensModel.binarylens._LimbDarkeningForMagnification
,MulensModel.binarylens._FiniteSource
Binary lens finite source magnification calculated using Adaptive Contouring method by Dominik 2007 MNRAS, 377, 1679
See also AdaptiveContouring website by Martin Dominik
For coordinate system convention see
point_source_magnification()
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- gamma: float, optional
Linear limb-darkening coefficient in gamma convention.
- u_limb_darkening: float
Linear limb-darkening coefficient in u convention. Note that either gamma or u_limb_darkening can be set. If neither of them is provided then limb darkening is ignored.
- accuracy: float, optional
Requested accuracy of the result defined as the sum of the area of the squares that determine the contour line and the estimated total enclosed area (see sec. 4 of the paper). As M. Dominik states: “this vastly overestimates the fractional error, and a suitable value should be chosen by testing how its variation affects the final results - I recommend starting at acc = 0.1.” It significantly affects execution time.
- ld_accuracy: float, optional
Requested limb-darkening accuracy. As M. Dominik states: ” Fractional uncertainty for the adaptive Simpson integration of the limb-darkening profile-related function during application of Green’s theorem.” It does not add execution time so can be set to very small value.
- trajectory:
- class MulensModel.BinaryLensPointSourceWithShearWM95Magnification(convergence_K=None, shear_G=None, **kwargs)¶
Bases:
MulensModel.binarylens.BinaryLensPointSourceWM95Magnification
The binary lens with shear and convergence: solutions, images, parities, magnifications, etc.
Uses Witt & Mao 1995 for the binary lens magnification calculation. See py:class:binarylens.BinaryLensPointSourceWM95Magnification
The binary lens with shear and convergence equation is the 9th order complex polynomial.
- Attributes :
- mass_1: float
mass of the primary (left-hand object) as a fraction of the total mass.
- mass_2: float
mass of the secondary (right-hand object) as a fraction of the total mass.
- separation: float
separation between the two bodies as a fraction of the Einstein ring.
- convergence_K: float
External mass sheet convergence.
- shear_G: complex
External mass sheat shear.
Note: mass_1 and mass_2 may be defined as a fraction of some other mass than the total mass. This is possible but not recommended - make sure you know what you’re doing before you start using this possibility.
If you’re using this class, then please cite Peirson et al. (2022; ApJ 927, 24).
- class MulensModel.BinaryLensPointSourceWithShearVBBLMagnification(**kwargs)¶
Bases:
MulensModel.binarylenswithshear.BinaryLensPointSourceWithShearWM95Magnification
The binary lens with shear and convergence: solutions, images, parities, magnifications, etc.
Uses VBBL.
The binary lens with shear and convergence equation is the 9th order complex polynomial.
- Attributes :
- mass_1: float
mass of the primary (left-hand object) as a fraction of the total mass.
- mass_2: float
mass of the secondary (right-hand object) as a fraction of the total mass.
- separation: float
separation between the two bodies as a fraction of the Einstein ring.
- convergence_K: float
External mass sheet convergence.
- shear_G: complex
External mass sheat shear.
Note: mass_1 and mass_2 may be defined as a fraction of some other mass than the total mass. This is possible but not recommended - make sure you know what you’re doing before you start using this possibility.
- class MulensModel.CausticsBinary(q, s)¶
Bases:
object
Class for the caustic structure corresponding to a given (q, s), i.e. mass ratio and separation. Implemented for 2-body lenses only.
- Attributes :
- q: float
mass ratio between the 2 bodies; always <= 1
- s: float
separation between the 2 bodies (as a fraction of the Einstein ring)
- plot(n_points=5000, **kwargs)¶
Plots the caustics using
matplotlib.pyplot.scatter()
.- Parameters :
- n_points: int, optional
The number of points to calculate along the caustic. Defaults to 5000.
**kwargs
:keywords accepted by
matplotlib.pyplot.scatter()
Note that default scaling of axis may not be equal on both axis. To mitigate this, use:
plt.gca().set_aspect('equal')
orplt.axis('equal')
(the other possible options are'scaled'
or'square'
).
- get_caustics(n_points=5000)¶
Returns x and y vectors corresponding to the outlines of the caustics. Origin is center of mass and larger mass is on the left (for q < 1).
- Parameters:
- n_pointsint, optional
The number of points to calculate along the caustic.
- Returns:
- x, ylist
Two lists of length n_points giving the x, y coordinates of the caustic points.
- property critical_curve¶
Critical curve stored as
CriticalCurve
object, read-only
- class MulensModel.CausticsPointWithShear(convergence_K, shear_G)¶
Bases:
MulensModel.causticsbinary.CausticsBinary
Class for the caustic structure produced by the Chang-Refsdal lens, i.e. single mass with shear and convergence.
- Attributes :
- convergence_K: float
convergence of the lens system
- shear_G: complex
shear of the the lens system
- class MulensModel.CausticsBinaryWithShear(q, s, convergence_K=0.0, shear_G=0j)¶
Bases:
MulensModel.causticsbinary.CausticsBinary
Caustics structure for a binary lens with shear and convergence.
- Attributes :
- q: float
mass ratio between the 2 bodies; always <= 1
- s: float
separation between the 2 bodies (as a fraction of the Einstein ring)
- convergence_K: float
convergence of the lens system
- shear_G: complex
shear of the the lens system
- class MulensModel.Coordinates(*args, **kwargs)¶
Bases:
astropy.coordinates.sky_coordinate.SkyCoord
A class for the coordinates (RA, Dec) of an event. Inherits from astropy.SkyCoord.
May be set as a str, pair of str, or SkyCoord object, e.g.
from astropy.coordinates import SkyCoord from astropy import units as u Coordinates('18:00:00 -30:00:00') Coordinates('18h00m00s', '-30d00m00s') Coordinates(SkyCoord('18:00:00 -30:00:00', unit=(u.hourangle, u.deg))) Coordinates(SkyCoord(270.000, -30.000, unit=u.deg))
If the unit keyword is not specified, defaults to unit=(u.hourangle, u.deg) where u is defined by “import astropy.units as u”.
You can print an instance of this class.
- property galactic_l¶
astropy.coordinates.angles.Longitude
Galactic longitude. Note that for connivance, the values l > 180 degrees are represented as 360-l.
- property galactic_b¶
astropy.coordinates.angles.Latitude
Galactic latitude calculated from (RA, Dec)
- property ecliptic_lon¶
astropy.coordinates.angles.Longitude
ecliptic longitude calculated from (RA, Dec)
- property ecliptic_lat¶
astropy.coordinates.angles.Latitude
ecliptic latitude calculated from (RA, Dec)
- property north_projected¶
np.array
North projected on the plane of the sky.
- property east_projected¶
np.array
East projected on the plane of the sky.
- v_Earth_projected(full_BJD)¶
Earth velocity at full_BJD projected on the plane of sky towards given coordinates.
- Parameters :
- full_BJD: float
Epoch for which projected velocity is requested. In most cases it is
t_0_par
- Returns :
- v_Earth_perp_N: float
North component of Earth’s projected velocity in km/s.
- v_Earth_perp_E: float
East component of Earth’s projected velocity in km/s.
- class MulensModel.Event(datasets=None, model=None, coords=None, fix_blend_flux=None, fix_source_flux=None, fix_source_flux_ratio=None, data_ref=0)¶
Bases:
object
Combines a microlensing model with data. Allows calculating chi^2 and making a number of plots.
- Arguments :
- datasets:
MulensData
object or a list of such objects Datasets that will be linked to the event. These datasets will be used for chi^2 calculation, plotting etc.
- model:
Model
Microlensing model that will be linked to the event. In order to get chi^2 for different sets of model parameters you should keep a single
Model
instance and change parameters for this model (i.e., do not provide separateModel
instances).- coords: str,
Coordinates
, or astropy.SkyCoord Coordinates of the event. If str, then needs format accepted by astropy.SkyCoord e.g.,
'18:00:00 -30:00:00'
.- fix_blend_flux, fix_source_flux: dict
Used to fix the source flux(es) or blend flux for a particular dataset. The dataset is the key, and the value to be fixed is the value. For example, to fix the blending of some dataset my_data to zero set fix_blend_flux={my_data: 0.}. See also
FitData
.- fix_source_flux_ratio: dict
Used to fix the flux ratio for a given band or dataset. The keys should be either
MulensData
objects or str. If aMulensData
object is specified, it will take precedence over a band.- data_ref: int or
MulensData
Reference dataset. If int then gives index of reference dataset in
datasets
. Default is the first dataset.
- datasets:
The datasets can be in magnitude or flux spaces. When we calculate chi^2 we do it in magnitude or flux space depending on value of
chi2_fmt
attribute. If dataset is in magnitude space and model results in negative flux, then we calculate chi^2 in flux space but only for the epochs with negative model flux.You can print an instance of this class. Information on model and datasets will be provided.
- plot(t_range=None, residuals=True, show_errorbars=None, show_bad=None, legend=True, trajectory=None, title=None, subtract_2450000=True, subtract_2460000=False, data_ref=None)¶
Basic plotting. Default is to plot the light curve with the residuals. If the model has 2 or more lenses, also plot the source trajectory with caustics. For more detailed control over the plotting, see
plot_model()
,plot_data()
, andplot_trajectory()
.- Keywords:
- t_range: list, tuple
Time range over which to show the light curve plot.
- residuals: bool
Whether or not to plot the residuals with the light curve.
- show_errorbars: bool
Whether or not to show the errorbars on the data. Default is None (see
plot_data()
).- show_bad: bool
Whether or not to show data points marked as “bad” (see
bad
). Default is None (seeplot_data()
).- legend: bool
Whether or not to show a legend for the datasets.
- trajectory: bool
Whether or not to plot the source trajectory. Defaults to False for single lens events and True for binary lenses.
- title: str
Title for the plot. Same title is used for trajectory plot, if applicable.
- subtract_2450000, subtract_2460000: bool
see
plot_data()
.- data_ref: int or MulensData
see
plot_data()
.
- plot_model(data_ref=None, **kwargs)¶
Plot the model light curve in magnitudes. See
MulensModel.model.Model.plot_lc()
for details.- Keywords :
- data_ref: int or MulensData
If data_ref is not specified, uses
data_ref
.**kwargs
:Keywords passed to
MulensModel.model.Model.plot_lc()
. You can use them to set time range plotted, fluxes, limb-darkening etc.
- plot_data(phot_fmt='mag', data_ref=None, show_errorbars=None, show_bad=None, subtract_2450000=False, subtract_2460000=False, **kwargs)¶
Plot the data scaled to the model.
- Keywords (all optional):
- phot_fmt: string (‘mag’, ‘flux’)
Whether to plot the data in magnitudes or in flux. Default is ‘mag’.
- data_ref: int or MulensData
If data_ref is not specified, uses
data_ref
.- show_errorbars: boolean or None
Do you want errorbars to be shown for all datasets? Default is None, which means the option is taken from each dataset plotting properties (for which default is True). If True, then data are plotted using matplotlib.errorbar(). If False, then data are plotted using matplotlib.scatter().
- show_bad: boolean or None
Do you want data marked as bad to be shown? Default is None, which means the option is taken from each dataset plotting properties (for which default is False). If bad data are shown, then they are plotted with ‘x’ marker.
- subtract_2450000, subtract_2460000: boolean
If True, subtracts 2450000 or 2460000 from the time axis to get more human-scale numbers. If using, make sure to also set the same settings for all other plotting calls (e.g.
plot_lc()
).**kwargs
:Passed to matplotlib plotting functions. Contrary to previous behavior,
**kwargs
are no longer remembered.
- plot_residuals(show_errorbars=None, data_ref=None, subtract_2450000=False, subtract_2460000=False, show_bad=None, **kwargs)¶
Plot the residuals (in magnitudes) to the model.
- Keywords:
For explanation of keywords, see doctrings in
plot_data()
. Note different order of keywords.
- plot_trajectory(**kwargs)¶
Plot the trajectory of the source. See
MulensModel.model.Model.plot_trajectory()
for details.
- plot_source_for_datasets(**kwargs)¶
Plot source positions for all linked datasets. See
MulensModel.model.Model.plot_source()
for details.Note: plots all points in datasets (including ones flagged as bad) using the same marker.
- get_flux_for_dataset(dataset)¶
Get the source and blend flux for a given dataset.
- Parameters :
dataset:
MulensData
or int If int should be the index (starting at 0) of the appropriate dataset in thedatasets
list.- Returns :
- source_flux: np.ndarray
flux of sources. see
source_fluxes
- blend_flux: float
blending flux. see
blend_flux
NOTE: This function does not recalculate fits or fluxes. If the data haven’t yet been fit to the model (i.e. self.fits = None), it will run
fit_fluxes()
. Otherwise, it just accesses the existing values. So if you change something inmodel
or some fit parameter (e.g.,fix_blend_flux
), be sure to runfit_fluxes()
first.
- get_ref_fluxes()¶
Get source and blending fluxes for the reference dataset. See
get_flux_for_dataset()
. If the reference dataset is not set, uses the first dataset as default. Seedata_ref
.
- get_chi2()¶
Calculates chi^2 of current model by fitting for source and blending fluxes.
- Returns :
- chi2: float
Chi^2 value
- get_chi2_for_dataset(index_dataset)¶
Calculates chi^2 for a single dataset
- Parameters :
- index_dataset: int
index that specifies for which dataset the chi^2 is requested
- Returns :
- chi2: float
chi2 for dataset[index_dataset].
- get_chi2_per_point(bad=False)¶
Calculates chi^2 for each data point of the current model by fitting for source and blending fluxes.
- Parameters :
- bad: bool
Should chi2 be also caclulated for points marked as bad in MulensData? If False (default), then bad epochs have chi2 of np.nan.
- Returns :
- chi2: list of np.ndarray
Chi^2 contribution from each data point, e.g.
chi2[data_num][k]
returns the chi2 contribution from the k-th point of dataset data_num.
- Example :
Assuming
event
is instance of Event class to get chi2 for 10-th point point of 0-th dataset.chi2 = event.get_chi2_per_point() print(chi2[0][10])
- get_chi2_gradient(parameters)¶
Fit for fluxes and calculate chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\).
- Parameters :
- parameters: str or list, required
Parameters with respect to which gradient is calculated. Currently accepted parameters are:
t_0
,u_0
,t_eff
,t_E
,pi_E_N
, andpi_E_E
. The parameters for which you request gradient must be defined in py:attr:~model.
- Returns :
- gradient: float or np.ndarray
chi^2 gradient
- calculate_chi2_gradient(parameters)¶
Calculate chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\).
- Parameters :
- parameters: str or list, required
Parameters with respect to which gradient is calculated. Currently accepted parameters are:
t_0
,u_0
,t_eff
,t_E
,pi_E_N
, andpi_E_E
. The parameters for which you request gradient must be defined in py:attr:~model.
- Returns :
- gradient: float or np.ndarray
chi^2 gradient
NOTE: Because this is not a ‘get’ function, it ASSUMES you have ALREADY fit for the fluxes, e.g. by calling get_chi2().
- fit_fluxes(bad=False)¶
Fit for the optimal fluxes for each dataset (and its chi2)
- property coords¶
see
Coordinates
- property datasets¶
a list of
MulensData
instances.
- property data_ref¶
Reference data set for scaling the model fluxes to (for plotting). May be set as a
MulensData
object or an index (int). Default is the first data set.- Returns :
index (int) of the relevant dataset.
- property chi2¶
float
Chi^2 value. Note this is a static property. It is only updated when
fit_fluxes()
orget_chi2()
is run. So, if you change one of the settings be sure to run one of those functions to update the chi2.
- property chi2_gradient¶
Return previously calculated chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\). See
get_chi2_gradient()
andcalculate_chi2_gradient()
.- Returns :
- gradient: float or np.ndarray
chi^2 gradient. Will return None if the chi2 gradient was not previously calculated using one of the functions mentioned above.
- property fits¶
list of
FitData
objectsThere is one
FitData
object for each dataset containing the information for fitting the model to that dataset, e.g. fitted fluxes, chi2 (for that dataset).
- property fluxes¶
list
An array giving the fitted source and blend flux(es) for each dataset.
- property source_fluxes¶
list
An array giving the fitted source flux(es) for each dataset.
- property blend_fluxes¶
list
An array giving the fitted blend flux for each dataset.
- property sum_function¶
str
Function used for adding chi^2 contributions. Can be either ‘numpy.sum’ (default value) or ‘math.fsum’. The latter is slightly slower and more accurate, which may be important for large datasets.
- class MulensModel.FitData(model=None, dataset=None, fix_blend_flux=False, fix_source_flux=False, fix_source_flux_ratio=False)¶
Bases:
object
Performs a least squares linear fit for given dataset and model to determine the source flux(es) and (optionally) blend flux. After creating the object, you must run
update()
to perform the linear fit for the fluxes and calculate the chi2. To perform the linear fit without calculating chi2, you can runfit_fluxes()
. If you change anything in the object, e.g. the model parameters, you must re-runupdate()
orfit_fluxes()
.- Arguments :
- model:
Model
object The model to fit to the data.
- dataset:
MulensData
object A single photometric dataset to be fitted.
- fix_blend_flux: False or float, optional
Default is False, i.e. allow the blend flux to be a free parameter. If set to a float, it will fix the blend value to that value.
- fix_source_flux: False, float, or list, optional
Default is False, i.e. allow the source flux to be a free parameter. If set to a float, it will fix the source value to that value. For binary source models, a list should be used to set the fluxes of the individual sources or fix one and not the other, e.g. [2.3, False] would fix source_flux_0 to 2.3 but allow a free fit to source_flux_1.
- fix_source_flux_ratio: False or float, optional
For binary source models, source_flux_ratio is the flux ratio between two components, i.e., source_flux_ratio = source_flux_1 / source_flux_0 Default is False, i.e. allow the source flux to be a free parameter. If set to a float, it will fix the source value to that value.
- model:
- update(bad=False)¶
Calculate the best-fit source and blend fluxes as well as the chi2.
- Parameters :
- bad: bool
Default is False. If True recalculates the data magnification for each point to ensure that there are values even for bad datapoints.
No returns.
- fit_fluxes()¶
Execute the linear least squares fit to determine the fitted fluxes. Sets the values of
source_fluxes
,blend_flux
, and (if applicable)source_flux
.Does not calculate chi2. To fit for the fluxes and calculate chi2, run
update()
.No parameters.
No returns.
- get_data_magnification(bad=False)¶
Calculates the model magnification for each data point.
- Parameters :
- bad: boolean
If True, calculates the magnification for all points. If False, only calculates the magnification for good data points. Values for bad data points are set to 0. Default is False.
- Returns :
- data_magnification: np.ndarray
The model magnification evaluated for each datapoint. If there is more than one source, the magnification of each source is reported separately.
- get_model_fluxes(bad=False)¶
Calculate model in flux space.
- Parameters :
- bad: bool
Default is False. If True recalculates the data magnification for each point to ensure that the values for bad datapoints are calculated (otherwise, they are set to the magnitude of the blend).
- Returns :
- model_flux: np.ndarray
The model flux evaluated for each datapoint.
- get_model_magnitudes(**kwargs)¶
Calculate model in magnitude space
- Parameters :
**kwargs
:
- Returns :
- model_mag: np.ndarray
The model magnitude evaluated for each datapoint.
- scale_fluxes(source_flux, blend_flux)¶
- Rescale the data fluxes to an arbitrary flux scale:
flux = source_flux_0 * (data.flux - blend_flux) / source_flux flux += blend_flux_0 err_flux = source_flux_0 * data.err_flux / source_flux
- Parameters :
- source_flux: float, list, np.array
Flux of the source in the desired system. If n_sources > 1 and source_flux has more than one element, the elements are summed to produce the overall scaling flux.
- blend_flux: float
Flux of the blend in the desired system
- Returns :
- flux: np.ndarray
Fluxes from the data rescaled to the desired system.
- err_flux: np.ndarray
Uncertainties of fluxes from the data rescaled to the desired system.
- get_residuals(phot_fmt=None, source_flux=None, blend_flux=None, bad=False)¶
Calculate the residuals for each datapoint relative to the model.
- Parameters :
- phot_fmt: str, optional
specify whether the residuals should be returned in magnitudes (‘mag’) or in flux (‘flux’). Default is ‘mag’. If ‘scaled’, will return the residuals in magnitudes scaled to source_flux, blend_flux.
- source_flux, blend_flux: float
reference source and blend fluxes for scaling the residuals
- bad: bool
Default is False. If True recalculates the data magnification for each point to ensure that there are values even for bad datapoints.
- Returns :
- residuals: np.ndarray
the residuals for the corresponding dataset.
- errorbars: np.ndarray
the scaled errorbars for each point. For plotting errorbars for the residuals.
- get_chi2_gradient(parameters)¶
Fits fluxes and calculates chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\).
- Parameters :
- parameters: str or list, required
Parameters with respect to which gradient is calculated. Currently accepted parameters are:
t_0
,u_0
,t_eff
,t_E
,pi_E_N
, andpi_E_E
. The parameters for which you request gradient must be defined in py:attr:~model.
- Returns :
- gradient: float or np.ndarray
chi^2 gradient
- calculate_chi2_gradient(parameters)¶
Calculates chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\) WITHOUT refitting for the fluxes. Saves computations if, e.g., you want to retrieve both py:attr:~chi2 and py:attr:~chi2_gradient.
- Parameters :
- parameters: str or list, required
Parameters with respect to which gradient is calculated. Currently accepted parameters are:
t_0
,u_0
,t_eff
,t_E
,pi_E_N
, andpi_E_E
. The parameters for which you request gradient must be defined in py:attr:~model.
- Returns :
- gradient: float or np.ndarray
chi^2 gradient
- get_d_A_d_params_for_point_lens_model(parameters)¶
Calculate d A / d parameters for a point lens model.
- Parameters :
- parameters: list
List of the parameters to take derivatives with respect to.
- Returns :
- dA_dparam: dict
Keys are parameter names from parameters argument above. Values are the partial derivatives for that parameter evaluated at each data point.
- get_d_A_d_rho()¶
Calculate d A / d rho for a point lens model.
No Inputs
- Returns :
- dA_drho: np.array
Values are the partial derivatives for rho evaluated at each data point.
- get_dataset_trajectory()¶
Retrieve a
Trajectory
object. If thedataset
has an ephemerides_file, apply it to the Trajectory, even if it is not part of themodel
.No parameters.
- Returns :
- trajectory:
Trajectory
Trajectory for given dataset.
- trajectory:
- get_d_A_d_u_for_PSPL_model()¶
- get_d_A_d_u_for_FSPL_model()¶
- get_d_A_d_u_for_point_lens_model()¶
- property chi2_gradient¶
float or np.ndarray
Previously calculated chi^2 gradient (also called Jacobian), i.e., \(d chi^2/d parameter\). See
get_chi2_gradient()
andcalculate_chi2_gradient()
.Gives None if the chi2 gradient was not previously calculated using one of the functions mentioned above.
- property chi2¶
float The total chi2 for the fitted dataset. Good points only. See
good
.If None, you need to run
update()
to execute the linear fit and calculate the chi2.
- property chi2_per_point¶
np.ndarray
The chi^2 contribution from each data point, e.g.,
chi2_per_point[k]
returns the chi2 contribution from the k-th point ofdataset
. Includes bad datapoints.If None, you need to run
update()
to execute the linear fit and calculate the chi2.
- property source_flux¶
float
The fitted source flux. Only defined for models with a single source. See also
source_fluxes
If None, you need to run
fit_fluxes()
orupdate()
to execute the linear fit.
- property source_fluxes¶
np.array
The fitted source flux(es).
If None, you need to run
fit_fluxes()
orupdate()
to execute the linear fit.
- property blend_flux¶
float
The fitted blend flux or the value set by fix_blend_flux (see Keywords).
If None, you need to run
fit_fluxes()
orupdate()
to execute the linear fit.
- property source_flux_ratio¶
float
source_flux_ratio = source_flux_1 / source_flux_0
i.e., the ratio of the fitted source fluxes or the value set by fix_source_flux_ratio (see Keywords).
If None, you need to run
fit_fluxes()
orupdate()
to execute the linear fit.
- property dataset¶
-
A single photometric dataset to be fitted.
- property data_magnification¶
Returns previously calculated magnifications. To calculate the magnifications (e.g., if something changed in the model), use
get_data_magnification()
.- Returns :
- data_magnification: np.ndarray
The model magnification evaluated for each datapoint. If there is more than one source, the magnification of each source is reported separately.
- property magnification_curve¶
Returns previously calculated magnification curve.
- property magnification_curves¶
Returns previously calculated magnification curves.
- Returns :
tuple of :py:class:`~MulensModel.magnification.MagnificationCurve objects, i.e., the model magnification curve evaluated for each datapoint.
- property gamma¶
float
Limb-darkening coefficient for this fit. Set by
bandpass
andget_limb_coeff_gamma()
.
- class MulensModel.Horizons(file_name)¶
Bases:
object
An Object to read and hold the standard JPL Horizons output, i.e. satellite ephemerides.
- Arguments :
- file_name: str
output from JPL Horizons file name
For info on preparing JPL Horizons file for import, see instructions.
- property time¶
np.ndarray
return times in TDB (reference frame depends on Horizons input)
- property xyz¶
Astropy.CartesianRepresentation
return X,Y,Z positions based on RA, DEC and distance
- class MulensModel.LimbDarkeningCoeffs¶
Bases:
object
Linear limb-darkening parameters. Both gamma and u conventions can be used. The u convention is more frequently used in studies other than microlensing. It has fixed flux at the center. An et al. 2002 (ApJ 572, 521) introduced the gamma convention:
gamma = (2 * u) / (3 - u)
u = 3 * gamma / (2 + gamma)
Note that the gamma convention has fixed total flux.
You can print an instance of this class.
- set_limb_coeff_gamma(bandpass, gamma)¶
Remembers limb darkening gamma coefficient for given band.
- Parameters :
- bandpass: str
Name of the filter.
- gamma: float
The value of gamma coefficient.
- get_limb_coeff_gamma(bandpass)¶
Gives limb darkening gamma coefficient for given band.
- Parameters :
- bandpass: str
Name of the filter.
- Returns :
- gamma: float
The value of gamma coefficient.
- set_limb_coeff_u(bandpass, u)¶
Remembers limb darkening u coefficient for given band
- Parameters :
- bandpass: str
Name of the filter.
- u: float
The value of u coefficient.
- get_limb_coeff_u(bandpass)¶
Gives limb darkening u coefficient for given band.
- Parameters :
- bandpass: str
Name of the filter.
- Returns :
- u: float
The value of u coefficient.
- get_weighted_limb_coeff_gamma(weights)¶
Get weighted limb darkening coefficient in gamma space.
- Parameters :
- weights: dict
A dictionary that for every band (keys; str type) gives its relative weight (value; float type), e.g.,
weights = {'I': 1.5, 'V': 1.}
will return gamma coefficient in the case when I band contributes 1.5 more than V band. Note that for each band used you have to first set to coefficient.
- Returns :
- gamma: float
The value of weighted gamma.
- class MulensModel.MagnificationCurve(times, parameters, parallax=None, coords=None, satellite_skycoord=None, gamma=0.0)¶
Bases:
object
The magnification curve calculated from the model light curve.
The key function is
set_magnification_methods()
, which specifies the method used to calculate the finite source magnification and when to use it..- Arguments :
- times: iterable of floats
the times at which to generate the magnification curve
- parameters:
ModelParameters
specifies the microlensing parameters
- parallax: dict, optional
dictionary specifying what parallax effects should be used, e.g.,
{'earth_orbital': True, 'satellite': False, 'topocentric': False}
- coords:
Coordinates
, optional sky coordinates of the event
- satellite_skycoord: Astropy.coordinates.SkyCoord, optional
sky coordinates of the satellite specified by the ephemerides file. See
MulensModel.mulensdata.MulensData.satellite_skycoord
.- gamma: float, optional
limb darkening coefficient in gamma convention; defaults to 0
- Attributes :
- trajectory:
trajectory
Trajectory used to calculate positions of the source that are used to calculate magnification values.
- trajectory:
- set_magnification_methods(methods, default_method)¶
Sets methods used for magnification calculation.
- For available methods, see:
get_point_lens_magnification()
and
- Parameters :
- methods: list
List that specifies which methods (str) should be used when (float values for Julian dates). Given method will be used for times between the times between which it is on the list, e.g.,
methods = [ 2455746., 'Quadrupole', 2455746.6, 'Hexadecapole', 2455746.7, 'VBBL', 2455747., 'Hexadecapole', 2455747.15, 'Quadrupole', 2455748.]
- default_method: str
Name of the method to be used for epochs outside the ranges specified in methods.
- set_magnification_methods_parameters(methods_parameters)¶
Set additional parameters for magnification calculation methods.
- Parameters :
- methods_parameters: dict
Dictionary that for method names (keys) returns dictionary in the form of
**kwargs
that are passed to given method, e.g.,{'VBBL': {'accuracy': 0.005}}
.
- get_magnification()¶
Calculate magnification.
- Returns :
- magnification: np.ndarray
Vector of magnifications.
- get_point_lens_magnification()¶
Calculate the Point Lens magnification.
- Allowed magnification methods (set by
set_magnification_methods()
) : point_source
:standard Paczynski equation for a point source/point lens.
finite_source_uniform_Gould94
:Uses the Gould 1994 ApJ, 421L, 71 prescription assuming a uniform (and circular) source. This method interpolates pre-computed tables. The relative interpolation errors are smaller than 10^-4.
finite_source_uniform_Gould94_direct
:Same as
finite_source_uniform_Gould94
, but calculates the underlying functions directly (i.e., without interpolation).finite_source_uniform_WittMao94
:Uses the Witt and Mao 1994 ApJ, 430, 505 method assuming a uniform (and circular) source. This method interpolates pre-computed tables. The relative interpolation errors are smaller than 10^-4.
finite_source_LD_WittMao94
:Uses the Witt and Mao 1994 ApJ, 430, 505 method and integrates multiple annuli to obtain magnification for a circular source including limb-darkening. For description of integration of multiple annuli see, e.g., Bozza et al. 2018 MNRAS, 479, 5157. This method interpolates pre-computed tables. The relative interpolation errors are smaller than 10^-4.
finite_source_LD_Yoo04
:Uses the Yoo et al. 2004 ApJ, 603, 139 prescription for a circular source including limb-darkening. This method interpolates pre-computed tables. The relative interpolation errors are smaller than 10^-4.
finite_source_LD_Yoo04_direct
:Same as
finite_source_LD_Yoo04
, but calculates the underlying functions directly (i.e., without interpolation), hence can be slow.finite_source_uniform_Lee09
:Uses the Lee et al. 2009 ApJ, 695, 200 method for a circular and uniform source. This method works well for large sources (rho ~ 1).
finite_source_LD_Lee09
:Uses the Lee et al. 2009 ApJ, 695, 200 method for a circular source including limb-darkening. This method works well for large sources (rho ~ 1) but can be slow compared to other methods.
- Returns :
- magnification: np.ndarray
Vector of magnifications.
- Allowed magnification methods (set by
- get_binary_lens_magnification()¶
Calculate the binary lens magnification. If the shear or convergence are set, then they are used.
- Allowed magnification methods (set by
set_magnification_methods()
) : point_source
:standard point source magnification calculation.
quadrupole
:From Gould 2008 ApJ, 681, 1593. See
hexadecapole_magnification()
hexadecapole
:From Gould 2008 ApJ, 681, 1593 See
hexadecapole_magnification()
VBBL
:Uses VBBinaryLensing (a Stokes theorem/contour integration code) by Valerio Bozza (Bozza 2010 MNRAS, 408, 2188). See
vbbl_magnification()
Adaptive_Contouring
:Uses AdaptiveContouring (a Stokes theorem/contour integration code) by Martin Dominik (Dominik 2007 MNRAS, 377, 1679). See
adaptive_contouring_magnification()
Note that it doesn’t work if shear or convergence are set.
point_source_point_lens
:Uses point-source _point_-_lens_ approximation; useful when you consider binary lens but need magnification very far from the lens (e.g. at separation u = 100).
- Returns :
- magnification: np.ndarray
Vector of magnifications.
- Allowed magnification methods (set by
- get_d_A_d_params(parameters)¶
Calculate d A / d parameters for a point lens model.
- Parameters :
- parameters: list
List of the parameters to take derivatives with respect to.
- Returns :
- dA_dparam: dict
Keys are parameter names from parameters argument above. Values are the partial derivatives for that parameter evaluated at each epoch.
- get_d_A_d_rho()¶
Calculate d A / d rho for a point lens model.
No Inputs
- Returns :
- dA_drho: np.array
Values are the partial derivatives for rho evaluated at each data point.
- property methods_for_epochs¶
list
for each epoch, decide which methods should be used to calculate magnification, but don’t run the calculations
- property methods_indices¶
dict
Keys are the magnification methods. Values are a boolean index that indicate which epochs should be calculated with each method.
- class MulensModel.Model(parameters=None, coords=None, ra=None, dec=None, ephemerides_file=None)¶
Bases:
object
A Model for a microlensing event with the specified parameters.
- Arguments :
parameters: dictionary,
ModelParameters
coords
: str, astropy.SkyCoords, MulensModel.Coordinates, optionalSky coordinates of the event. If type is str, then it is assumed that the units are hour angle and degrees for RA and Dec, respectively.
- ra, dec: str, optional
Sky coordinates of the event.
- ephemerides_file: str, optional
Specify name of the file with satellite ephemerides. See
MulensData
for more details. Note that if you provide file name here, then it will affect all calculations for this model. In most cases, you want to combine ground-based and satellite data and in those cases set ephemerides_file for specificMulensData
instance to pass satellite information.
- Attributes :
- ephemerides_file: str
Name of file with satellite ephemerides.
Default values for parallax are all True. Use
parallax()
to turn different parallax effects ON/OFF. If using satellite parallax, you may also specify an ephemerides_file (seeMulensData
).Note that you can print an instance of Model, which shows you parameters in a nice way, e.g.,
model = Model(parameters={'t_0': 2456789.0, ....}) print(model)
This will provide information on parameter values, coordinates, methods used for magnification calculations, and limb-darkening coefficients.
- plot_magnification(times=None, t_range=None, t_start=None, t_stop=None, dt=None, n_epochs=None, subtract_2450000=False, subtract_2460000=False, satellite_skycoord=None, gamma=None, source_flux_ratio=None, **kwargs)¶
Plot the model magnification curve.
- Keywords :
see
plot_lc()
- gamma:
- satellite_skycoord:
- source_flux_ratio: float
If the model has two sources, source_flux_ratio is the ratio of source_flux_2 / source_flux_1
**kwargs
:any arguments accepted by
matplotlib.pyplot.plot()
.
- get_lc(times=None, t_range=None, t_start=None, t_stop=None, dt=None, n_epochs=None, source_flux=None, blend_flux=None, source_flux_ratio=None, gamma=None, bandpass=None)¶
Calculate model light curve in magnitudes.
- Keywords :
- times: [float, list, numpy.ndarray]
a list of times at which to plot the magnifications
t_range, t_start, t_stop, dt, n_epochs: see
set_times()
- source_flux: float or list
Explicitly specify the source flux(es) in a system where flux = 1 corresponds to
MulensModel.utils.MAG_ZEROPOINT
(= 22 mag). If the model has n_source > 1, source_flux may be specified as a list: one value for each source. Alternatively, if source_flux is specified as a float, source_flux_ratio should also be specified. Then, source_flux is taken to be the flux of the first source, and the other source fluxes are derived using source_flux_ratio.- blend_flux: float
Explicitly specify the blend flux in a system where flux = 1 corresponds to
MulensModel.utils.MAG_ZEROPOINT
(= 22 mag).- source_flux_ratio: float, Optional
If the model has two sources, source_flux_ratio is the ratio of source_flux_2 / source_flux_1.
- gamma, bandpass:
- Returns :
- magnitudes: numpy.ndarray
Magnitude values for each epoch.
- plot_lc(times=None, t_range=None, t_start=None, t_stop=None, dt=None, n_epochs=None, source_flux=None, blend_flux=None, source_flux_ratio=None, gamma=None, bandpass=None, subtract_2450000=False, subtract_2460000=False, phot_fmt='mag', **kwargs)¶
Plot the model light curve in magnitudes.
- Keywords :
- times: [float, list, numpy.ndarray]
a list of times at which to plot the magnifications
- t_range, t_start, t_stop, dt, n_epochs:
see
set_times()
- source_flux, blend_flux, source_flux_ratio:
see
get_lc()
- gamma, bandpass:
- subtract_2450000, subtract_2460000: boolean, optional
If True, subtracts 2450000 or 2460000 from the time axis to get more human-scale numbers. If using, make sure to also set the same settings for all other plotting calls (e.g.
plot_data()
)- phot_fmt: str
Specifies whether the photometry is plotted in magnitude or flux space. Accepts either ‘mag’ or ‘flux’. Default = ‘mag’.
**kwargs
:any arguments accepted by
matplotlib.pyplot.plot()
.
- plot_caustics(n_points=5000, epoch=None, **kwargs)¶
Plot the caustic structure. See
MulensModel.caustics.Caustics.plot()
for binary lenses. For a single lens it just marks (0, 0) point and the first two parameters are ignored.- Additional parameters :
- n_points: int, optional
The number of points to calculate along the caustic. Defaults to 5000.
- epoch: float, optional
Epoch for which separation s will be used. Important for models with orbital motion. Defaults to t_0_kep, which defaults to t_0.
**kwargs
:keywords accepted by
matplotlib.pyplot.scatter()
- update_caustics(epoch=None)¶
Updates
caustics
property for given epoch.- Parameters :
- epoch: float
For orbital motion models, epoch for which separation s is calculated to calculate
caustics
. Defaults to t_0_kep, which defaults to t_0.
- plot_trajectory(times=None, t_range=None, t_start=None, t_stop=None, dt=None, n_epochs=None, caustics=False, arrow=True, satellite_skycoord=None, arrow_kwargs=None, **kwargs)¶
Plot the source trajectory.
Keywords (all optional) :
- times, t_range, t_start, t_stop, dt, n_epochs:
May all be used to specify exactly when to plot the source trajectory. See also
plot_lc()
andset_times()
.- caustics: boolean
plot the caustic structure in addition to the source trajectory. default=False (off). For finer control of plotting features, e.g. color, use
plot_caustics()
instead.- arrow: boolean
Show the direction of the source motion. Default is True.
satellite_skycoord: astropy.SkyCoord or
MulensModel.satelliteskycoord.SatelliteSkyCoord
Allows the user to specify that the trajectory is calculated for a satellite. If astropy.SkyCoord object is provided, then these are satellite positions for all epochs. See also
get_satellite_coords()
- arrow_kwargs: dict
Kwargs that are passed to
pyplot.arrow()
. If no color is given here, then we use one specified in**kwargs
and if nothing is there, then we use black. The size of the arrow is determined based on limits of current axis. If those are not adequate, then change the size by specifying width keyword and maybe other as well. Note that arrow_kwargs are of dict type and are different than**kwargs
.**kwargs
Controls plotting features of the trajectory. It’s passed to
pyplot.plot()
.
Note that in order to have equal scaling of both axis (i.e., make circles look circular), you have to call appropriate pyplot command. This can be one of these commands:
pyplot.axis('equal') pyplot.axis('scaled') pyplot.axis('square') pyplot.gca().set_aspect('equal')
They have slightly different behavior.
- plot_source(times=None, **kwargs)¶
Plot source: circles of the radius rho at positions corresponding to source positions at times. When the rho is not defined, then X symbols are plotted.
- Parameters:
- times: float or np.ndarray
epochs for which source positions will be plotted
**kwargs
:Keyword arguments passed to matplotlib.Circle. Examples:
color='red'
,fill=False
,linewidth=3
,alpha=0.5
. When the rho is not defined, then keyword arguments are passed to matplotlib.plot.
Note that it is likely that with default axis scaling, the circles may be plotted as ellipses. To mitigate it, use:
plt.gca().set_aspect('equal')
orplt.axis('equal')
(the other possible options are'scaled'
or'square'
).
- get_trajectory(times, satellite_skycoord=None)¶
Get the trajectory of the source for the given set of times. For multi-source models, one trajectory per source is returned.
- Parameters :
- times: np.ndarray, list of floats, or float
Epochs for which source positions are requested.
- satellite_skycoord: astropy.SkyCoord
Allows the user to specify that the trajectory is calculated for a satellite. If astropy.SkyCoord object is provided, then these are satellite positions for all epochs. See also
get_satellite_coords()
- Returns :
- trajectories: :py:class:`~MulensModel.trajectory.Trajectory object or a list of them
Single object for single source model, a list otherwise.
- set_times(t_range=None, t_start=None, t_stop=None, dt=None, n_epochs=1000)¶
Return a list of times. If no keywords are specified, default is 1000 epochs from [\(t_0 - 1.5 * t_E\), \(t_0 + 1.5 * t_E\)] range.
For multi-source models, respectively, minimum and maximum of t_0_N values are used.
- Parameters (all optional):
- t_range: [list, tuple]
A range of times of the form [t_start, t_stop]
- t_start, t_stop: float
a start or stop time.
- dt: float
the interval spacing between successive points
- n_epochs: int
the number of epochs (evenly spaced)
- Returns :
- times: np.ndarray
Vector of epochs.
- set_magnification_methods(methods, source=None)¶
Sets methods used for magnification calculation. See
MagnificationCurve
for a list of implemented methods.- Parameters :
- methods: list
List that specifies which methods (str) should be used when (float values for Julian dates). Given method will be used for times between the times between which it is on the list, e.g.,
methods = [ 2455746., 'Quadrupole', 2455746.6, 'Hexadecapole', 2455746.7, 'VBBL', 2455747., 'Hexadecapole', 2455747.15, 'Quadrupole', 2455748.]
- source: int or None, optional
Which source do the given methods apply to? Accepts 1, 2, or None (i.e., all sources). Default is None
- get_magnification_methods(source=None)¶
Gets methods used for magnification calculation. See
set_magnification_methods()
- Parameters :
- source: int or None, optional
Which source do the given methods apply to? Accepts 1, 2, or None (i.e., all sources). Default is None.
- property methods¶
list of methods used for magnification calculation or dict of lists if there are multiple sources.
- property default_magnification_method¶
Stores information on method to be used, when no method is directly specified. See
MagnificationCurve
for a list of implemented methods.- Parameters:
- method: str
Name of the method to be used.
- set_magnification_methods_parameters(methods_parameters)¶
Set additional parameters for magnification calculation methods.
- Parameters :
- methods_parameters: dict
Dictionary that for method names (keys) returns dictionary in the form of
**kwargs
that are passed to given method, e.g.,{'VBBL': {'accuracy': 0.005}}
.
- get_magnification_methods_parameters(method)¶
Get additional parameters for a specific magnification calculation method or methods.
- Parameters :
- method: str, list
Name of method or a list of the names for which parameters will be returned.
- Returns :
- method_parameters: dict
- set_limb_coeff_gamma(bandpass, coeff)¶
Store gamma limb darkening coefficient for given band. See also
LimbDarkeningCoeffs
.- Parameters :
- bandpass: str
Bandpass for the coefficient you provide.
- coeff: float
Value of the coefficient.
- get_limb_coeff_gamma(bandpass)¶
Get gamma limb darkening coefficient for given band.
- Parameters :
- bandpass: str
Bandpass for which coefficient will be provided.
- Returns :
- gamma: float
limb darkening coefficient
- set_limb_coeff_u(bandpass, coeff)¶
Store u limb darkening coefficient for given band. See also
MulensModel.limbdarkeningcoeffs.LimbDarkeningCoeffs
.- Parameters :
- bandpass: str
Bandpass for which coefficient you provide.
- coeff: float
Value of the coefficient.
- get_limb_coeff_u(bandpass)¶
Get u limb darkening coefficient for given band.
- Parameters :
- bandpass: str
Bandpass for which coefficient will be provided.
- Returns :
- u: float
limb darkening coefficient
- parallax(earth_orbital=None, satellite=None, topocentric=None)¶
Specifies the types of the microlensing parallax that will be included in calculations.
Parameters :
- earth_orbital: boolean, optional
Do you want to include the effect of Earth motion about the Sun? Default is False.
- satellite: boolean, optional
Do you want to include the effect of difference due to separation between the Earth and satellite? Note that this separation changes over time. Default is False.
- topocentric: boolean, optional
Do you want to include the effect of different positions of observatories on the Earth? Default is False. Note that this is significant only for very high magnification events and if high quality datasets are analyzed. Hence, this effect is rarely needed. Not Implemented yet.
- get_parallax()¶
Returns dict that specifies the types of the microlensing parallax that are included in calculations.
- Returns :
- parallax: dict
For keys
'earth_orbital'
,'satellite'
, and'topocentric'
returns bool.
- get_satellite_coords(times)¶
Get astropy.SkyCoord object that gives satellite positions for given times. see also
MulensModel.satelliteskycoord.SatelliteSkyCoord
- Parameters :
- times: np.ndarray or list
Epochs for which satellite position is requested.
- Returns :
- satellite_skycoord: astropy.SkyCoord
SkyCoord giving satellite positions. The parameter representation is set to ‘spherical’. If ephemerides_file is not set, returns None.
- get_magnification(time, satellite_skycoord=None, gamma=None, bandpass=None, source_flux_ratio=None, separate=None)¶
Calculate the model magnification for the given time(s).
- Parameters :
- time: np.ndarray, list of floats, or float
Times for which magnification values are requested.
- satellite_skycoord: astropy.coordinates.SkyCoord, optional
SkyCoord object that gives satellite positions. Must be the same length as time parameter. Use only for satellite parallax calculations.
- gamma: float, optional
The limb-darkening coefficient in gamma convention. Default is 0 which means no limb darkening effect.
- bandpass: str, optional
The bandpass for setting the limb-darkening coefficient. Expects that you have used
set_limb_coeff_gamma()
orset_limb_coeff_u()
. Only ONE of ‘gamma’ or ‘bandpass’ may be specified.- source_flux_ratio: float, list
If the model has 2 sources, source_flux_ratio is the ratio of source_flux_2 / source_flux_1. For N sources, source_flux_ratio a list of the ratios of source_flux_i / source_flux_1 where i = 2, N.
- separate: boolean, optional
For multi source models, return magnification of each source separately. Defaults to True if source_flux_ratio is provided and False otherwise (then only effective magnification is returned).
- Returns :
- magnification: np.ndarray
A vector of calculated magnification values. For binary source models, the effective magnification is returned (unless separate=True).
- get_magnification_curve(time, satellite_skycoord, gamma)¶
Create a
MagnificationCurve
object for a given set of times.- Parameters :
- time: np.ndarray, list of floats, or float
Times for which magnification values are requested.
- satellite_skycoord: astropy.coordinates.SkyCoord, optional
SkyCoord object that gives satellite positions. Must be the same length as time parameter. Use only for satellite parallax calculations.
- gamma: float, optional
The limb-darkening coefficient in gamma convention. Default is 0 which means no limb darkening effect.
- Return:
py:class:~MulensModel.magnificationcurve.MagnificationCurve
- get_magnification_curves(time, satellite_skycoord, gamma)¶
Create a list of
MagnificationCurve
objects for multiple sources, given a set of times.- Parameters :
- time: np.ndarray, list of floats, or float
Times for which magnification values are requested.
- satellite_skycoord: astropy.coordinates.SkyCoord, optional
SkyCoord object that gives satellite positions. Must be the same length as time parameter. Use only for satellite parallax calculations.
- gamma: float, optional
The limb-darkening coefficient in gamma convention. Default is 0 which means no limb darkening effect.
- Return:
list of py:class:~MulensModel.magnificationcurve.MagnificationCurve
- property caustics¶
Caustics
Caustics for given model
- property parameters¶
-
Model parameters.
- property n_lenses¶
int
number of objects in the lens system
- property n_sources¶
int
number of luminous sources; it’s possible to be 1 for xallarap model
- is_static()¶
- property coords¶
see
Coordinates
- property bandpasses¶
list
List of all bandpasses for which limb darkening coefficients are set.
- class MulensModel.ModelParameters(parameters)¶
Bases:
object
A class for the basic microlensing model parameters (t_0, u_0, t_E, s, q, alpha, etc.). Can handle point lens or binary lens. The pi_E assumes NE coordinates (Parallel, Perpendicular coordinates are not supported).
- Arguments :
- parameters: dictionary
A dictionary of parameters and their values. See
which_parameters()
for valid parameter combinations.
- Attributes :
- parameters: dictionary
A dictionary of parameters and their values. Do not use it to change parameter values, instead use e.g.:
model_parameters.u_0 = 0.1
orsetattr(model_parameters, 'u_0', 0.1)
.
- Example:
- Define a point lens model:
params = ModelParameters({'t_0': 2450000., 'u_0': 0.3, 't_E': 35.})
- Then you can print the parameters:
print(params)
- property t_0¶
float
The time of minimum projected separation between the source and the lens center of mass.
- property u_0¶
float
The minimum projected separation between the source and the lens center of mass.
- property t_star¶
float
t_star = rho * t_E = source radius crossing time in days
- property t_eff¶
float
t_eff = u_0 * t_E = effective timescale in days
- property t_E¶
float
The Einstein timescale in days.
- property rho¶
float
source size as a fraction of the Einstein radius
- property alpha¶
float
The angle of the source trajectory relative to the binary lens axis (or primary-secondary axis). Measured counterclockwise, i.e., according to convention advocated by Skowron et al. 2011 (ApJ, 738, 87), but shifted by 180 deg.
- property q¶
float
mass ratio of the two lens components. Only 2 bodies allowed.
- property convergence_K¶
float
Convergence of external mass sheet.
- property shear_G¶
complex
Shear of external mass sheet.
- property s¶
float
separation of the two lens components relative to Einstein ring size
- property pi_E_N¶
float
The North component of the microlensing parallax vector.
- property pi_E_E¶
float
The East component of the microlensing parallax vector.
- property t_0_par¶
float
The reference time for the calculation of parallax. If not set explicitly, then it is assumed t_0_par = t_0. If there are multiple sources, t_0_1 is used.
Note that this is a reference value and not the fitting parameter. It is best to fix it at the begin of calculations.
- property pi_E_mag¶
float
The magnitude of the microlensing parallax vector.
- property x_caustic_in¶
float
Curvelinear coordinate (in Cassan (2008) parameterization) of caustic entrance for a static binary lens model. See
UniformCausticSampling
.
- property x_caustic_out¶
float
Curvelinear coordinate (in Cassan (2008) parameterization) of caustic exit for a static binary lens model. See
UniformCausticSampling
.
- property t_caustic_in¶
float
Epoch of caustic entrance for a static binary lens model in Cassan (2008) parameterization) See
UniformCausticSampling
.
- property t_caustic_out¶
float
Epoch of caustic exit for a static binary lens model in Cassan (2008) parameterization) See
UniformCausticSampling
.
- property t_0_kep¶
float
The reference time for the calculation of lens orbital motion. If not set explicitly, then it is assumed t_0_kep = t_0 (or t_0_1 for multi-source models).
Note that this is a reference value and not the fitting parameter. It is best to fix it at the begin of calculations.
- property xi_period¶
float
Orbital period of the source system (xallarap) in days.
- property xi_semimajor_axis¶
float
Semi-major axis of the source orbit (xallarap) in the theta_E units.
- property xi_Omega_node¶
float
The longitude of the ascending node of the xallarap orbit, i.e., the angle from relative lens-source proper motion direction to the ascending node direction. The units are degrees.
- property xi_inclination¶
float
The inclination of the xallarap orbit, i.e., the angle between source-orbit plane and the sky plane. The units are degrees.
- property xi_argument_of_latitude_reference¶
float
The argument of latitude for the xallarap orbit at
t_0_xi
. The argument of latitude is a sum of the true anomaly and the argument of periapsis. In standard notation: u = nu + omega. This parameter is internally used to calculate perihelion passage (T_0 in standard notation). The units are degrees.
- property xi_eccentricity¶
float
The eccentricity of the xallarap orbit. Has to be in [0, 1) range.
- property xi_omega_periapsis¶
float
The argument of periapsis of the xallarap orbit, i.e., the angle between the ascending node and periapsis measured in the direction of motion. The units are degrees.
- property t_0_xi¶
float
Reference epoch for xallarap orbit. If not provided, then it defaults to
t_0
.
- property q_source¶
float
The mass ratio of the second and the first source. This is value must be positive and can be > 1. Defined only for xallarap binary-source models because it does not affect the magnification for binary-source models without xallarap.
- property xallarap_reference_position¶
np.ndarray of shape (2, 1)
The position of the first source at
t_0_xi
relative to the source center of mass. It is a 2D vector that is subtracted from the source position along the orbit in order to calculate the shift caused by xallarap.
- property t_0_1¶
float
The time of minimum projected separation between the source no. 1 and the lens center of mass.
- property t_0_2¶
float
The time of minimum projected separation between the source no. 2 and the lens center of mass.
- property u_0_1¶
float
The minimum projected separation between the source no. 1 and the lens center of mass.
- property u_0_2¶
float
The minimum projected separation between the source no. 2 and the lens center of mass.
- property t_star_1¶
float
t_star_1 = rho_1 * t_E_1 = source no. 1 radius crossing time in days
- property t_star_2¶
float
t_star_2 = rho_2 * t_E_2 = source no. 2 radius crossing time in days.
- property rho_1¶
float
source no. 1 size as a fraction of the Einstein radius
- property rho_2¶
float
source no. 2 size as a fraction of the Einstein radius
- get_s(epoch)¶
Returns the value of separation
s
at a given epoch or epochs (if orbital motion parameters are set).- Arguments :
- epoch: float, list, np.ndarray
The time(s) at which to calculate
s
.
- Returns :
- separation: float or np.ndarray
Value(s) of separation for given epochs.
- get_alpha(epoch)¶
Returns the value of angle
alpha
at a given epoch or epochs (if orbital motion parameters are set).- Arguments :
- epoch: float, list, np.ndarray
The time(s) at which to calculate
alpha
.
- Returns :
- angle: float
Value(s) of angle for given epochs in degrees
- property gamma_parallel¶
float
Parallel component of instantaneous velocity of the secondary relative to the primary in 1/year. It is parallel to the primary-secondary axis. Equals
ds_dt
/s
. Cannot be set.
- property gamma_perp¶
float
Perpendicular component of instantaneous velocity of the secondary relative to the primary. It is perpendicular to the primary-secondary axis. It has sign opposite to
dalpha_dt
and is in rad/yr, not deg/yr. Cannot be set.
- property gamma¶
float
Instantaneous velocity of the secondary relative to the primary in 1/year. Cannot be set.
- is_finite_source()¶
Checks if model has finite source. For binary source models it checks if either of the sources is finite.
- Returns:
- is_finite_source: boolean
True if at least one source has finite size.
- is_static()¶
Checks if model is static, i.e., orbital motion parameters are not set.
- Returns :
- is_static: boolean
True if dalpha_dt or ds_dt are set.
- property n_lenses¶
int
Number of objects in the lens system.
- property n_sources¶
int
Number of luminous sources. It can be be 1 for a xallarap model.
- property is_external_mass_sheet¶
bool
Whether an external mass sheet is included in the model
- property is_external_mass_sheet_with_shear¶
bool
Whether an external mass sheet is included in the model with non-zero shear
- property is_xallarap¶
bool
Whether the parameters include the xallarap or not.
- property source_1_parameters¶
-
Parameters of source 1 in multi-source model.
Do not change returned values. To change parameters of the source 1, simply change the parameters of double source instance.
- property source_2_parameters¶
-
Parameters of source 2 in multi-source model.
Do not change returned values. To change parameters of the source 1, simply change the parameters of double source instance.
- property uniform_caustic_sampling¶
-
An instance of the class
UniformCausticSampling
that is used to calculate standard parameters based on the curvelinear coordinates. The main usage is access to the jacobian() function. In most cases, you do not need to access this property directly.
- as_dict()¶
Give parameters as a dict.
- Returns :
- dictionary: dict
The dictionary of model parameters.
- class MulensModel.MulensData(data_list=None, file_name=None, phot_fmt='mag', chi2_fmt='flux', ephemerides_file=None, add_2450000=False, add_2460000=False, bandpass=None, bad=None, good=None, plot_properties=None, **kwargs)¶
Bases:
object
A set of photometric measurements for a microlensing event.
- Examples of how to define a MulensData object:
data = MulensData(file_name=SAMPLE_FILE_01)
data = MulensData(data_list=[[Dates], [Magnitudes], [Errors]])
Parallax calculations assume that the dates supplied are BJD_TDB. See
Trajectory
. If you aren’t using parallax, the time system shouldn’t matter as long as it is consistent across all MulensData and Model objects. If you have multiple datasets, then you also need multiple instances of MulensData class.- Keywords :
- data_list: [list of lists, numpy.ndarray], optional
The list that contains three lists or numpy.ndarrays that specify: time, magnitude or flux, and its uncertainty (in that order). The lengths of these three objects must be the same.
- file_name: str, optional
The path to a file with columns: Date, Magnitude/Flux, Err. Loaded using
numpy.loadtxt()
. See**kwargs
.
Either data_list or file_name is required.
- phot_fmt: str
Specifies whether the photometry is provided in magnitude or flux space. Accepts either ‘mag’ or ‘flux’. Default = ‘mag’.
- chi2_fmt: str
Specifies whether the format used for chi^2 calculation should be done in Magnitude or Flux spaces. Accepts either ‘mag’ or ‘flux’. Default is ‘flux’ because almost always the errors are gaussian in flux space.
- ephemerides_file: str, optional
Specify the ephemerides of a satellite over the period when the data were taken. You may want to extend the time range to get nicer plots. Will be interpolated as necessary to model the satellite parallax effect. See instructions on getting satellite positions. Note that there is no check on time format (e.g., BJD TBD vs. HJD) and it should be the same as in data_list or file_name.
- add_2450000: boolean, optional
Adds 2450000 to the input dates. Useful if the dates are supplied as HJD-2450000.
- add_2460000: boolean, optional
Adds 2460000 to the input dates. Useful if the dates are supplied as HJD-2460000.
bandpass: see
bandpass
- bad: boolean np.ndarray, optional
Flags for bad data (data to exclude from fitting and plotting). Should be the same length as the number of data points.
- good: boolean np.ndarray, optional
Flags for good data, should be the same length as the number of data points.
- plot_properties: dict, optional
Specify properties for plotting, e.g.
color
,marker
,label
,alpha
,zorder
,markersize
,visible
, and also theshow_bad
andshow_errorbars
properties.Note: pyplot functions errorbar() and scatter() are used to plot data with errorbars and without them, respectively. The type and size of marker are specified using different keywords: (‘fmt’, ‘markersize’) for errorbar() and (‘marker’, ‘size’) for scatter(). You can use either convention in
plot_properties
and they will be translated to appropriate keywords. If there are similar problems with other keywords, then they won’t be translated unless you contact code authors.- Other special keys :
- show_errorbars: boolean, optional
Whether or not to show the errorbars for this dataset.
- show_bad: boolean, optional
Whether or not to plot data points flagged as bad.
**kwargs
:Kwargs passed to np.loadtxt(). Works only if
file_name
is set.
You can print an instance of this class, which always provides label and information on the total number of epochs and the number of bad epochs. If applicable, additional information is provided: bandpass, ephemerides file, color used for plotting, and errorbar scaling.
- property plot_properties¶
dict
Settings that specify how the photometry should be plotted.
The keys in this dict could be either special keys introduced here (i.e.,
show_bad
andshow_errorbars
) or keys accepted by matplotlib.pyplot plotting functions. The latter could be for examplecolor
,marker
,label
,alpha
,zorder
,markersize
, orvisible
.See
MulensData
for more information.
- plot(phot_fmt=None, show_errorbars=None, show_bad=None, subtract_2450000=False, subtract_2460000=False, **kwargs)¶
Plot the data.
Uses
plot_properties
for label, color, etc. This settings can be changed by setting**kwargs
.You can plot in either flux or magnitude space.
- Keywords:
- phot_fmt: string (‘mag’, ‘flux’)
Whether to plot the data in magnitudes or in flux. Default is the same as
input_fmt
.- show_errorbars: boolean
If show_errorbars is True (default), plots with matplotlib.errorbar(). If False, plots with matplotlib.scatter().
- show_bad: boolean
If False, bad data are suppressed (default). If True, shows points marked as bad (
mulensdata.MulensData.bad
) as ‘x’- subtract_2450000, subtract_2460000: boolean
If True, subtracts 2450000 or 2460000 from the time axis to get more human-scale numbers. If using it, make sure to also set the same settings for all other plotting calls (e.g.
plot_lc()
).**kwargs
:passed to matplotlib plotting functions.
- set_limb_darkening_weights(weights)¶
Save a dictionary of weights that will be used to evaluate the limb darkening coefficient. See also
LimbDarkeningCoeffs
- Parameters :
- weights: dict
A dictionary that specifies weight for each bandpass. Keys are str and values are float, e.g.,
{'I': 1.5, 'V': 1.}
if the I-band gamma limb-darkening coefficient is 1.5-times larger than the V-band.
- property time¶
np.ndarray
vector of dates
- property mag¶
np.ndarray
magnitude vector
- property err_mag¶
np.ndarray
vector of magnitude errors
- property flux¶
numpy.ndarray
Vector of the measured brightness in flux units.
- property err_flux¶
np.ndarray
Vector of uncertainties of flux values.
- property bad¶
np.ndarray boolean
flags marking bad data
- property n_epochs¶
int
give total number of epochs (including bad data)
- data_and_err_in_input_fmt()¶
Gives photometry in input format (mag or flux).
- Returns :
- data: np.ndarray
Magnitudes or fluxes
- data_err: np.ndarray
Uncertainties of magnitudes or of fluxes
- data_and_err_in_chi2_fmt()¶
Gives photometry in format used for chi2 calculation (flux in most cases, but magnitude possible).
- Returns :
- data: np.ndarray
Magnitudes or fluxes
- data_err: np.ndarray
Uncertainties of magnitudes or of fluxes
- property bandpass¶
String
Bandpass of given dataset (primary usage is limb darkening), e.g. ‘I’ or ‘V’. Returns None if not set.
- property satellite_skycoord¶
Astropy.SkyCoord object for satellite positions at epochs covered by the dataset
- Returns :
- skycoord: astropy.coordinates.SkyCoord
satellite positions at epochs covered by the dataset
- property input_fmt¶
str (‘mag’ or ‘flux’)
Input format - same as phot_fmt keyword in __init__().
- property chi2_fmt¶
str (‘mag’ or ‘flux’)
Photometry format used for chi^2 calculations. Default is ‘flux’.
- property ephemerides_file¶
str
File with satellite ephemeris.
- copy()¶
Returns a copy of given instance with settings copied
- Returns :
- mulens_data:
MulensData
Copy of self.
- mulens_data:
- scale_errorbars(factor=None, minimum=None)¶
Scale magnitude errorbars by multiplying by factor and then adding minimum in squares.
- Parameters :
- factor: float
Multiplication factor. Typically used for faint events.
- minimum: float
Floor of magnitude errorbars. Typically used for bright events.
- property errorbars_scale_factors¶
list of two floats
Parameters used by
scale_errorbars()
. The first one is multiplication factor, while the second is the minimum scatter added in squares.
- property errorbars_scaling_equation¶
str
Equation and parameters used by
scale_errorbars()
.
- class MulensModel.Lens(total_mass=None, mass=None, mass_1=None, mass_2=None, a_proj=None, distance=None, q=None, s=None, epsilon=None)¶
Bases:
object
A mass or system of masses.
Standard parameter combinations for defining a Lens object:
2+ body lens:
Note that s, q, and epsilon may be single values, lists, or numpy ndarrays.
If units are not specified for a given mass, it is assumed the value given is in Solar Masses.
If units are not specified for distance it is assumed the value is given in kpc.
- property epsilon¶
[float, list, numpy.ndarray]
An array of mass fractions for each lens components: m_i/total_mass. Stored as a numpy.ndarray.
- property q¶
[float, list, numpy.ndarray]
mass ratio for companions relative to the primary
Array of mass ratios defined relative to the primary (m_i/m_1). Size is number of components -1. Set as a list, np.ndarray, or single value.
Note: if
total_mass
is defined before q, it is assumed this is the mass of the primary. If you want this to actually be the total mass, define it after defining q.
- property s¶
[float, list, numpy.ndarray]
Separation between the components of the lens as a fraction of the Einstein ring. A np.ndarray or single value.
Definitions for more than 2 lens bodies TBD
- property total_mass¶
astropy.Quantity
The total mass of the lens (sum of all components). An astropy.Quantity with mass units. If set as a float, units are assumed to be solMass.
- property mass¶
astropy.Quantity
The mass of a point lens –> total mass. An astropy.Quantity with mass units. May be set as a float (in which case solMass is assumed).
- property mass_1¶
astropy.Quantity
The mass of the primary. Defined as total_mass * epsilon[0]. An astropy.Quantity with mass units. If set as a float, units are assumed to be solMass.
- property mass_2¶
astropy.Quantity
The mass of the secondary. Defined as total_mass * epsilon[1]. An astropy.Quantity with mass units. If set as a float, units are assumed to be solMass.
Note that if total_mass is defined before mass_2, and there is no epsilon corresponding to mass_2, mass_2 is added to the total_mass.
- property mass_3¶
astropy.Quantity
The mass of the tertiary. Defined as total_mass * epsilon[2]. An astropy.Quantity with mass units. If set as a float, units are assumed to be solMass.
Note that if total_mass is defined before mass_3, and there is no epsilon corresponding to mass_3, mass_3 is added to the total_mass.
- property n_masses¶
int, read-only
number of masses in the system.
- property distance¶
astropy.Quantity
The distance to the lens.
May be set as a float. If no unit is given, the value is assumed to be kpc.
- property pi_L¶
astropy.Quantity
The parallax to the lens in milliarcseconds. May be set as a float, in which case units are assumed to be milliarcseconds.
- property a_proj¶
astropy.Quantity
Projected separation between the components of the lens in AU. An astropy.Quantity with distance units. If set as float (without units), AU is assumed.
- property caustics¶
A
Caustics
object, read-only
- plot_caustics(n_points=5000, **kwargs)¶
A function to plot the x,y coordinates (scaled to the Einstein ring) of the caustics. Pyplot scatter is used for plotting. See
MulensModel.caustics.Caustics.plot()
.- Parameters :
- n_points: int
Number of points be plotted.
**kwargs
:Keyword arguments passed to Pyplot scatter
- class MulensModel.Source(distance=None, pi_S=None, angular_radius=None, limb_darkening=None)¶
Bases:
object
Physical properties of a source (background) star.
- Attributes :
limb_darkening:
LimbDarkeningCoeffs
Limb darkening coefficients of the source.
- property distance¶
astropy.Quantity
The distance to the source. May be set as a float. The distance should either be given in pc, or if no unit is given, the value is assumed to be kpc (u.kpc).
- property pi_S¶
astropy.Quantity
The parallax to the source in millarcseconds. May be set as a float. If no units are specified, assumes milliarcseconds (u.mas).
- property angular_radius¶
astropy.Quantity
Angular radius of the source. May be set as a float. If units are not specified, assumed to be microarcseconds (u.mas).
- class MulensModel.MulensSystem(lens=None, source=None, mu_rel=None)¶
Bases:
object
A microlensing system consisting of a lens and a source.
- property lens¶
A
Lens
object. Physical properties of the lens. Note: lens mass must be in solMasses.
- property mu_rel¶
astropy.Quantity
Relative proper motion between the source and lens stars. If set as a float, units are assumed to be mas/yr.
- property t_E¶
astropy.Quantity
The Einstein crossing time (in days). If set as a float, assumes units are in days.
- property pi_rel¶
astropy.Quantity, read-only
The source-lens relative parallax in milliarcseconds.
- property pi_E¶
float, read-only
The Einstein ring radius. It’s equal to pi_rel / theta_E. Dimensionless.
- property theta_E¶
astropy.Quantity, read-only
The angular Einstein Radius in milliarcseconds.
- property r_E¶
astropy.Quantity, read-only
The physical size of the Einstein Radius in the Lens plane (in AU).
- property r_E_tilde¶
astropy.Quantity, read-only
The physical size of the Einstein Radius projected onto the Observer plane (in AU).
- plot_magnification(u_0=None, alpha=None, **kwargs)¶
Plot the magnification curve for the lens. u_0 must always be specified. If the lens has more than one body, alpha must also be specified.
- Parameters :
- u_0: float
Impact parameter between the source and the lens (as a fraction of the Einstein ring)
- alpha: astropy.Quantity, float
If float then degrees are assumed as a unit. See
MulensModel.modelparameters.ModelParameters.alpha
**kwargs
:
- plot_caustics(n_points=5000, **kwargs)¶
Plot the caustics structure using Pyplot scatter. See
MulensModel.caustics.Caustics.plot()
- Parameters :
- n_points: int
Number of points be plotted.
**kwargs
:Keyword arguments passed to Pyplot scatter
- class MulensModel.PointSourcePointLensMagnification(trajectory)¶
Bases:
MulensModel.pointlens._PointLensMagnification
Equations for calculating point-source–point-lens magnification and its derivatives.
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- trajectory:
- get_magnification()¶
Calculate the magnification
Parameters : None
- Returns :
- magnification: float or np.ndarray
The magnification for each point specified by u in
trajectory
.
- get_d_A_d_u()¶
Calculate dA/du for PSPL point-source–point-lens model.
No parameters.
- Returns :
- dA_du: np.ndarray
Derivative dA/du evaluated at each epoch.
- class MulensModel.FiniteSourceUniformGould94Magnification(direct=False, **kwargs)¶
Bases:
MulensModel.pointlens._PointLensMagnification
Equations for calculating finite-source–point-lens magnification and its derivatives following the Gould 1994 ApJ, 421L, 71 prescription assuming a uniform (and circular) source.
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- direct: bool
Use direct calculation (slow) instead of interpolation. Default is False.
- trajectory:
- get_d_A_d_u()¶
Calculate dA/du for USPL uniform-source–point-lens model.
No parameters.
- Returns :
- dA_du: np.ndarray
Derivative dA/du evaluated at each epoch.
- get_magnification()¶
Calculate magnification for point lens and finite source (for a uniform source). The approximation was proposed by: Gould A. 1994 ApJ 421L, 71 “Proper motions of MACHOs” and later the integral calculation was simplified by: Yoo J. et al. 2004 ApJ 603, 139 “OGLE-2003-BLG-262: Finite-Source Effects from a Point-Mass Lens” This approach assumes rho is small (rho < 0.1). For larger sources use
FiniteSourceUniformLee09Magnification
.- Returns :
- magnification: float, np.array
The finite-source source magnification for each epoch.
- get_d_A_d_params(parameters)¶
Return the gradient of the magnification with respect to the FSPL parameters for each epoch.
- Parameters :
- parameters: list
List of the parameters to take derivatives with respect to.
- Returns:
- dA_dparams: dict
Keys are parameter names from parameters argument above. Values are the partial derivatives for that parameter evaluated at each epoch.
- get_d_A_d_rho()¶
Return the derivative of the magnification with respect to rho for each epoch.
No parameters.
- Returns :
- dA_drho: np.ndarray
Derivative dA/drho evaluated at each epoch.
- property z_¶
np.ndarray
Magnitude of lens-source separation scaled to rho for each epoch (i.e., z = u/rho).
- property b0¶
np.ndarray
The value of the B_0(z) function for each epoch.
- property db0¶
np.ndarray
Derivative of the B_0(z) function for each epoch.
- class MulensModel.FiniteSourceLDYoo04Magnification(gamma=None, **kwargs)¶
Bases:
MulensModel.pointlens.FiniteSourceUniformGould94Magnification
Equations for calculating finite-source–point-lens magnification and its derivatives following the Gould 1994 ApJ, 421L, 71 prescription assuming a limb-darkened (and circular) source.
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- direct: bool
Use direct calculation (slow) instead of interpolation. Default is False.
- gamma: float
The limb-darkening coefficient. See also
LimbDarkeningCoeffs
- trajectory:
- get_magnification()¶
Calculate magnification for point lens and finite source (for a uniform source). The approximation was proposed by: Gould A. 1994 ApJ 421L, 71 “Proper motions of MACHOs” and later the integral calculation was simplified by: Yoo J. et al. 2004 ApJ 603, 139 “OGLE-2003-BLG-262: Finite-Source Effects from a Point-Mass Lens” This approach assumes rho is small (rho < 0.1). For larger sources use
FiniteSourceLDLee09Magnification
.- Returns :
- magnification: float, np.array
The finite-source source magnification for each epoch.
- get_d_A_d_u()¶
Finite-source modification to derivatives of other parameters. = Factor due to rho for multiplying derivatives due to non-rho parameters.
- get_d_A_d_rho()¶
Return the derivative of the magnification with respect to rho for each epoch.
No parameters.
- Returns :
- dA_drho: np.ndarray
Derivative dA/drho evaluated at each epoch.
- property b1¶
np.ndarray
Value of the B_1(z) function for each epoch.
- property db1¶
np.ndarray
Derivative of the B_1(z) function for each epoch.
- property gamma¶
float
Limb-darkening gamma coefficient.
- class MulensModel.PointSourcePointLensWithShearMagnification(trajectory)¶
Bases:
MulensModel.pointlens._PointLensMagnification
Lens of point mass plus shear and convergence, i.e., Chang-Refsdal lens.
This lens was first considered by Chang and Refsdal (1979; Nature, 282, 561).
- Arguments :
- trajectory:
Trajectory
Including trajectory.parameters =
ModelParameters
- trajectory:
- get_magnification()¶
Calculate point source magnification for the lens composed of a single mass plus a mass sheet.
- Returns :
- pspl_magnification: np.ndarray
The point-source–point-lens magnification for each point specified by trajectory.
- class MulensModel.B0B1Utils¶
Bases:
object
Data and methods used for interpolation for finite-source point-lens magnification calculations.
There are no parameters of __init__(). In general, this class is not directly used by a user.
- interpolate_B0(z)¶
Interpolate B_0(z) function.
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- B0: np.ndarray
Values of B_0.
- interpolate_B0minusB1(z)¶
Interpolate B_0(z)-B_1(z) function.
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- B0_minus_B1: np.ndarray
Values of B_0(z)-B_1(z).
- interpolate_B1(z)¶
Interpolate B_1(z) function.
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- B1: np.ndarray
Values of B_1(z).
- interpolate_B0prime(z)¶
Interpolate derivative of B_0(z), i.e., d B_0(z)/d z.
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- dB0_dz: np.ndarray
Values of d B_0(z)/d z.
- interpolate_B1prime(z)¶
Interpolate derivative of B_1(z), i.e., d B_1(z)/d z.
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- dB1_dz: np.ndarray
Values of d B_1(z)/d z.
- get_interpolation_mask(z)¶
Get mask of z values for which interpolation of B_0(z), B_1(z) etc. is allowed
- Parameters :
- z: np.ndarray
Values of u/rho.
- Returns :
- mask: np.ndarray
Mask indicating for which z values interpolation is allowed.
- property z_max_interpolation¶
Maximum value of z for which interpolation is allowed.
float
- class MulensModel.EllipUtils¶
Bases:
object
Data and methods used for interpolation for finite-source point-lens magnification calculations.
There are no parameters of __init__(). In general, this class is not directly used by a user.
- class MulensModel.SatelliteSkyCoord(ephemerides_file, satellite=None)¶
Bases:
object
An object that gives the Astropy.SkyCoord of satellite for a given epoch based on an ephemerides file.
- Keywords :
- ephemerides_file: str
path to file with satellite ephemerides from JPL horizons, for examples see data/ephemeris_files/Spitzer_ephemeris_01.dat or data/ephemeris_files/K2_ephemeris_01.dat
- satellite: str, optional
Just the name of the satellite.
- Attributes :
- satellite: str
name of the satellite
- get_satellite_coords(times)¶
Calculate the coordinates of the satellite for given times using cubic interpolation. The epheris provided is extrapolated up to 3 minutes beyond the range provided.
- Parameters :
- times: np.ndarray or list of floats
Epochs for which satellite coordinates will be calculated.
- Returns :
- satellite_skycoord: Astropy.coordinates.SkyCoord
SkyCoord for satellite at epochs times.
- class MulensModel.Trajectory(times=None, parameters=None, x=None, y=None, parallax=None, coords=None, satellite_skycoord=None, earth_coords=None)¶
Bases:
object
The (dimensionless) X, Y trajectory of the source in the source plane. This class includes internal functions that calculate how microlensing parallax affects the trajectory.
For binary lens, the origin of the coordinate system is at the center of mass with higher mass (assuming q < 1) at negative X and Y=0.
This class follows the conventions defined in Appendix A of Skowron et al. (2011) except the definition of alpha, which is shifted by 180 deg.
- Arguments :
- times: [float, list, np.ndarray]
the times at which to generate the source trajectory, e.g. a vector. Either times OR (x AND y) are required.
parameters: instance of
ModelParameters
, requireda ModelParameters object specifying the microlensing parameters
- x, y: [float, list, np.ndarray]
Dimensionless X, Y coordinates of trajectory. Either times OR (x AND y) are required.
- parallax: boolean dictionary, optional
specifies what parallax effects should be used. Default is False for each of ‘earth_orbital’, ‘satellite’, and ‘topocentric’. (differs from
Model
which defaults to True)
coords: str, or
Coordinates
, Astropy.coordinates.SkyCoord, optionalsky coordinates of the event; required for parallax calculations
- satellite_skycoord: Astropy.coordinates.SkyCoord, optional
sky coordinates of the satellite specified by the ephemerides file. See
satellite_skycoord
.
- Attributes :
- parameters:
ModelParameters
input
ModelParameters
- satellite_skycoord: Astropy.coordinates.SkyCoord
input
satellite_skycoord
- parallax: dict
specifies which types of microlensing parallax will be taken into account; boolean dict with keys:
earth_orbital
,satellite
, andtopocentric
(default values are False)- coords:
Coordinates
event coordinates
- parameters:
- property x¶
np.ndarray
Dimensionless X coordinates of trajectory.
- property y¶
np.ndarray
Dimensionless Y coordinates of trajectory.
- property times¶
np.ndarray or None
Epochs for which trajectory is calculated. None if the class was initialized without epochs provided.
- property parallax_delta_N_E¶
dict
Net North (key=’N’) and East (key=’E’) components of the parallax offset calculated for each time stamp (so sum of the offsets from all parallax types).
- property d_perp¶
np.ndarray
Projection of the Earth-Satellite separation vector (in AU) on the sky.
- class MulensModel.UniformCausticSampling(s, q, n_points=10000)¶
Bases:
object
Uniform sampling of a binary lens caustic. Note that calculations take some time for given (s, q). Keep that in mind, when optimizing your fitting routine.
- Arguments :
- s: float
Separation of the two lens components relative to the Einstein ring size.
- q: float
Mass ratio of the two lens components.
- n_points: int
Number of points used for internal integration. Default value should work fine.
Instead of standard parameters (t_0, u_0, t_E, alpha), here we use four other parameters: two epochs of caustic crossing (t_caustic_in, t_caustic_out) and two curvelinear coordinates of caustic crossing (x_caustic_in, x_caustic_out).
The curvelinear coordinate, x_caustic, is defined so that going from 0 to 1 draws all caustics for given separation and mass ratio. We use 0-1 range, which is a different convention than in the papers cited below (we also use different symbols for epochs of caustic crossing and curvelinear coordinates). For a wide topology (i.e., 2 caustics), there is a value between 0 and 1 (called
x_caustic_sep
) which separates the caustics and a trajectory exists only if x_caustic_in and x_caustic_out correspond to the same caustic, i.e., both are smaller thanx_caustic_sep
or both are larger thanx_caustic_sep
. For a close topology (i.e., 3 caustics), there are two such separating values.For description of the curvelinear coordinates, see:
Cassan A. et al. 2010 A&A 515, 52 “Bayesian analysis of caustic-crossing microlensing events”
In order to visualize the curvelinear coordinates, you can run a code like:
import matplotlib.pyplot as plt import numpy as np sampling = UniformCausticSampling(s=1.1, q=0.3) color = np.linspace(0., 1., 200) points = [sampling.caustic_point(c) for c in color] x = [p.real for p in points] y = [p.imag for p in points] plt.scatter(x, y, c=color) plt.axis('equal') plt.colorbar() plt.show()
This will show an intermediate topology. Change s=1.1 to s=2. to plot a wide topology, or to s=0.7 to plot a close topology.
To be specific, the central caustics are plotted counter-clockwise and x_caustic=0. corresponds to right-hand point where the caustic crosses the X-axis. For a wide topology, the planetary caustic is plotted in a similar way. For a close topology, the lower planetary caustic is plotted counter-clockwise and the upper planetary caustic is symmetric, thus plotted clockwise. For planetary caustics in a close topology, the zero-point of x_caustic values is defined in a very complicated way, however it is a smooth function of s and q.
For more advanced fitting of binary lens events see:
- get_standard_parameters(x_caustic_in, x_caustic_out, t_caustic_in, t_caustic_out)¶
Get standard binary lens parameters (i.e.,
t_0
,u_0
,t_E
,alpha
; seeModelParameters
) based on provided curvelinear parameters.Note that this function quite frequently raises
ValueError
exception. This is because not all (s
,q
,x_caustic_in
, andx_caustic_out
) correspond to real trajectories. The returned values are in conventions used byModel
.- Keywords :
- x_caustic_in: float
Curvelinear coordinate of caustic entrance. Must be in (0, 1) range.
- x_caustic_out: float
Curvelinear coordinate of caustic exit. Must be in (0, 1) range.
- t_caustic_in: float
Epoch of caustic entrance.
- t_caustic_out: float
Epoch of caustic exit.
- Returns :
- parameters: dict
Dictionary with standard binary parameters, i.e, keys are
t_0
,u_0
,t_E
, andalpha
.
- get_x_in_x_out(u_0, alpha)¶
Calculate where given trajectory crosses the caustic.
- Parameters :
- u_0: float
The parameter u_0 of source trajectory, i.e., impact parameter.
- alpha: float
Angle defining the source trajectory.
- Returns :
- x_caustic_points: list of float
Caustic coordinates of points where given trajectory crosses the caustic. The length is 0, 2, 4, or 6. Note that if there are 4 or 6 points, then only some pairs will produce real trajectories.
- get_uniform_sampling(n_points, n_min_for_caustic=10, caustic=None)¶
Sample uniformly (x_caustic_in, x_caustic_out) space according to Jacobian and requirement that x_caustic_in corresponds to caustic entrance, and x_caustic_out corresponds to caustic exit. The Jacobian is defined by Eq. 23 of:
Cassan A. et al. 2010 A&A 515, 52 “Bayesian analysis of caustic-crossing microlensing events”
and the above requirement is defined under Eq. 27 of that paper.
Relative number of points per caustic is not yet specified. Points do not repeat. This function is useful for sampling starting distribution for model fitting. For example sampling see Cassan et al. (2010) bottom panel of Fig. 1.
- Parameters :
- n_points: int
number of points to be returned
- n_min_for_caustic: int
minimum number of points in each caustic
- caustic: int or None
Select which caustic will be sampled. None means all caustics. Can be 1, 2, or 3 but has to be <=
n_caustics
.
- Returns :
- x_caustic_in: np.ndarray
Randomly drawn entrance points.
- x_caustic_out: np.ndarray
Corresponding randomly drawn exit points.
- jacobian(x_caustic_in, x_caustic_out)¶
Evaluates Eq. 23 from Cassan et al. (2010) with condition under Eq. 27.
- Parameters :
- x_caustic_in: float
Point of caustic entrance.
- x_caustic_out: float
Point of caustic exit.
- Returns :
- jacobian: float
Value of Jacobian. Returns 0. if trajectory does not exist.
- check_valid_trajectory(x_caustic_in, x_caustic_out)¶
Check if given (x_caustic_in, x_caustic_out) define an existing trajectory. An obvious case, when they don’t is when both caustic points are on the same fold, but other cases exists.
- Parameters :
- x_caustic_in: float
Coordinate of putative caustic entrance.
- x_caustic_out: float
Coordinate of putative caustic exit.
- Returns :
- check: bool
True if input defines a trajectory, False if it does not.
- caustic_point(x_caustic)¶
Calculate caustic position corresponding to given x_caustic.
- Keywords :
- x_caustic: float
Curvelinear coordinate of the point considered. Has to be in 0-1 range.
- Returns :
- point: numpy.complex128
Caustic point in complex coordinates.
- which_caustic(x_caustic)¶
Indicates on which caustic given point is.
- Keywords :
- x_caustic: float
Curvelinear coordinate to be checked
- Returns :
- i_caustic: int
Number indicating the caustic:
1
- central caustic,2
- planetary caustic; for close configuration it is the lower of the two planetary caustics,3
- upper planetary caustic.
- property n_caustics¶
int
Number of caustics: 1 for resonant topology, 2 for wide topology, or 3 for close topology.
- property s¶
float
separation of the two lens components relative to Einstein ring size
- property q¶
float
Mass ratio.
- class MulensModel.Utils¶
Bases:
object
A number of small functions used in different places
- static get_flux_from_mag(mag, zeropoint=None)¶
Transform magnitudes into fluxes.
- Parameters :
- mag: np.ndarray or float
Values to be transformed.
- zeropoint: float, optional
Zeropoint of magnitude scale. Defaults to 22. - double check if you want to change this.
- Returns :
- flux: np.ndarray or float
Calculated fluxes. Type is the same as mag parameter.
- static get_flux_and_err_from_mag(mag, err_mag, zeropoint=None)¶
Transform magnitudes and their uncertainties into flux space.
- Parameters :
- mag: np.ndarray or float
Magnitude values to be transformed.
- err_mag: np.ndarray or float
Uncertainties of magnitudes to be transformed.
- zeropoint: float
Zeropoint of magnitude scale. Defaults to 22. - double check if you want to change this.
- Returns :
- flux: np.ndarray or float
Calculated fluxes. Type is the same as mag parameter.
- err_flux: np.ndarray or float
Calculated flux uncertainties. Type is float if both mag and err_mag are floats and np.ndarray otherwise.
- static get_mag_from_flux(flux, zeropoint=None)¶
Transform fluxes into magnitudes.
- Parameters :
- flux: np.ndarray or float
Values to be transformed.
- zeropoint: float
Zeropoint of magnitude scale. Defaults to 22. - double check if you want to change this.
- Returns :
- mag: np.ndarray or float
Calculated fluxes. Type is the same as flux parameter.
- static get_mag_and_err_from_flux(flux, err_flux, zeropoint=None)¶
Transform fluxes and their uncertainties into magnitude space.
- Parameters :
- flux: np.ndarray or float
Flux values to be transformed.
- err_flux: np.ndarray or float
Uncertainties of fluxes to be transformed.
- zeropoint: float
Zeropoint of magnitude scale. Defaults to 22. - double check if you want to change this.
- Returns :
- mag: np.ndarray or float
Calculated fluxes. Type is the same as mag parameter.
- err_mag: np.ndarray or float
Calculated magnitude uncertainties. Type is float if both flux and err_flux are floats and np.ndarray otherwise.
- static gamma_to_u(gamma)¶
Transform gamma limb darkening coefficient to u in convention introduced by An et al. 2008 (ApJ 681, 1593).
- Parameters :
- gamma: float
Limb darkening coefficient in gamma convention.
- Returns :
- u: float
Limb darkening coefficient in u convention.
- static u_to_gamma(u)¶
Transform u limb darkening coefficient to gamma in convention introduced by An et al. 2008 (ApJ 681, 1593).
- Parameters :
- u: float
Limb darkening coefficient in u convention.
- Returns :
- gamma: float
Limb darkening coefficient in gamma convention.
- static get_n_caustics(s, q)¶
Find number of caustics for binary lens.
- Parameters :
- s: float
Separation of binary lens relative to theta_E.
- q: float
Mass ratio (<1).
- Returns :
- n_caustics: int
Number of caustics: 1, 2, or 3.
- static velocity_of_Earth(full_BJD)¶
Calculate 3D velocity of Earth for given epoch.
If you need velocity projected on the plane of the sky, then use
v_Earth_projected()
- Parameters :
- full_BJD: float
Barycentric Julian Data. Full means it should begin with 245… or 246…
- Returns :
- velocity: np.ndarray (float, size of (3,))
3D velocity in km/s. The frame follows Astropy conventions.
- static complex_fsum(arguments)¶
Accurate floating points sum of complex numbers in iterable.
- Parameters :
arguments: iterable (e.g., list or np.ndarray)
- Returns :
- sum: complex
Sum of arguments.
- static vector_product_normalized(vector_1, vector_2)¶
Get vector that is perpendicular to the 2 input ones and is normalized.
- Parameters :
- vector_1: np.ndarray (3,)
First vector.
- vector_2: np.ndarray (3,)
Second vector.
- Returns :
- vector_product_norm: np.ndarray (3,)
Normalized vector product of the inputs.
- static dot(cartesian, vector)¶
Dot product of 2 vectors represented by astropy.CartesianRepresentation and np.ndarray.
- Parameters :
- cartesian: astropy.CartesianRepresentation
First vector.
- vector: np.ndarray
Second vector (length 3).
- Returns :
- dot_product: astropy.Quantity
Dot product of inputs.
- static date_change(text)¶
Change format of month in date, e.g. ‘2015-Oct-30 12:00’ -> ‘2015-10-30 12:00’
- Parameters :
- text: str
Date to be changed.
- Returns :
- out_text: str
Changed text.
- static astropy_version_check(minimum)¶
Check if astropy package is installed at given or later version.
- Parameters :
- minimum: str
Minimum version, e.g., ‘3.1.2’.
- Returns :
- out: bool
Is the installed version later?