primpy package

Calculations for the primordial Universe.

primpy subpackages

primpy.bigbang module

General setup for equations for standard Big Bang cosmology.

primpy.bigbang.get_H0(h, units='planck')[source]

Get present-day Hubble parameter from little Hubble h.

primpy.bigbang.get_a0(h, Omega_K0, units='planck')[source]

Get present-day scale factor from curvature density parameter.

primpy.bigbang.get_N_BBN(h, Omega_K0)[source]

Get the epoch of Big Bang nucleosynthesis in terms of e-folds (in Planck units).

primpy.bigbang.get_w_delta_reh(N_end, N_reh, log_cHH_end, log_cHH_reh)[source]

Get the e.o.s. parameter for reheating w_reh from e-folds and comoving Hubble horizon.

primpy.bigbang.get_rho_crit_kg_im3(h)[source]

Get present-day critical density from little Hubble h.

primpy.bigbang.get_Omega_r0(h)[source]

Get present-day radiation density parameter from little Hubble h.

primpy.bigbang.Hubble_parameter(N, Omega_m0, Omega_K0, h, units='planck')[source]

Hubble parameter (in reduced Planck units) at N=ln(a) during standard Big Bang.

Parameters:
Nfloat, np.ndarray

e-folds of the scale factor N=ln(a) during standard Big Bang evolution, where the scale factor would be given in reduced Planck units (same as output from primpy).

Omega_m0float

matter density parameter today

Omega_K0float

curvature density parameter today

hfloat

dimensionless Hubble parameter today, “little h”

unitsstr

Output units, can be any of {‘planck’, ‘H0’, ‘SI’} returning units of 1/tp, km/s/Mpc or 1/s respectively.

Returns:
Hfloat

Hubble parameter during standard Big Bang evolution of the Universe. In reduced Planck units [tp^-1].

Notes

Omega_r0 is derived from the Hubble parameter using Planck’s law. Omega_L0 is derived from the other density parameters to sum to one.

primpy.bigbang.no_Big_Bang_line(Omega_m0)[source]

Return Omega_L0 for dividing line between universes with/without Big Bang.

Parameters:
Omega_m0float

matter density parameter today

Returns:
Omega_L0float

Density parameter of cosmological constant Lambda along the dividing line between a Big Bang evolution (for smaller Omega_L0) and universes without a Big Bang (for larger Omega_L0).

primpy.bigbang.expand_recollapse_line(Omega_m0)[source]

Return Omega_L0 for dividing line between expanding/recollapsing universes.

Parameters:
Omega_m0float

matter density parameter today

Returns:
Omega_L0float

Density parameter of cosmological constant Lambda along the dividing line between expanding (for larger Omega_L0) and recollapsing (for smaller Omega_L0) universes.

primpy.bigbang.comoving_Hubble_horizon(N, Omega_m0, Omega_K0, h, units='planck')[source]

Comoving Hubble horizon at N=ln(a) during standard Big Bang.

Parameters:
Nfloat, np.ndarray

e-folds of the scale factor N=ln(a) during standard Big Bang evolution, where the scale factor would be given in reduced Planck units (same as output from primpy).

Omega_m0float

matter density parameter today

Omega_K0float

curvature density parameter today

hfloat

dimensionless Hubble parameter today, “little h”

unitsstr

Output units, can be any of {‘planck’, ‘Mpc’, ‘SI’} returning units of lp, Mpc, or m respectively.

Returns:
cHHfloat

Comoving Hubble horizon during standard Big Bang evolution of the Universe.

Notes

Omega_r0 is derived from the Hubble parameter using Planck’s law. Omega_L0 is derived from the other density parameters to sum to one.

primpy.bigbang.conformal_time(N_start, N, Omega_m0, Omega_K0, h)[source]

Conformal time during standard Big Bang evolution from N_start to N.

Parameters:
N_startfloat

e-folds of the scale factor N=ln(a) during standard Big Bang evolution at lower integration limit (e.g. at end of inflation), where the scale factor would be given in reduced Planck units (same as output from primpy).

Nfloat, np.ndarray

e-folds of the scale factor N=ln(a) during standard Big Bang evolution at upper integration limit (e.g. at end of inflation), where the scale factor would be given in reduced Planck units (same as output from primpy).

Omega_m0float

matter density parameter today

Omega_K0float

curvature density parameter today

hfloat

dimensionless Hubble parameter today, “little h”

Returns:
etafloat, np.ndarray

conformal time passing between a_start and a during standard Big Bang evolution of the Universe. Same shape as N.

Notes

Omega_r0 is derived from the Hubble parameter using Planck’s law. Omega_L0 is derived from the other density parameters to sum to one.

primpy.bigbang.conformal_time_ratio(Omega_m0, Omega_K0, h, b_forward, b_backward=None)[source]

Conformal time ratio before to after the end of inflation (until today).

Parameters:
Omega_m0float

matter density parameter today

Omega_K0float

curvature density parameter today

hfloat

dimensionless Hubble parameter today, “little h”

b_forwardBunch object same as returned by scipy.integrate.solve_ivp()

Solution returned by primpy.solver.solve(). Needs to have been run with track_eta=True.

b_backwardBunch object same as returned by scipy.integrate.solve_ivp(), default: None

Additional solution returned by primpy.solver.solve(). This second solution is assumed to be an integration from inflation start backwards in time.

Returns:
ratiofloat

Ratio of conformal time before (during and before inflation) to after (from the end of inflation until today). Needs to be >1 in order to solve the horizon problem.

primpy.equations module

General setup for ODEs.

class primpy.equations.Equations[source]

Bases: ABC

Base class for equations.

Allows one to compute derivatives and derived variables. Most of the other classes take ‘equations’ as an object.

Attributes:
idxdict

dictionary mapping variable names to indices in the solution vector y

independent_variablestring

name of independent variable

Methods

__call__(x, y)

Vector of derivatives.

add_variable(*args)

Add dependent variables to the equations.

sol(sol, **kwargs)

Post-processing of scipy.integrate.solve_ivp() solution.

__call__(x, y)[source]

Vector of derivatives.

Parameters:
xfloat

independent variable

ynp.ndarray

dependent variables

Returns:
dynp.ndarray

Vector of derivatives

sol(sol, **kwargs)[source]

Post-processing of scipy.integrate.solve_ivp() solution.

add_variable(*args)[source]

Add dependent variables to the equations.

  • creates an index for the location of variable in y

  • creates a class method of the same name with signature name(self, x, y) that should be used to extract the variable value in an index-independent manner.

