EXOSIMS.PlanetPopulation package

Submodules

EXOSIMS.PlanetPopulation.AlbedoByRadius module

class EXOSIMS.PlanetPopulation.AlbedoByRadius.AlbedoByRadius(SAG13coeffs=[[0.38, -0.19, 0.26, 0.0], [0.73, -1.18, 0.59, 3.4]], SAG13starMass=1.0, Rprange=[0.6666666666666666, 17.0859375], arange=[0.09084645, 1.45354324], ps=[0.2, 0.5], Rb=[1.4], **specs)[source]

Bases: SAG13

Planet Population module based on SAG13 occurrence rates.

NOTE: This assigns constant albedo based on radius ranges.

SAG13coeffs

Coefficients used by the SAG13 broken power law. The 4 lines correspond to Gamma, alpha, beta, and the minimum radius.

Type:

float 4x2 ndarray

Gamma

Gamma coefficients used by SAG13 broken power law.

Type:

float ndarray

alpha

Alpha coefficients used by SAG13 broken power law.

Type:

float ndarray

beta

Beta coefficients used by SAG13 broken power law.

Type:

float ndarray

Rplim

Minimum radius used by SAG13 broken power law.

Type:

float ndarray

SAG13starMass

Assumed stellar mass corresponding to the given set of coefficients.

Type:

astropy Quantity

mu

Gravitational parameter associated with SAG13starMass.

Type:

astropy Quantity

Ca

Constants used for sampling.

Type:

float 2x1 ndarray

ps

Constant geometric albedo values.

Type:

float nx1 ndarray

Rb

Planetary radius break points for albedos in earthRad.

Type:

float (n-1)x1 ndarray

Rbs

Planetary radius break points with 0 padded on left and np.inf padded on right

Type:

float (n+1)x1 ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and planetary radius are jointly distributed. Eccentricity is a Rayleigh distribution. Albedo is a constant value based on planetary radius.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

get_p_from_Rp(Rp)[source]

Generate constant albedos for radius ranges

Parameters:

Rp (astropy Quantity array) – Planetary radius with units of earthRad

Returns:

Albedo values

Return type:

float ndarray

EXOSIMS.PlanetPopulation.AlbedoByRadiusDulzPlavchan module

class EXOSIMS.PlanetPopulation.AlbedoByRadiusDulzPlavchan.AlbedoByRadiusDulzPlavchan(starMass=1.0, occDataPath=None, esigma=0.13962979814050144, ps=[0.2, 0.5], Rb=[1.4], **specs)[source]

Bases: DulzPlavchan

Planet Population module based on occurrence rate tables from Shannon Dulz and Peter Plavchan.

NOTE: This assigns constant albedo based on radius ranges.

SAG13coeffs

Coefficients used by the SAG13 broken power law. The 4 lines correspond to Gamma, alpha, beta, and the minimum radius.

Type:

float 4x2 ndarray

Gamma

Gamma coefficients used by SAG13 broken power law.

Type:

float ndarray

alpha

Alpha coefficients used by SAG13 broken power law.

Type:

float ndarray

beta

Beta coefficients used by SAG13 broken power law.

Type:

float ndarray

Rplim

Minimum radius used by SAG13 broken power law.

Type:

float ndarray

SAG13starMass

Assumed stellar mass corresponding to the given set of coefficients.

Type:

astropy Quantity

mu

Gravitational parameter associated with SAG13starMass.

Type:

astropy Quantity

Ca

Constants used for sampling.

Type:

float 2x1 ndarray

ps

Constant geometric albedo values.

Type:

float nx1 ndarray

Rb

Planetary radius break points for albedos in earthRad.

Type:

float (n-1)x1 ndarray

Rbs

Planetary radius break points with 0 padded on left and np.inf padded on right

Type:

float (n+1)x1 ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and planetary radius are jointly distributed. Eccentricity is a Rayleigh distribution. Albedo is a constant value based on planetary radius.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

get_p_from_Rp(Rp)[source]

