MulensModel.modelparameters module

MulensModel.modelparameters.which_parameters(*args)

Prints information on valid parameter combinations that can be used to define a model or a particular effect. May be called with no arguments (returns information on many types of models) or with one argument referring to a specific model (e.g. PSPL) or effect (e.g. parallax).

Valid arguments: str

Model types: ‘PSPL’, ‘FSPL’, ‘PSBL’, ‘FSBL’

Effects: ‘point lens’, ‘binary lens’, ‘finite source’, ‘parallax’, ‘lens orbital motion’

class MulensModel.modelparameters.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 or setattr(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

“day” is the default unit. Can be set as float or astropy.Quantity, but always returns float in units of days.

property t_eff

float

t_eff = u_0 * t_E = effective timescale

“day” is the default unit. Can be set as float or astropy.Quantity, but always returns float in units of days.

property t_E

float

The Einstein timescale. “day” is the default unit. Can be set as float or astropy.Quantity, but always returns float in units of days.

property rho

float

source size as a fraction of the Einstein radius

property alpha

astropy.Quantity

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. May be set as a float –> assumes “deg” is the default unit. Regardless of input value, returns value in degrees.

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

list of floats

The microlensing parallax vector. Must be set as a vector/list (i.e. [pi_E_N, pi_E_E]). To get the magnitude of pi_E, use pi_E_mag

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.

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 ds_dt

astropy.Quantity

Change rate of separation s in 1/year. Can be set as AstroPy.Quantity or as float (1/year is assumed default unit). Regardless of input value, returns value in 1/year.

property dalpha_dt

astropy.Quantity

Change rate of angle alpha in deg/year. Can be set as AstroPy.Quantity or as float (deg/year is assumed default unit). Regardless of input value, returns value in deg/year.

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.

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 xallrap 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 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

“day” is the default unit. Can be set as float or astropy.Quantity, but always returns float in units of days.

property t_star_2

float

t_star_2 = rho_2 * t_E_2 = source no. 2 radius crossing time

“day” is the default unit. Can be set as float or astropy.Quantity, but always returns float in units of 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 :
separation: astropy.Quantity

Value(s) of angle for given epochs in degrees

property gamma_parallel

astropy.Quantity

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

astropy.Quantity

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

astropy.Quantity

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

ModelParameters

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

ModelParameters

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

UniformCausticSampling

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.