Parameters:
*argsstr

Name of the dependent variables

primpy.events module

Setup for event tracking in scipy.integrate.solve_ivp().

class primpy.events.Event(equations, direction=0, terminal=False, value=0)[source]

Bases: object

Base class for event tracking.

Gives a more usable wrapper to callable event to be passed to scipy.integrate.solve_ivp().

Parameters:
equations: Equations

The equations for computing derived variables.

direction: int, default: 0

The direction of the root finding (if any), one of {-1, 0, +1}.

terminal: bool, default: False

Whether to stop at this root or continue integrating.

value: float, default: 0

Offset to root.

Methods

__call__(x, y)

Vector of derivatives.

__call__(x, y)[source]

Vector of derivatives.

Parameters:
xfloat

independent variable

ynp.ndarray

dependent variables

Returns:
rootfloat

event occurs when this is zero from a given direction

class primpy.events.UntilTEvent(equations, value, direction=0, terminal=True)[source]

Bases: Event

Stop after a given amount of time t has passed.

Methods

__call__(x, y)

Root of t - value.

__call__(x, y)[source]

Root of t - value.

class primpy.events.UntilNEvent(equations, value, direction=0, terminal=True)[source]

Bases: Event

Stop after a given number of e-folds _N.

Methods

__call__(x, y)

Root of _N - value.

__call__(x, y)[source]

Root of _N - value.

class primpy.events.InflationEvent(equations, direction=0, terminal=False, value=0, **kwargs)[source]

Bases: Event

Track inflation start/end.

Methods

__call__(x, y)

Root of V - dphidt**2.

__call__(x, y)[source]

Root of V - dphidt**2.

class primpy.events.AfterInflationEndEvent(equations, direction=1, terminal=True, value=0)[source]

Bases: Event

Go a bit past the end of inflation.

Methods

__call__(x, y)

Root of w - value.

__call__(x, y)[source]

Root of w - value.

class primpy.events.CollapseEvent(equations, direction=0, terminal=True, value=0)[source]

Bases: Event

Stop if Universe collapses, i.e. test whether H**2 turns negative.

Methods

__call__(x, y)

Root of H2 - value.

__call__(x, y)[source]

Root of H2 - value.

class primpy.events.Phi0Event(equations, direction=0, terminal=True, value=0)[source]

Bases: Event

Track zero crossings of inflaton phi.

Methods

__call__(x, y)

Root of phi - value.

__call__(x, y)[source]

Root of phi - value.

class primpy.events.ModeExitEvent(equations, value, direction=0, terminal=True)[source]

Bases: Event

Track when mode exits the horizon aH.

Methods

__call__(x, y)

Root of _logaH - log(value).

__call__(x, y)[source]

Root of _logaH - log(value).

primpy.exceptionhandling module

Custom exceptions and warnings.

exception primpy.exceptionhandling.PrimpyError[source]

Bases: Exception

Base class for exceptions in primpy.

exception primpy.exceptionhandling.InflationStartError(message, *args, **kwargs)[source]

Bases: PrimpyError

Exception when the inflation start condition for open or closed universes is violated.

Parameters:
messagestr

Explanation of the error.

Other Parameters:
geometrystr, default: “all types of”

Should be either ‘open’ or ‘closed’ to relate to the respective condition at inflation start.

exception primpy.exceptionhandling.StepSizeError(message, *args)[source]

Bases: PrimpyError

Warning when the scipy integrator failed because of a too small step size.

Parameters:
messagestr

Explanation of the warning.

exception primpy.exceptionhandling.InsufficientInflationError(message, *args)[source]

Bases: PrimpyError

Exception when there are insufficient number of e-folds for inflation.

Parameters:
messagestr

Explanation of the error.

exception primpy.exceptionhandling.BigBangError(message, *args)[source]

Bases: PrimpyError

Exceptions for the standard Big Bang evolution.

Parameters:
messagestr

Explanation of the error.

exception primpy.exceptionhandling.PrimpyWarning[source]

Bases: UserWarning

Base class for warnings in primpy.

exception primpy.exceptionhandling.InflationWarning(message, *args)[source]

Bases: PrimpyWarning

Warnings for the inflationary background evolution.

Parameters:
messagestr

Explanation of the warning.

exception primpy.exceptionhandling.CollapseWarning(message, *args)[source]

Bases: InflationWarning

Warning when the Universe has collapsed.

Parameters:
messagestr

Explanation of the warning.

exception primpy.exceptionhandling.InflationStartWarning(message, *args, **kwargs)[source]

Bases: InflationWarning

Warnings when the start of inflation could not be determined.

Parameters:
messagestr

Explanation of the warning.

Other Parameters:
eventsdict, default: None

Dictionary of events captured. Can be any of N_events, t_events, phi_events.

exception primpy.exceptionhandling.InflationEndWarning(message, *args, **kwargs)[source]

Bases: InflationWarning

Warnings when the end of inflation could not be determined.

Parameters:
messagestr

Explanation of the warning.

Other Parameters:
eventsdict, default: None

Dictionary of events captured. Can be any of N_events, t_events, phi_events.

solBunch object same as returned by scipy.integrate.solve_ivp(), default: None

Bunch object returned by primpy.solver.solve().

exception primpy.exceptionhandling.BigBangWarning(message, *args)[source]

Bases: PrimpyWarning

Warnings for the standard Big Bang evolution.

Parameters:
messagestr

Explanation of the warning.

primpy.inflation module

General setup for equations for cosmic inflation.

class primpy.inflation.InflationEquations(K, potential, verbose=False)[source]

Bases: Equations, ABC

Base class for inflation equations.

Methods

H(x, y)

Hubble parameter.

H2(x, y)

Hubble parameter squared.

V(x, y)

Inflationary Potential.

d2Vdphi2(x, y)

Second derivative of inflationary potential.

dVdphi(x, y)

First derivative of inflationary potential.

get_H2(N, dphi, V, K)

Get the Hubble parameter squared from the background equations.

get_d2H(N, H, dH, dphi, d2phi, K)

Get the 2nd time derivative of the Hubble parameter from the background equations.

get_d2phi(H2, dH_H, dphi, dV)

Get the 2nd time derivative of the inflaton field from the background equations.

get_d3H(N, H, dH, d2H, dphi, d2phi, d3phi, K)

Get the 3rd time derivative of the Hubble parameter from the background equations.

get_d3phi(H, dH, d2H, dphi, d2phi, dV, d2V)