Generate constant albedos for radius ranges

Parameters:

Rp (astropy Quantity array) – Planetary radius with units of earthRad

Returns:

Albedo values

Return type:

float ndarray

EXOSIMS.PlanetPopulation.Brown2005EarthLike module

class EXOSIMS.PlanetPopulation.Brown2005EarthLike.Brown2005EarthLike(eta=1, arange=[0.6377303505401009, 1.3665650368716449], erange=[0.0, 0.35], constrainOrbits=False, **specs)[source]

Bases: PlanetPopulation

Population of Earth-Like Planets from Brown 2005 paper

This implementation is intended to enforce this population regardless of JSON inputs. The only inputs that will not be disregarded are erange and constrainOrbits.

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and eccentricity are uniformly distributed with all other parameters constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

gen_radius_nonorm(n)[source]

Generate planetary radius values in Earth radius. This one just generates a bunch of EarthRad

Parameters:

n (integer) – Number of target systems. Total number of samples generated will be, on average, n*self.eta

Returns:

Planet radius values in units of Earth radius

Return type:

astropy Quantity array

EXOSIMS.PlanetPopulation.DulzPlavchan module

class EXOSIMS.PlanetPopulation.DulzPlavchan.DulzPlavchan(starMass=1.0, occDataPath=None, esigma=0.13962979814050144, **specs)[source]

Bases: PlanetPopulation

Population based on occurrence rate tables from Shannon Dulz and Peter Plavchan.

The data comes as either Period-Radius or semi-major axis-mass pairs. If occDataPath is not specified, the nominal Period-Radius table is loaded.

Parameters:

specs – user specified values

starMass

stellar mass in M_sun used to convert period to semi-major axis

occDataPath

path on local disk to occurrence rate table

esigma

Sigma value of Rayleigh distribution for eccentricity.

Type:

float

Notes: 1. Mass/Radius and semi-major axis are specified in occurrence rate tables. User specified values will be ignored. 2. Albedo is sampled as in KeplerLike1 and KeplerLike2. 3. Eccentricity is Rayleigh distributed with user defined sigma parameter.

MfromRp(Rp)[source]

Converts mass to radius using Chen and Kipping

Parameters:

M (astropy Quantity array) – Planet mass in units of Earth mass

Returns:

Planet radius in units of Earth radius

Return type:

astropy Quantity array

RpfromM(M)[source]

Converts mass to radius using Chen and Kipping

Parameters:

M (astropy Quantity array) – Planet mass in units of Earth mass

Returns:

Planet radius in units of Earth radius

Return type:

Rp (astropy Quantity array)

dist_albedo(p)[source]

Probability density function for albedo

Parameters:

p (float ndarray) – Albedo value(s)

Returns:

Albedo probability density

Return type:

float ndarray

dist_eccen(e)[source]

Probability density function for eccentricity

Parameters:

e (float ndarray) – Eccentricity value(s)

Returns:

Eccentricity probability density

Return type:

float ndarray

dist_eccen_from_sma(e, a)[source]

Probability density function for eccentricity constrained by semi-major axis, such that orbital radius always falls within the provided sma range.

This provides a Rayleigh distribution between the minimum and maximum allowable values.

Parameters:
  • e (float ndarray) – Eccentricity values

  • a (float ndarray) – Semi-major axis value in AU. Not an astropy quantity.

Returns:

Probability density of eccentricity constrained by semi-major axis

Return type:

float ndarray

dist_radius(Rp)[source]

Probability density function for planetary radius.

Note that this is a marginalized distribution.

Parameters:

Rp (float ndarray) – Planetary radius value(s)

Returns:

Planetary radius probability density

Return type:

float ndarray

dist_sma(a)[source]

Probability density function for semi-major axis.

Note that this is a marginalized distribution.

Parameters:

a (float ndarray) – Semi-major axis value(s)

Returns:

Semi-major axis probability density

Return type:

float ndarray

gen_albedo(n)[source]

Generate geometric albedo values

