MulensModel.mulensobjects package¶
Submodules¶
Module contents¶
- class MulensModel.mulensobjects.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.mulensobjects.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.mulensobjects.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