Get the 3rd time derivative of the inflaton field from the background equations.

get_d4phi(H, dH, d2H, d3H, dphi, d2phi, ...)

Get the 4th time derivative of the inflaton field from the background equations.

get_dH(N, H, dphi, K)

Get the 1st time derivative of the Hubble parameter from the background equations.

get_dH_H(H2, dphi, K)

Get the 1st time derivative of the Hubble parameter normalised by the Hubble parameter..

get_delta_1(H, dH, dphi, d2phi)

Get the 2nd slow-roll parameter for n=1.

get_delta_2(H, dH, d2H, dphi, d2phi, d3phi)

Get the 2nd slow-roll parameter for n=2.

get_delta_3(H, dH, d2H, d3H, dphi, d2phi, ...)

Get the 2nd slow-roll parameter for n=3.

get_epsilon_1H(H, dH)

Get the 1st Hubble flow parameter.

get_epsilon_2H(H, dH, d2H[, kind])

Get the 2nd Hubble flow parameter.

get_epsilon_3H(H, dH, d2H, d3H[, kind])

Get the 3rd Hubble flow parameter.

inflating(x, y)

Inflation diagnostic for event tracking.

postprocessing_inflation_end(sol)

Extract end point of inflation from event tracking.

postprocessing_inflation_start(sol)

Extract starting point of inflation from event tracking.

sol(sol, **kwargs)

Post-processing of scipy.integrate.solve_ivp() solution.

w(x, y)

Equation of state parameter.

static get_H2(N, dphi, V, K)[source]

Get the Hubble parameter squared from the background equations.

Parameters:
Nfloat or array_like

Number of e-folds N=ln(a).

dphifloat or array_like

1st derivative of inflaton field.

Vfloat or array_like

Inflation potential at phi.

Kint

Curvature parameter.

Returns:
H2float or array_like

Hubble parameter squared.

static get_dH(N, H, dphi, K)[source]

Get the 1st time derivative of the Hubble parameter from the background equations.

Parameters:
Nfloat or array_like

Number of e-folds N=ln(a).

Hfloat or array_like

Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

Kint

Curvature parameter.

Returns:
dHfloat or array_like

1st derivative of Hubble parameter.

get_dH_H(H2, dphi, K)[source]

Get the 1st time derivative of the Hubble parameter normalised by the Hubble parameter..

Parameters:
Nfloat or array_like

Number of e-folds N=ln(a).

H2float or array_like

Hubble parameter squared.

dphifloat or array_like

1st derivative of inflaton field.

Kint

Curvature parameter.

Returns:
dH_Hfloat or array_like

1st derivative of Hubble parameter normalised by Hubble parameter: dH/H.

static get_d2H(N, H, dH, dphi, d2phi, K)[source]

Get the 2nd time derivative of the Hubble parameter from the background equations.

Parameters:
Nfloat or array_like

Number of e-folds N=ln(a).

Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

Kint

Curvature parameter.

Returns:
d2Hfloat or array_like

2nd derivative of Hubble parameter.

static get_d3H(N, H, dH, d2H, dphi, d2phi, d3phi, K)[source]

Get the 3rd time derivative of the Hubble parameter from the background equations.

Parameters:
Nfloat or array_like

Number of e-folds N=ln(a).

Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

d3phifloat or array_like

3rd derivative of inflaton field.

Kint

Curvature parameter.

Returns:
d3Hfloat or array_like

3rd derivative of Hubble parameter.

static get_d2phi(H2, dH_H, dphi, dV)[source]

Get the 2nd time derivative of the inflaton field from the background equations.

Parameters:
H2float or array_like

Hubble parameter squared.

dH_Hfloat or array_like

1st derivative of Hubble parameter normalised by Hubble parameter dH/H.

dphifloat or array_like

1st derivative of inflaton field.

dVfloat or array_like

1st derivative of inflation potential at phi.

Returns:
d2phifloat or array_like

2nd derivative of inflaton field.

static get_d3phi(H, dH, d2H, dphi, d2phi, dV, d2V)[source]

Get the 3rd time derivative of the inflaton field from the background equations.

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

dVfloat or array_like

1st derivative of inflation potential at phi.

d2Vfloat or array_like

2nd derivative of inflation potential at phi.

Returns:
d3phifloat or array_like

3rd derivative of inflaton field.

static get_d4phi(H, dH, d2H, d3H, dphi, d2phi, d3phi, dV, d2V, d3V)[source]

Get the 4th time derivative of the inflaton field from the background equations.

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

d3Hfloat or array_like

3rd derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

d3phifloat or array_like

3rd derivative of inflaton field.

dVfloat or array_like

1st derivative of inflation potential at phi.

d2Vfloat or array_like

2nd derivative of inflation potential at phi.

d3Vfloat or array_like

3rd derivative of inflation potential at phi.

Returns:
d4phifloat or array_like

4th derivative of inflaton field.

static get_epsilon_1H(H, dH)[source]

Get the 1st Hubble flow parameter.

This definition of the 1st Hubble flow parameter was suggested e.g. by Stewart & Lyth (1993) in eq. (23). https://arxiv.org/abs/gr-qc/9302019

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

Returns:
epsilon_1Hfloat or array_like

1st Hubble flow parameter.

static get_epsilon_2H(H, dH, d2H, kind=None)[source]

Get the 2nd Hubble flow parameter.

This definition of the 2nd Hubble flow parameter was suggested e.g. by Leach, Liddle, Martin & Schwarz (2003) in eq. (15). https://arxiv.org/abs/astro-ph/0202094

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

kindstr or None, default=None

Switch to alternative formulation given by Gong (2004).

Returns:
epsilon_2Hfloat or array_like

2nd Hubble flow parameter.

static get_epsilon_3H(H, dH, d2H, d3H, kind=None)[source]

Get the 3rd Hubble flow parameter.

This definition of the 3rd Hubble flow parameter was suggested e.g. by Leach, Liddle, Martin & Schwarz (2003) in eq. (15). https://arxiv.org/abs/astro-ph/0202094

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

d3Hfloat or array_like

3rd derivative of Hubble parameter.

kindstr or None, default=None

Switch to alternative formulation given by Gong (2004).

Returns:
epsilon_3Hfloat or array_like

3rd Hubble flow parameter.

static get_delta_1(H, dH, dphi, d2phi)[source]

Get the 2nd slow-roll parameter for n=1.