The albedo is determined by sampling the semi-major axis distribution, and then calculating the albedo from the physical model.

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet albedo values

Return type:

float ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and planetary radius come from the occurrence rate tables and are assumed to be log-uniformly distributed within bins. Eccentricity is a Rayleigh distribution. Albedo is dependent on the PlanetPhysicalModel but is calculated such that it is independent of other parameters.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

gen_sma_radius(n)[source]

Generates semi-major axis and planetary radius samples

Parameters:

n (int) – number of samples

Returns:

a (astropy Quantity array):

Semi-major axis samples in units of AU

Rp (astropy Quantity array):

Planetary radius samples in units of Earth radius

Return type:

tuple

EXOSIMS.PlanetPopulation.EarthTwinHabZone1 module

class EXOSIMS.PlanetPopulation.EarthTwinHabZone1.EarthTwinHabZone1(eta=0.1, **specs)[source]

Bases: PlanetPopulation

Population of Earth twins (1 R_Earth, 1 M_Eearth, 1 p_Earth) On circular Habitable zone orbits (0.7 to 1.5 AU)

Note that these values may not be overwritten by user inputs. This implementation is intended to enforce this population regardless of JSON inputs.

dist_sma(a)[source]

Probability density function for uniform semi-major axis distribution in AU

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

float ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis is uniformly distributed and all other parameters are constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

EXOSIMS.PlanetPopulation.EarthTwinHabZone1SDET module

class EXOSIMS.PlanetPopulation.EarthTwinHabZone1SDET.EarthTwinHabZone1SDET(eta=0.1, **specs)[source]

Bases: PlanetPopulation

Population of Earth twins (1 R_Earth, 1 M_Eearth, 1 p_Earth) On circular Habitable zone orbits (0.7 to 1.5 AU)

Warning

Note that these values may not be overwritten by user inputs. This implementation is intended to enforce this population regardless of the Input Specification.

Parameters:
dist_sma(a)[source]

Probability density function for uniform semi-major axis distribution in AU

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

float numpy.ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis is uniformly distributed and all other parameters are constant.

Parameters:

n (int) – Number of samples to generate

Returns:

a (astropy.units.Quantity numpy.ndarray):

Semi-major axis in units of AU

e (float numpy.ndarray):

Eccentricity

p (float numpy.ndarray):

Geometric albedo

Rp (astropy.units.Quantity numpy.ndarray):

Planetary radius in units of earthRad

Return type:

tuple

EXOSIMS.PlanetPopulation.EarthTwinHabZone2 module

class EXOSIMS.PlanetPopulation.EarthTwinHabZone2.EarthTwinHabZone2(eta=0.1, erange=[0.0, 0.9], constrainOrbits=True, **specs)[source]

Bases: EarthTwinHabZone1

Population of Earth twins (1 R_Earth, 1 M_Eearth, 1 p_Earth) On eccentric habitable zone orbits (0.7 to 1.5 AU).

This implementation is intended to enforce this population regardless of JSON inputs. The only inputs that will not be disregarded are erange and constrainOrbits.

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and eccentricity are uniformly distributed with all other parameters constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

EXOSIMS.PlanetPopulation.EarthTwinHabZone3 module

class EXOSIMS.PlanetPopulation.EarthTwinHabZone3.EarthTwinHabZone3(eta=0.1, **specs)[source]

Bases: PlanetPopulation

Population of Earth twins (1 R_Earth, 1 M_Eearth, 1 p_Earth) On circular Habitable zone orbits (0.75 to 1.77 AU)

Note that these values may not be overwritten by user inputs. This implementation is intended to enforce this population regardless of JSON inputs.

dist_sma(a)[source]

Probability density function for uniform semi-major axis distribution in AU

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

f (float ndarray)

gen_sma(n)[source]

Generate semi-major axis values in AU

The Earth-twin population assumes a uniform distribution between the minimum and maximum values.

Parameters:

n (integer) – Number of samples to generate

Returns:

Semi-major axis values in units of AU

Return type:

a (astropy Quantity array)

EXOSIMS.PlanetPopulation.EarthTwinHabZoneSDET module

class EXOSIMS.PlanetPopulation.EarthTwinHabZoneSDET.EarthTwinHabZoneSDET(eta=0.1, **specs)[source]

Bases: PlanetPopulation

Population of Earth twins (1 R_Earth, 1 M_Eearth, 1 p_Earth) On circular Habitable zone orbits (0.7 to 1.5 AU)

Note that these values may not be overwritten by user inputs. This implementation is intended to enforce this population regardless of JSON inputs.

dist_radius(Rp)[source]

Probability density function for planetary radius in Earth radius

The prototype provides a log-uniform distribution between the minimum and maximum values.

Parameters:

Rp (float ndarray) – Planetary radius value(s) in Earth radius. Not an astropy quantity.

Returns:

Planetary radius probability density

Return type:

float ndarray

dist_sma(a)[source]

Probability density function for uniform semi-major axis distribution in AU

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

f (float ndarray)

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis is uniformly distributed and all other parameters are constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy.units.Quantity numpy.ndarray):

Semi-major axis in units of AU

e (float numpy.ndarray):

Eccentricity

p (float numpy.ndarray):

Geometric albedo

Rp (astropy.units.Quantity numpy.ndarray):

Planetary radius in units of earthRad

Return type:

tuple

EXOSIMS.PlanetPopulation.Guimond2019 module

class EXOSIMS.PlanetPopulation.Guimond2019.Guimond2019(eta=1, arange=[0.1, 50.0], erange=[0.0, 0.999], prange=[0.434, 0.434], constrainOrbits=False, **specs)[source]

Bases: PlanetPopulation

Population of Earth-Like Planets from Brown 2005 paper

This implementation is intended to enforce this population regardless of JSON inputs. The only inputs that will not be disregarded are erange and constrainOrbits.

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and eccentricity are uniformly distributed with all other parameters constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

gen_radius_nonorm(n)[source]

Generate planetary radius values in Earth radius. This one just generates a bunch of EarthRad

Parameters:

n (integer) – Number of target systems. Total number of samples generated will be, on average, n*self.eta

Returns:

Planet radius values in units of Earth radius

Return type:

astropy Quantity array

loguniform(low=0.001, high=1, size=None)[source]
Parameters:
  • float

  • float – upper end of samples to generate

  • integer – number of values to sample

Returns:

of logarithmically distributed values

Return type:

ndarray

EXOSIMS.PlanetPopulation.JupiterTwin module

class EXOSIMS.PlanetPopulation.JupiterTwin.JupiterTwin(eta=1, erange=[0.0, 0.048], constrainOrbits=True, **specs)[source]

Bases: PlanetPopulation

Population of Jupiter twins (11.209 R_Earth, 317.83 M_Eearth, 1 p_Earth) On eccentric orbits (0.7 to 1.5 AU)*5.204. Numbers pulled from nssdc.gsfc.nasa.gov/planetary/factsheet/jupiterfact.html

This implementation is intended to enforce this population regardless of JSON inputs. The only inputs that will not be disregarded are erange and constrainOrbits.

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and eccentricity are uniformly distributed with all other parameters constant.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

EXOSIMS.PlanetPopulation.KeplerLike1 module

class EXOSIMS.PlanetPopulation.KeplerLike1.KeplerLike1(smaknee=30, esigma=0.13962979814050144, prange=[0.083, 0.882], Rprange=[1, 22.6], **specs)[source]

Bases: PlanetPopulation

Population based on Kepler radius distribution with RV-like semi-major axis distribution with exponential decay.

Parameters:

**specs – user specified values

smaknee

Location (in AU) of semi-major axis decay point (knee). Not an astropy quantity.

Type:

float

esigma

Sigma value of Rayleigh distribution for eccentricity.

Type:

float