This definition of the 2nd slow roll parameter was suggested e.g. by Stewart & Lyth (1993) in eq. (23). https://arxiv.org/abs/gr-qc/9302019

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

Returns:
delta_1float or array_like
static get_delta_2(H, dH, d2H, dphi, d2phi, d3phi)[source]

Get the 2nd slow-roll parameter for n=2.

This definition of higher orders of the 2nd slow roll parameter was suggested e.g. by Stewart & Gong (2001) in eq. (22). https://arxiv.org/abs/astro-ph/0101225

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

d3phifloat or array_like

3rd derivative of inflaton field.

Returns:
delta_2float or array_like
static get_delta_3(H, dH, d2H, d3H, dphi, d2phi, d3phi, d4phi)[source]

Get the 2nd slow-roll parameter for n=3.

This definition of higher orders of the 2nd slow roll parameter was suggested e.g. by Stewart & Gong (2001) in eq. (22). https://arxiv.org/abs/astro-ph/0101225

Parameters:
Hfloat or array_like

Hubble parameter.

dHfloat or array_like

1st derivative of Hubble parameter.

d2Hfloat or array_like

2nd derivative of Hubble parameter.

d3Hfloat or array_like

3rd derivative of Hubble parameter.

dphifloat or array_like

1st derivative of inflaton field.

d2phifloat or array_like

2nd derivative of inflaton field.

d3phifloat or array_like

3rd derivative of inflaton field.

d4phifloat or array_like

4th derivative of inflaton field.

Returns:
delta_3float or array_like
H(x, y)[source]

Hubble parameter.

H2(x, y)[source]

Hubble parameter squared.

V(x, y)[source]

Inflationary Potential.

dVdphi(x, y)[source]

First derivative of inflationary potential.

d2Vdphi2(x, y)[source]

Second derivative of inflationary potential.

w(x, y)[source]

Equation of state parameter.

inflating(x, y)[source]

Inflation diagnostic for event tracking.

postprocessing_inflation_start(sol)[source]

Extract starting point of inflation from event tracking.

postprocessing_inflation_end(sol)[source]

Extract end point of inflation from event tracking.

sol(sol, **kwargs)[source]

Post-processing of scipy.integrate.solve_ivp() solution.

primpy.initialconditions module

Initial conditions for inflation.

class primpy.initialconditions.InitialConditions(equations, **kwargs)[source]

Bases: ABC

Base class for initial conditions.

Methods

__call__(y0, **ivp_kwargs)

Initialise background equations of inflation.

abstract __call__(y0, **ivp_kwargs)[source]

Initialise background equations of inflation.

class primpy.initialconditions.SlowRollIC(equations, phi_i, t_i=None, eta_i=None, x_end=1e+300, **kwargs)[source]

Bases: InitialConditions

Slow-roll initial conditions given phi_i and either of N_i or Omega_Ki.

Class for setting up initial conditions during slow-roll inflation where the potential energy dominates over the kinetic energy.

Methods

__call__(y0, **ivp_kwargs)

Set background equations of inflation for _N, phi and dphi.

__call__(y0, **ivp_kwargs)[source]

Set background equations of inflation for _N, phi and dphi.

class primpy.initialconditions.InflationStartIC(equations, phi_i, t_i=None, eta_i=None, x_end=1e+300, **kwargs)[source]

Bases: InitialConditions

Inflation start initial conditions given phi_i and either of N_i or Omega_Ki.

Class for setting up initial conditions at the start of inflation, when the curvature density parameter was maximal after kinetic dominance.

Methods

__call__(y0, **ivp_kwargs)

Set background equations of inflation for _N, phi and dphi.

__call__(y0, **ivp_kwargs)[source]

Set background equations of inflation for _N, phi and dphi.

class primpy.initialconditions.ISIC_Nt(equations, N_tot, phi_i_bracket, t_i=None, eta_i=None, x_end=1e+300, verbose=False, **kwargs)[source]

Bases: InflationStartIC

Inflation start initial conditions given N_tot and either of N_i or Omega_Ki.

Methods

__call__(y0, **ivp_kwargs)

Set background equations of inflation optimizing for N_tot.

__call__(y0, **ivp_kwargs)[source]

Set background equations of inflation optimizing for N_tot.

class primpy.initialconditions.ISIC_NsOk(equations, N_star, Omega_K0, h, phi_i_bracket, t_i=None, eta_i=None, x_end=1e+300, verbose=False, **kwargs)[source]

Bases: InflationStartIC

Inflation start initial conditions given N_star, Omega_K0, h.

Additionally either N_i or Omega_Ki need to be specified.

Methods

__call__(y0, **ivp_kwargs)

Set background equations of inflation optimizing for N_star.

__call__(y0, **ivp_kwargs)[source]

Set background equations of inflation optimizing for N_star.

primpy.oscode_solver module

Setup for running pyoscode.solve().

primpy.oscode_solver.solve_oscode(background, k, vacuum=None, drop_closed_large_scales=True, fac_beg=0, fac_end=100, rtol=5e-05, even_grid=False, num_eval=0, **kwargs)[source]

Run pyoscode.solve() and store information for post-processing.

This is a wrapper around pyoscode.solve() to calculate the solution to the Mukhanov-Sasaki equation.

Parameters:
backgroundBunch object same as returned by scipy.integrate.solve_ivp()

Bunch object as returned by primpy.solver.solve(). Solution to the inflationary background equations used to calculate the frequency and damping term passed to oscode.

kint, float, np.ndarray

Comoving wavenumber used to evolve the Mukhanov-Sasaki equation.

vacuumtuple

Set of vacuum initial conditions to be computed. Choose any of (‘k’, ‘HD’, ‘RST’). default : (‘RST’, )

drop_closed_large_scalesbool

If true, this will set the PPS for closed universes on comoving scales of k < 1 to close to zero (1e-30). Strictly speaking, the PPS for closed universes is only defined for rational numbers k > 2. default : True

fac_begint, float

Integration of the mode evolution starts when the considered scale 1/k is within a factor of fac_beg of the comoving Hubble horizon, i.e. when 1/k > 1/aH / fac_beg. fac_beg=0 starts integration immediately. default : 0

fac_endint, float

Integration of the mode evolution stops when the considered scale 1/k exceeds the comoving Hubble horizon by a factor of fac_end, i.e. when 1/k > 1/aH * fac_end. default : 100

rtolfloat

Tolerance passed to pyoscode. default : 5e-5

even_gridbool

Set this to True if the grid of the independent variable is equally spaced. default : False