Notes: 1. The gen_mass function samples the Radius and calculates the mass from there. Any user-set mass limits are ignored. 2. The gen_albedo function samples the sma, and then calculates the albedos from there. Any user-set albedo limits are ignored. 3. The Rprange is fixed to (1,22.6) R_Earth and cannot be overwritten by user settings (the JSON input will be ignored) 4. The radius piece-wise distribution (from Fressin et al 2012) provides the normalization required to get the proper overall eta. The gen_radius method provided here normalizes in order to return exactly the number of samples requested. A second method (gen_radius_nonorm) is provided for generating the simulated universe population. The latter assumes a poisson distribution for occurences in each bin. 5. Eccentricity is assumed to be Rayleigh distributed with a user-settable sigma parameter (defaults to value from Fressin et al 2012).

dist_albedo(p)[source]

Probability density function for albedo

Parameters:

p (float ndarray) – Albedo value(s)

Returns:

Albedo probability density

Return type:

float ndarray

dist_eccen(e)[source]

Probability density function for eccentricity

Parameters:

e (float ndarray) – Eccentricity value(s)

Returns:

Eccentricity probability density

Return type:

float ndarray

dist_eccen_from_sma(e, a)[source]

Probability density function for eccentricity constrained by semi-major axis, such that orbital radius always falls within the provided sma range.

This provides a Rayleigh distribution between the minimum and maximum allowable values.

Parameters:
  • e (float ndarray) – Eccentricity values

  • a (float ndarray) – Semi-major axis value in AU. Not an astropy quantity.

Returns:

Probability density of eccentricity constrained by semi-major axis

Return type:

float ndarray

dist_radius(Rp)[source]

Probability density function for planetary radius in Earth radius

Parameters:

Rp (float ndarray) – Planetary radius value(s) in Earth radius. Not an astropy quantity.

Returns:

Planetary radius probability density

Return type:

float ndarray

dist_sma(a)[source]

Probability density function for semi-major axis in AU

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

float ndarray

gen_albedo(n)[source]

Generate geometric albedo values

The albedo is determined by sampling the semi-major axis distribution, and then calculating the albedo from the physical model.

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet albedo values

Return type:

float ndarray

gen_mass(n)[source]

Generate planetary mass values in Earth Mass

The mass is determined by sampling the radius and then calculating the mass from the physical model.

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet mass values in units of Earth mass

Return type:

astropy Quantity array

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis is distributed RV like with exponential decay. Eccentricity is a Rayleigh distribution. Albedo is dependent on the PlanetPhysicalModel but is calculated such that it is independent of other parameters. Planetary radius comes from the Kepler observations.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

gen_radius(n)[source]

Generate planetary radius values in Earth radius

Samples a radius distribution defined as log-uniform in each of 9 radius bins with fixed occurrence rates.

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet radius values in units of Earth radius

Return type:

astropy Quantity array

gen_radius_nonorm(n)[source]

Generate planetary radius values in Earth radius.

Samples a radius distribution defined as log-uniform in each of 9 radius bins with fixed occurrence rates. The rates in the bins determine the overall occurrence rates of all planets.

Parameters:

n (integer) – Number of target systems. Total number of samples generated will be, on average, n*self.eta

Returns:

Planet radius values in units of Earth radius

Return type:

astropy Quantity array

gen_sma(n)[source]

Generate semi-major axis values in AU

Samples a power law distribution with exponential turn-off determined by class attribute smaknee

Parameters:

n (integer) – Number of samples to generate

Returns:

Semi-major axis values in units of AU

Return type:

astropy Quantity array

EXOSIMS.PlanetPopulation.KeplerLike2 module

class EXOSIMS.PlanetPopulation.KeplerLike2.KeplerLike2(smaknee=30, esigma=0.25, **specs)[source]

Bases: KeplerLike1

Population based on Kepler radius distribution with RV-like semi-major axis distribution with exponential decay.

NOTE: This is an exact clone of KeplerLike1, but uses (approximate) inverse transform sampling instead of simple rejection sampling for performance improvements.

Parameters:

**specs – user specified values

smaknee

Location (in AU) of semi-major axis decay point (knee). Not an astropy quantity.

Type:

float

esigma

Sigma value of Rayleigh distribution for eccentricity.

Type:

float

Notes: 1. The gen_mass function samples the Radius and calculates the mass from there. Any user-set mass limits are ignored. 2. The gen_albedo function samples the sma, and then calculates the albedos from there. Any user-set albedo limits are ignored. 3. The Rprange is fixed to (1,22.6) R_Earth and cannot be overwritten by user settings (the JSON input will be ignored) 4. The radius piece-wise distribution provides the normalization required to get the proper overall eta. The gen_radius method provided here normalizes in order to return exactly the number of samples requested. A second method (gen_radius_nonorm) is provided for generating the simulated universe population. The latter assumes a poisson distribution for occurences in each bin. 5. Eccentricity is assumed to be Rayleigh distributed with a user-settable sigma parameter (defaults to 0.25).

gen_sma(n)[source]

Generate semi-major axis values in AU

Samples a power law distribution with exponential turn-off determined by class attribute smaknee

Parameters:

n (integer) – Number of samples to generate

Returns:

Semi-major axis in units of AU

Return type:

a (astropy Quantity array)

EXOSIMS.PlanetPopulation.KnownRVPlanets module

class EXOSIMS.PlanetPopulation.KnownRVPlanets.KnownRVPlanets(smaknee=30, esigma=0.25, rvplanetfilepath=None, planetfile='planets_2019.05.31_11.18.02.votable', **specs)[source]

Bases: KeplerLike1

Population consisting only of known RV planets. Eccentricity and sma distributions are taken from KeplerLike1 (Rayleigh and power law with exponential decay, respectively). Mass is sampled from power law and radius is assumed to be calculated from mass via the physical model.

The data file read in by this class also provides all of the information about the target stars, and so no StarCatalog object is needed (only the KnownRvPlanetsTargetList implementation).

To download a new copy of the data file:

  1. Navigate to the IPAC exoplanet archive at http://exoplanetarchive.ipac.caltech.edu/cgi-bin/TblView/nph-tblView?app=ExoTbls&config=planets

  2. Type ‘radial’ (minus quotes) in the ‘Discovery Method’ search box and hit enter.

  3. In the ‘Download Table’ menu select ‘VOTable Format’, ‘Download all Columns’ and ‘Download Currently Filtered Rows’.

  4. In the ‘Download Table’ menu click ‘Download Table’.

Parameters:

**specs – user specified values

smaknee

Location (in AU) of semi-major axis decay point (knee). Not an astropy quantity.

Type:

float

esigma

Sigma value of Rayleigh distribution for eccentricity.

Type:

float

rvplanetfilepath

Full path to RV planet votable file from IPAC. If None, assumes default file in PlanetPopulation directory of EXOSIMS.

Type:

string

period

Orbital period in units of day. Error in perioderr.

Type:

astropy Quantity array

tper

Periastron time in units of jd. Error in tpererr.

Type:

astropy Time

Notes:

gen_mass(n)[source]

Generate planetary mass values in Earth mass

The mass is determined by sampling the RV mass distribution from Cumming et al. 2010

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet mass values in units of Earth mass

Return type:

Mp (astropy Quantity array)

gen_radius(n)[source]

Generate planetary radius values in Earth radius

Samples the mass distribution and then converts to radius using the physical model.

Parameters:

n (integer) – Number of samples to generate

Returns:

Planet radius values in units of Earth radius

Return type:

Rp (astropy Quantity array)

EXOSIMS.PlanetPopulation.SAG13 module

class EXOSIMS.PlanetPopulation.SAG13.SAG13(SAG13coeffs=[[0.38, -0.19, 0.26, 0.0], [0.73, -1.18, 0.59, 3.4]], SAG13starMass=1.0, Rprange=[0.6666666666666666, 17.0859375], arange=[0.09084645, 1.45354324], **specs)[source]

Bases: KeplerLike2