num_evalint

Number of interpolation points used for dense output. This number is applied dynamically to the range determined from the parameters fac_beg and fac_end. Zero means no dense output. default : 0

Returns:
solBunch object same as returned by scipy.integrate.solve_ivp()

Solution to the inverse value problem, containing the primordial power spectrum value corresponding to the wavenumber k. Monkey-patched version of the Bunch type usually returned by scipy.integrate.solve_ivp().

Other Parameters:
y0(float, float, float, float)

Initial values (y0_1, dy0_1, y0_2, dy0_2) of perturbations and their derivatives for two independent solutions. The perturbations (y0_1, y0_2) are scaled with k and their derivatives with k**2 in order to produce freeze-out values of about order(~1). default : determined by input inflationary potential

primpy.parameters module

Constants and parameters for primpy.

primpy.parameters.K_STAR = 0.05

wavenumber at pivot scale in units of 1/Mpc

primpy.parameters.T_CMB_K = 2.72548

+/- 0.00057, in Kelvin, arXiv:0911.1955

primpy.parameters.T_nu_TCMB = 0.7137658555036082

T_ncdm in CLASS

primpy.parameters.N_eff = 3.044

equivalent to N_ur in CLASS for massless neutrinos

primpy.parameters.z_BBN = 1000000000.0

rough estimate of redshift of Big Bang Nucleosynthesis

primpy.parameters.g0 = 3.937090909090909

effective number of relativistic degrees of freedom

primpy.parameters.rho_gamma0_kg_im3 = 4.644956134403369e-31

energy density of photons in kg/m^3

primpy.parameters.rho_nu0_kg_im3 = 3.211126340248152e-31

energy density of relativistic neutrinos in kg/m^3

primpy.parameters.rho_r0_kg_im3 = 7.85608247465152e-31

total radiation energy density in kg/m^3

primpy.perturbations module

Comoving curvature perturbations.

class primpy.perturbations.PrimordialPowerSpectrum(background, k, vacuum=None)[source]

Bases: object

Primordial Power spectrum of curvature perturbations.

class primpy.perturbations.Perturbation(background, k)[source]

Bases: ABC

Perturbation for wavenumber k.

Methods

oscode_postprocessing(oscode_sol[, vacuum])

Post-processing for pyoscode.solve() solution.

oscode_postprocessing(oscode_sol, vacuum=None)[source]

Post-processing for pyoscode.solve() solution.

Translate oscode dictionary output to solve_ivp output with attributes x and y.

Parameters:
oscode_sollist

List [scalar_1, scalar_2, tensor_1, tensor_2] of two independent solutions each for both scalar and tensor modes, where each element is a dictionary returned by pyoscode.solve().

class primpy.perturbations.Mode(background, k, **kwargs)[source]

Bases: Equations, ABC

Template for scalar or tensor modes.

Methods

get_vacuum_ic_HD()

Get initial conditions for HD vacuum.

get_vacuum_ic_RST()

Get initial conditions for RST vacuum.

get_vacuum_ic_k()

Get initial conditions for HD approximation.

mukhanov_sasaki_frequency_damping()

Frequency and damping term of the Mukhanov-Sasaki equations.

abstract mukhanov_sasaki_frequency_damping()[source]

Frequency and damping term of the Mukhanov-Sasaki equations.

abstract get_vacuum_ic_k()[source]

Get initial conditions for HD approximation.

abstract get_vacuum_ic_HD()[source]

Get initial conditions for HD vacuum.

abstract get_vacuum_ic_RST()[source]

Get initial conditions for RST vacuum.

class primpy.perturbations.ScalarMode(background, k, vacuum=None, **kwargs)[source]

Bases: Mode, ABC

Template for scalar modes.

class primpy.perturbations.TensorMode(background, k, vacuum=None, **kwargs)[source]

Bases: Mode, ABC

Template for tensor modes.

primpy.potentials module

Inflationary potentials.

class primpy.potentials.InflationaryPotential(**pot_kwargs)[source]

Bases: ABC

Base class for inflaton potential and derivatives.

Attributes:
name

Name of the inflationary potential.

perturbation_ic

Set of well scaling initial conditions for perturbation module.

phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

tag

3 letter tag identifying the type of inflationary potential.

tex

Tex string useful for labelling the inflationary potential.

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1(phi)

Approximation of 1st Hubble flow parameter with potential slow-roll parameters.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2(phi)

Approximation of 2nd Hubble flow parameter with potential slow-roll parameters.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3(phi)

Approximation of 3rd Hubble flow parameter with potential slow-roll parameters.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4(phi)

Approximation of 4th Hubble flow parameter with potential slow-roll parameters.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

sr_As2Lambda(A_s, phi_star, N_star, **pot_kwargs)

Get potential amplitude Lambda from PPS amplitude A_s using slow-roll approximation.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

abstract property tag

3 letter tag identifying the type of inflationary potential.

abstract property name

Name of the inflationary potential.

abstract property tex

Tex string useful for labelling the inflationary potential.

abstract property perturbation_ic

Set of well scaling initial conditions for perturbation module.

abstract V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

abstract dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

abstract d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

abstract d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

abstract d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

abstract inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

abstract get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

abstract get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

abstract get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

abstract get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

get_epsilon_1(phi)[source]

Approximation of 1st Hubble flow parameter with potential slow-roll parameters.

get_epsilon_2(phi)[source]

Approximation of 2nd Hubble flow parameter with potential slow-roll parameters.

get_epsilon_3(phi)[source]

Approximation of 3rd Hubble flow parameter with potential slow-roll parameters.

get_epsilon_4(phi)[source]

Approximation of 4th Hubble flow parameter with potential slow-roll parameters.

abstract property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

abstract sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

abstract sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_As2Lambda(A_s, phi_star, N_star, **pot_kwargs)[source]

Get potential amplitude Lambda from PPS amplitude A_s using slow-roll approximation.

Find the inflaton amplitude Lambda (4th root of potential amplitude) that produces the desired amplitude A_s of the primordial power spectrum using the slow-roll approximation.

Note that either phi_star or N_star have to be passed as None, and will be calculated subsequently from the respective other.

Parameters:
A_sfloat

Target amplitude A_s of scalar primordial power spectrum (at the pivot scale k=0.05 Mpc^{-1}).

phi_starfloat

Inflaton field value at the time of horizon crossing of the pivot scale.

N_starfloat

Number of e-folds of inflation remaining at the time of horizon crossing of the pivot scale.

Returns:
Lambdafloat

Potential amplitude parameter Lambda corresponding approximately to A_s.

phi_starfloat

Estimated inflaton field value at the time of horizon crossing of the pivot scale, inferred from N_star. (Exact if passed as input.)

N_starfloat

Estimated number of e-folds of inflation remaining at the time of horizon crossing of the pivot scale, inferred from phi_star. (Exact if passed as input.)

class primpy.potentials.MonomialPotential(**pot_kwargs)[source]

Bases: InflationaryPotential

Monomial potential: V(phi) = Lambda**4 * phi**p.

Attributes:
phi_end

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_Nstar2ns(N_star, **pot_kwargs)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star, **pot_kwargs)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s, **pot_kwargs)

Slow-roll approximation for inferring N_star from n_s.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_r2Nstar(r, **pot_kwargs)

Slow-roll approximation for inferring N_star from r.

tag = 'mnp'
name = 'MonomialPotential'
tex = '$\\phi^p$'
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

static sr_Nstar2ns(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from r.

class primpy.potentials.LinearPotential(**pot_kwargs)[source]

Bases: MonomialPotential

Linear potential: V(phi) = Lambda**4 * phi.

Methods

sr_Nstar2ns(N_star, **pot_params)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star, **pot_params)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s, **pot_params)

Slow-roll approximation for inferring N_star from n_s.

sr_r2Nstar(r, **pot_params)

Slow-roll approximation for inferring N_star from r.

tag = 'mn1'
name = 'LinearPotential'
tex = '$\\phi^1$'
perturbation_ic = (1, 0, 0, 1)
static sr_Nstar2ns(N_star, **pot_params)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s, **pot_params)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_params)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r, **pot_params)[source]

Slow-roll approximation for inferring N_star from r.

class primpy.potentials.QuadraticPotential(**pot_kwargs)[source]

Bases: MonomialPotential

Quadratic potential: V(phi) = Lambda**4 * phi**2.

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

inv_V(V)

Inverse function phi(V) with respect to potential V.

sr_Nstar2ns(N_star, **pot_kwargs)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star, **pot_kwargs)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s, **pot_kwargs)

Slow-roll approximation for inferring N_star from n_s.

sr_r2Nstar(r, **pot_kwargs)

Slow-roll approximation for inferring N_star from r.

tag = 'mn2'
name = 'QuadraticPotential'
tex = '$\\phi^2$'
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

static sr_Nstar2ns(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from r.

class primpy.potentials.CubicPotential(**pot_kwargs)[source]

Bases: MonomialPotential

Linear potential: V(phi) = Lambda**4 * phi.

Methods

sr_Nstar2ns(N_star, **pot_kwargs)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star, **pot_kwargs)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s, **pot_kwargs)

Slow-roll approximation for inferring N_star from n_s.

sr_r2Nstar(r, **pot_kwargs)

Slow-roll approximation for inferring N_star from r.

tag = 'mn3'
name = 'CubicPotential'
tex = '$\\phi^3$'
perturbation_ic = (1, 0, 0, 1)
static sr_Nstar2ns(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from r.

class primpy.potentials.QuarticPotential(**pot_kwargs)[source]

Bases: MonomialPotential

Linear potential: V(phi) = Lambda**4 * phi.

Methods

sr_Nstar2ns(N_star, **pot_kwargs)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star, **pot_kwargs)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s, **pot_kwargs)

Slow-roll approximation for inferring N_star from n_s.

sr_r2Nstar(r, **pot_kwargs)

Slow-roll approximation for inferring N_star from r.

tag = 'mn4'
name = 'QuarticPotential'
tex = '$\\phi^4$'
perturbation_ic = (1, 0, 0, 1)
static sr_Nstar2ns(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_kwargs)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from r.

class primpy.potentials.StarobinskyPotential(**pot_kwargs)[source]

Bases: InflationaryPotential

Starobinsky potential: V(phi) = Lambda**4 * (1 - exp(-sqrt(2/3) * phi))**2.

Attributes:
phi_end

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

phi2efolds(phi)

Get e-folds N from inflaton phi.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_Nstar2ns(N_star)

Slow-roll approximation for inferring n_s from N_star.

sr_Nstar2r(N_star)

Slow-roll approximation for inferring r from N_star.

sr_ns2Nstar(n_s)

Slow-roll approximation for inferring N_star from n_s.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_r2Nstar(r)

Slow-roll approximation for inferring N_star from r.

tag = 'stb'
name = 'StarobinskyPotential'
tex = 'Starobinsky'
gamma = np.float64(0.816496580927726)
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

static sr_Nstar2ns(N_star)[source]

Slow-roll approximation for inferring n_s from N_star.

static sr_ns2Nstar(n_s)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star)[source]

Slow-roll approximation for inferring r from N_star.

static sr_r2Nstar(r)[source]

Slow-roll approximation for inferring N_star from r.

static phi2efolds(phi)[source]

Get e-folds N from inflaton phi.

Find the number of e-folds N till end of inflation from inflaton phi using the slow-roll approximation.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Nfloat or np.ndarray

Number of e-folds N until end of inflation.

class primpy.potentials.NaturalPotential(**pot_kwargs)[source]

Bases: InflationaryPotential

Natural inflation potential: V(phi) = Lambda**4 * (1 - cos(pi*phi/phi0)).

Natural inflation with phi0 = pi * f where f is the standard parameter used in definitions of natural inflation. Here we use phi0 the position of the maximum and we have a minus in our definition such that the minimum is at zero instead of the maximum.

Attributes:
phi_end

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

phi2efolds(phi, phi0)

Get e-folds N from inflaton phi.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_Nstar2ns(N_star, **pot_kwargs)

Slow-roll approximation for the spectral index n_s.

sr_Nstar2r(N_star, **pot_kwargs)

Slow-roll approximation for the tensor-to-scalar ratio r.

sr_ns2Nstar(n_s, **pot_kwargs)

Slow-roll approximation for inferring N_star from n_s.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_r2Nstar(r, **pot_kwargs)

Slow-roll approximation for inferring N_star from r.

tag = 'nat'
name = 'NaturalPotential'
tex = 'Natural'
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

static sr_Nstar2ns(N_star, **pot_kwargs)[source]

Slow-roll approximation for the spectral index n_s.