Planet Population module based on SAG13 occurrence rates.

This is the current working model based on averaging multiple studies. These do not yet represent official scientific values.

SAG13coeffs

Coefficients used by the SAG13 broken power law. The 4 lines correspond to Gamma, alpha, beta, and the minimum radius.

Type:

float 4x2 ndarray

Gamma

Gamma coefficients used by SAG13 broken power law.

Type:

float ndarray

alpha

Alpha coefficients used by SAG13 broken power law.

Type:

float ndarray

beta

Beta coefficients used by SAG13 broken power law.

Type:

float ndarray

Rplim

Minimum radius used by SAG13 broken power law.

Type:

float ndarray

SAG13starMass

Assumed stellar mass corresponding to the given set of coefficients.

Type:

astropy Quantity

mu

Gravitational parameter associated with SAG13starMass.

Type:

astropy Quantity

Ca

Constants used for sampling.

Type:

float 2x1 ndarray

dist_radius(Rp)[source]

Marginalized probability density function for planetary radius in Earth radius.

Parameters:

Rp (float ndarray) – Planetary radius value(s) in Earth radius. Not an astropy quantity.

Returns:

Planetary radius probability density

Return type:

float ndarray

dist_sma(a)[source]

Marginalized probability density function for semi-major axis in AU.

Parameters:

a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity.

Returns:

Semi-major axis probability density

Return type:

float ndarray

dist_sma_given_radius(a, beta, m, C, smaknee)[source]

Conditional probability density function of semi-major axis given planetary radius.

Parameters:
  • a (float ndarray) – Semi-major axis value(s) in AU. Not an astropy quantity

  • beta (float) – Exponent for distribution

  • m (float) – Gravitational parameter (AU3/year2)

  • C (float) – Normalization for distribution

  • smaknee (float) – Coefficient for decay

Returns:

Probability density

Return type:

float ndarray

dist_sma_radius(a, R)[source]

Joint probability density function for semi-major axis (AU) and planetary radius in Earth radius.

This method performs a change of variables on the SAG13 broken power law (originally in planetary radius and period).

Parameters:
  • a (float ndarray) – Semi-major axis values in AU. Not an astropy quantity

  • R (float ndarray) – Planetary radius values in Earth radius. Not an astropy quantity

Returns:

Joint (semi-major axis and planetary radius) probability density matrix of shape (len(R),len(a))

Return type:

float ndarray

gen_plan_params(n)[source]

Generate semi-major axis (AU), eccentricity, geometric albedo, and planetary radius (earthRad)

Semi-major axis and planetary radius are jointly distributed. Eccentricity is a Rayleigh distribution. Albedo is dependent on the PlanetPhysicalModel but is calculated such that it is independent of other parameters.

Parameters:

n (integer) – Number of samples to generate

Returns:

a (astropy Quantity array):

Semi-major axis in units of AU

e (float ndarray):

Eccentricity

p (float ndarray):

Geometric albedo

Rp (astropy Quantity array):

Planetary radius in units of earthRad

Return type:

tuple

gen_radius_sma(n)[source]

Generate radius values in earth radius and semi-major axis values in AU.

This method is called by gen_radius and gen_sma.

Parameters:

n (integer) – Number of samples to generate

Returns:

Rp (astropy Quantity array):

Planet radius values in units of Earth radius

a (astropy Quantity array):

Semi-major axis values in units of AU

Return type:

tuple

EXOSIMS.PlanetPopulation.SolarSystem module

class EXOSIMS.PlanetPopulation.SolarSystem.SolarSystem(prange=[0.1, 0.7], Rprange=[0.01, 30.0], **specs)[source]

Bases: PlanetPopulation

Population of Earth-Like Planets from Brown 2005 paper

This implementation is intended to enforce this population regardless of JSON inputs. The only inputs that will not be disregarded are erange and constrainOrbits.

gen_plan_params(nPlans)[source]

Values taken From mallama2018PlantProperties Seidlemenn 1992 :param float: nPlans, the number of planets