static sr_ns2Nstar(n_s, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from n_s.

static sr_Nstar2r(N_star, **pot_kwargs)[source]

Slow-roll approximation for the tensor-to-scalar ratio r.

static sr_r2Nstar(r, **pot_kwargs)[source]

Slow-roll approximation for inferring N_star from r.

static phi2efolds(phi, phi0)[source]

Get e-folds N from inflaton phi.

Find the number of e-folds N till end of inflation from inflaton phi using the slow-roll approximation.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

phi0float

Inflaton distance between local maximum and minima.

Returns:
Nfloat or np.ndarray

Number of e-folds N until end of inflation.

class primpy.potentials.DoubleWellPotential(**pot_kwargs)[source]

Bases: InflationaryPotential

Double-Well potential: V(phi) = Lambda**4 * (1 - (phi/phi0)**p)**2.

Double-Well shifted such that left minimum is at zero: phi -> phi-phi0

Attributes:
phi_end

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

tag = 'dwp'
name = 'DoubleWellPotential'
tex = 'Double-Well (p)'
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

class primpy.potentials.DoubleWell2Potential(**pot_kwargs)[source]

Bases: DoubleWellPotential

Quadratic Double-Well potential: V(phi) = Lambda**4 * (1 - (phi/phi0)**2)**2.

Double-Well shifted such that left minimum is at zero: phi -> phi-phi0

Attributes:
phi_end

Methods

phi2efolds(phi_shifted, phi0)

Get e-folds N from inflaton phi.

tag = 'dw2'
name = 'DoubleWell2Potential'
tex = 'Double-Well (quadratic)'
perturbation_ic = (1, 0, 0, 1)
property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

static phi2efolds(phi_shifted, phi0)[source]

Get e-folds N from inflaton phi.

Find the number of e-folds N till end of inflation from inflaton phi using the slow-roll approximation.

Parameters:
phi_shiftedfloat or np.ndarray

Inflaton field phi shifted by phi0 such that left potential minimum is at zero.

phi0float

Inflaton distance between local maximum and minima.

Returns:
Nfloat or np.ndarray

Number of e-folds N until end of inflation.

class primpy.potentials.DoubleWell4Potential(**pot_kwargs)[source]

Bases: DoubleWellPotential

Quartic Double-Well potential: V(phi) = Lambda**4 * (1 - (phi/phi0)**4)**2.

Double-Well shifted such that left minimum is at zero: phi -> phi-phi0

Methods

phi2efolds(phi_shifted, phi0)

Get e-folds N from inflaton phi.

phi_end_squared(phi0)

Get inflaton at end of inflation using slow-roll.

tag = 'dw4'
name = 'DoubleWell4Potential'
tex = 'Double-Well (quartic)'
perturbation_ic = (1, 0, 0, 1)
static phi_end_squared(phi0)[source]

Get inflaton at end of inflation using slow-roll.

Parameters:
phi0float

Inflaton distance between local maximum and minima.

Returns:
phi_end2float

Inflaton phi squared at end of inflation. (unshifted!)

classmethod phi2efolds(phi_shifted, phi0)[source]

Get e-folds N from inflaton phi.

Find the number of e-folds N till end of inflation from inflaton phi using the slow-roll approximation.

Parameters:
phi_shiftedfloat or np.ndarray

Inflaton field phi shifted by phi0 such that left potential minimum is at zero.

phi0float

Inflaton distance between local maximum and minima.

Returns:
Nfloat or np.ndarray

Number of e-folds N until end of inflation.

class primpy.potentials.TmodelPotential(**pot_kwargs)[source]

Bases: InflationaryPotential

T-model potential: V(phi) = Lambda**4 * (tanh(phi / (sqrt(6) * alpha)))**(2*p).

Attributes:
phi_end

Methods

V(phi)

Inflationary potential V(phi).

d2V(phi)

Second derivative V''(phi) with respect to inflaton phi.

d3V(phi)

Third derivative V'''(phi) with respect to inflaton phi.

d4V(phi)

Fourth derivative V''''(phi) with respect to inflaton phi.

dV(phi)

First derivative V'(phi) with respect to inflaton phi.

get_epsilon_1V(phi)

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)

Get 4th potential slow-roll parameter.

inv_V(V)

Inverse function phi(V) with respect to potential V.

sr_N2phi(N)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_phi2N(phi)

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

tag = 'tmp'
name = 'TmodelPotential'
tex = 'T-model'
perturbation_ic = (1, 0, 0, 1)
V(phi)[source]

Inflationary potential V(phi).

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
Vfloat or np.ndarray

Inflationary potential V(phi).

dV(phi)[source]

First derivative V’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
dVfloat or np.ndarray

1st derivative of inflationary potential: V’(phi).

d2V(phi)[source]

Second derivative V’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d2Vfloat or np.ndarray

2nd derivative of inflationary potential: V’’(phi).

d3V(phi)[source]

Third derivative V’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d3Vfloat or np.ndarray

3rd derivative of inflationary potential: V’’’(phi).

d4V(phi)[source]

Fourth derivative V’’’’(phi) with respect to inflaton phi.

Parameters:
phifloat or np.ndarray

Inflaton field phi.

Returns:
d4Vfloat or np.ndarray

4th derivative of inflationary potential: V’’’’(phi).

inv_V(V)[source]

Inverse function phi(V) with respect to potential V.

Parameters:
Vfloat or np.ndarray

Inflationary potential V.

Returns:
phifloat or np.ndarray

Inflaton field phi.

get_epsilon_1V(phi)[source]

Get 1st potential slow-roll parameter.

get_epsilon_2V(phi)[source]

Get 2nd potential slow-roll parameter.

get_epsilon_3V(phi)[source]

Get 3rd potential slow-roll parameter.

get_epsilon_4V(phi)[source]

Get 4th potential slow-roll parameter.

property phi_end

Inflaton field value at the end of inflation.

Value inferred from epsilon_1V = 1 assuming a positive phi_end.

sr_phi2N(phi)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

sr_N2phi(N)[source]

Convert from inflaton field phi to e-folds N assuming slow-roll approximation.

class primpy.potentials.FeatureFunction[source]

Bases: ABC

Feature in the inflationary potential.

Methods

F(x, x0, a, b)

Feature function.

d2F(x, x0, a, b)

Feature function 2nd derivative.

d3F(x, x0, a, b)

Feature function 3rd derivative.

dF(x, x0, a, b)

Feature function derivative.

abstract static F(x, x0, a, b)[source]

Feature function.

abstract static dF(x, x0, a, b)[source]

Feature function derivative.

abstract static d2F(x, x0, a, b)[source]

Feature function 2nd derivative.

abstract static d3F(x, x0, a, b)[source]

Feature function 3rd derivative.

class primpy.potentials.GaussianDip[source]

Bases: FeatureFunction

Gaussian: F(x) = -a * exp(-(x-x0)**2 / (2*b**2)).

Methods

F(x, x0, a, b)

F(x) = -a * exp(-(x-x0)**2 / (2*b**2)).

d2F(x, x0, a, b)

F''(x) = a/b**4 * (b**2 - (x-x0)**2) * exp(-(x-x0)**2 / (2*b**2)).

d3F(x, x0, a, b)

F'''(x) = a/b**6 * (x-x0) * ((x-x0)**2 - 3*b**2) * exp(-(x-x0)**2 / (2*b**2)).

dF(x, x0, a, b)

F'(x) = a/b**2 * (x-x0) * exp(-(x-x0)**2 / (2*b**2)).

static F(x, x0, a, b)[source]

F(x) = -a * exp(-(x-x0)**2 / (2*b**2)).

static dF(x, x0, a, b)[source]

F’(x) = a/b**2 * (x-x0) * exp(-(x-x0)**2 / (2*b**2)).

static d2F(x, x0, a, b)[source]

F’’(x) = a/b**4 * (b**2 - (x-x0)**2) * exp(-(x-x0)**2 / (2*b**2)).

static d3F(x, x0, a, b)[source]

F’’’(x) = a/b**6 * (x-x0) * ((x-x0)**2 - 3*b**2) * exp(-(x-x0)**2 / (2*b**2)).

class primpy.potentials.TanhStep[source]

Bases: FeatureFunction

Tanh step function: F(x) = a * tanh((x - x0) / b).

Methods

F(x, x0, a, b)

F(x) = a * tanh((x-x0)/b).

d2F(x, x0, a, b)

F''(x) = -2*a/b**2 * tanh((x-x0)/b) * (1 - tanh((x-x0)/b)**2).

d3F(x, x0, a, b)

F'''(x) = -2*a/b**3 * (1 - 4*tanh((x-x0)/b)**2 + 3*tanh((x-x0)/b)**4).

dF(x, x0, a, b)

F'(x) = a/b * (1 - tanh((x-x0)/b)**2).

static F(x, x0, a, b)[source]

F(x) = a * tanh((x-x0)/b).

static dF(x, x0, a, b)[source]

F’(x) = a/b * (1 - tanh((x-x0)/b)**2).

static d2F(x, x0, a, b)[source]

F’’(x) = -2*a/b**2 * tanh((x-x0)/b) * (1 - tanh((x-x0)/b)**2).

static d3F(x, x0, a, b)[source]

F’’’(x) = -2*a/b**3 * (1 - 4*tanh((x-x0)/b)**2 + 3*tanh((x-x0)/b)**4).

class primpy.potentials.FeaturePotential(**pot_kwargs)[source]

Bases: InflationaryPotential, FeatureFunction

Inflationary potential with a feature: V(phi) = V0(phi) * (1+F(phi)).

Methods

V(phi)

Inflationary potential V0(phi) with a feature function F(phi).

d2V(phi)

Second derivative of the inflationary potential with a feature.

d3V(phi)

Third derivative of the inflationary potential with a feature.

dV(phi)

First derivative of the inflationary potential with a feature.

V(phi)[source]

Inflationary potential V0(phi) with a feature function F(phi).

V(phi) = V0(phi) * (1 + F(phi))

dV(phi)[source]

First derivative of the inflationary potential with a feature.

V’(phi) = V0’(phi) * (1 + F(phi)) + V0(phi) * F’(phi)

d2V(phi)[source]

Second derivative of the inflationary potential with a feature.

V’’(phi) = V0’’(phi) * (1 + F(phi)) + 2*V0’(phi)*F’(phi) + V0(phi)*F’’(phi)

d3V(phi)[source]

Third derivative of the inflationary potential with a feature.

\[V'''(\phi) = V_0'''(\phi) * (1 + F(\phi)) + 3 * V_0''(\phi) * F'(\phi) + 3 * V_0'(\phi) * F''(\phi) + V_0(\phi) * F'''(\phi)\]
class primpy.potentials.StarobinskyGaussianDipPotential(**pot_kwargs)[source]

Bases: FeaturePotential, StarobinskyPotential, GaussianDip

Starobinsky potential with a Gaussian dip.

tag = 'sgd'
name = 'StarobinskyGaussianDipPotential'
tex = 'Starobinsky with a Gaussian dip'
class primpy.potentials.StarobinskyTanhStepPotential(**pot_kwargs)[source]

Bases: FeaturePotential, StarobinskyPotential, TanhStep

Starobinsky potential with a hyperbolic tangent step.

tag = 'sts'
name = 'StarobinskyTanhStepPotential'
tex = 'Starobinsky with a tanh step'

primpy.solver module

General setup for running scipy.integrate.solve_ivp().

(Modified from “primordial” by Will Handley.)

primpy.solver.solve(ic, *args, **kwargs)[source]

Run scipy.integrate.solve_ivp() and store information in sol for post-processing.

This is a wrapper around scipy.integrate.solve_ivp(), with easier reusable objects for the equations and initial conditions.

Parameters:
icprimpy.initialconditions.InitialConditions

Initial conditions specifying relevant equations, variables, and initial numerical values.

*args, **kwargs

All other arguments are identical to scipy.integrate.solve_ivp().

Returns:
solBunch object same as returned by scipy.integrate.solve_ivp()

Solution to the inverse value problem.

primpy.units module

Units and constants for primpy.

primpy.units.mp_kg = np.float64(4.341358399139358e-09)

reduced Planck mass in kg

primpy.units.tp_s = np.float64(2.7027701565693675e-43)

reduced Planck time in s

primpy.units.lp_m = np.float64(8.102701086469756e-35)

reduced Planck length in m

primpy.units.Tp_K = np.float64(2.826075522438417e+31)

reduced Planck temperature in K

primpy.units.mp_GeV = np.float64(2.4353234600842885e+18)

reduced Planck mass in GeV

primpy.units.tp_iGeV = np.float64(4.106230717973657e-19)

reduced Planck time in GeV^-1

primpy.units.lp_iGeV = np.float64(4.106230717973658e-19)

reduced Planck length in GeV^-1

primpy.units.Mpc_m = 3.085677581491367e+22

conversion factor from Mpc to m

primpy.units.a_B = 7.565733250280007e-16

radiation constant (Planck’s law)