primpy API

Subpackages

primpy.potentials module

primpy.potentials: inflationary potentials.

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

Bases: ABC

Base class for inflaton potential and derivatives.

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:

phi (float or np.ndarray) – Inflaton field phi.

Returns:

V – Inflationary potential V(phi).

Return type:

float or np.ndarray

abstract dV(phi)[source]

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

Parameters:

phi (float or np.ndarray) – Inflaton field phi.

Returns:

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

Return type:

float or np.ndarray

abstract d2V(phi)[source]

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

Parameters:

phi (float or np.ndarray) – Inflaton field phi.

Returns:

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

Return type:

float or np.ndarray

abstract d3V(phi)[source]

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

Parameters:

phi (float or np.ndarray) – Inflaton field phi.

Returns:

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

Return type:

float or np.ndarray

abstract inv_V(V)[source]

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

Parameters:

V (float or np.ndarray) – Inflationary potential V.

Returns:

phi – Inflaton field phi.

Return type:

float or np.ndarray

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

Get potential amplitude Lambda from PPS amplitude A_s.

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

Bases: InflationaryPotential

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

V(phi)[source]

V(phi) = Lambda**4 * phi**p.

dV(phi)[source]

V(phi) = Lambda**4 * phi**(p-1) * p.

d2V(phi)[source]

V(phi) = Lambda**4 * phi**(p-2) * p * (p-1).

d3V(phi)[source]

V(phi) = Lambda**4 * phi**(p-3) * p * (p-1) * (p-2).

inv_V(V)[source]

phi(V) = (V / Lambda**4)**(1/p).

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.

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Monomial potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

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

Bases: MonomialPotential

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

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.

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Linear potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

V(phi)

V(phi) = Lambda**4 * phi**p.

d2V(phi)

V(phi) = Lambda**4 * phi**(p-2) * p * (p-1).

d3V(phi)

V(phi) = Lambda**4 * phi**(p-3) * p * (p-1) * (p-2).

dV(phi)

V(phi) = Lambda**4 * phi**(p-1) * p.

inv_V(V)

phi(V) = (V / Lambda**4)**(1/p).

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

Bases: MonomialPotential

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

V(phi)[source]

V(phi) = Lambda**4 * phi**2.

dV(phi)[source]

V’(phi) = 2 * Lambda**4 * phi.

d2V(phi)[source]

V’’(phi) = 2 * Lambda**4.

d3V(phi)[source]

V’’’(phi) = 0.

inv_V(V)[source]

phi(V) = sqrt(V) / Lambda**2.

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.

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

Get potential amplitude Lambda from PPS amplitude A_s.

Find the inflaton mass m (i.e. essentially the potential amplitude) that produces the desired amplitude A_s of the primordial power spectrum using the slow roll approximation.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Quadratic potential (Lambda**2 = mass).

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

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

Bases: MonomialPotential

Linear potential: V(phi) = Lambda**4 * 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.

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Linear potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

V(phi)

V(phi) = Lambda**4 * phi**p.

d2V(phi)

V(phi) = Lambda**4 * phi**(p-2) * p * (p-1).

d3V(phi)

V(phi) = Lambda**4 * phi**(p-3) * p * (p-1) * (p-2).

dV(phi)

V(phi) = Lambda**4 * phi**(p-1) * p.

inv_V(V)

phi(V) = (V / Lambda**4)**(1/p).

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

Bases: MonomialPotential

Linear potential: V(phi) = Lambda**4 * 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.

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Linear potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

V(phi)

V(phi) = Lambda**4 * phi**p.

d2V(phi)

V(phi) = Lambda**4 * phi**(p-2) * p * (p-1).

d3V(phi)

V(phi) = Lambda**4 * phi**(p-3) * p * (p-1) * (p-2).

dV(phi)

V(phi) = Lambda**4 * phi**(p-1) * p.

inv_V(V)

phi(V) = (V / Lambda**4)**(1/p).

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

Bases: InflationaryPotential

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

V(phi)[source]

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

dV(phi)[source]

V’(phi) = Lambda**4 * 2 * gamma * exp(-2 * gamma * phi) * (-1 + exp(gamma * phi)).

d2V(phi)[source]

V’’(phi) = Lambda**4 * 2 * gamma**2 * exp(-2*gamma*phi) * (2 - exp(gamma*phi)).

d3V(phi)[source]

V’’’(phi) = Lambda**4 * 2 * gamma**3 * exp(-2*gamma*phi) * (-4 + exp(gamma*phi)).

inv_V(V)[source]

phi(V) = -np.log(1 - np.sqrt(V) / Lambda**2) / gamma.

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:

phi (float or np.ndarray) – Inflaton field phi.

Returns:

N – Number of e-folds N until end of inflation.

Return type:

float or np.ndarray

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Starobinsky potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the 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.

V(phi)[source]

V(phi) = Lambda**4 * (1 - cos(pi*phi/phi0)).

dV(phi)[source]

V(phi) = Lambda**4 * sin(pi*phi/phi0) * pi / phi0.

d2V(phi)[source]

V(phi) = Lambda**4 * cos(pi*phi/phi0) * (pi / phi0)**2.

d3V(phi)[source]

V(phi) = -Lambda**4 * sin(pi*phi/phi0) * (pi / phi0)**3.

inv_V(V)[source]

phi(V) = arccos(1 - V / Lambda**4) * phi0 / pi.

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:
  • phi (float or np.ndarray) – Inflaton field phi.

  • phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

N – Number of e-folds N until end of inflation.

Return type:

float or np.ndarray

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Keyword Arguments:

phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Natural inflation potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the 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

V(phi)[source]

V(phi) = Lambda**4 * (1 - (phi/phi0)**p)**2.

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

dV(phi)[source]

V’(phi) = 2p*Lambda**4 * (-1 + (phi / phi0)**p) * phi**(p - 1) / phi0**p.

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

d2V(phi)[source]

V’’(phi) = 2p*Lambda**4 * (1-p+(2*p-1)*(phi/phi0)**p) * phi**(p-2) / phi0**p.

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

d3V(phi)[source]

V’’’(phi) = 2p(p-1)Lambda**4 * (2-p+(4*p-2)*(phi/phi0)**p) * phi**(p-3) / phi0**p.

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

inv_V(V)[source]

phi(V) = phi0 * (1 - sqrt(V) / Lambda**2)**(1/p).

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

Get potential amplitude Lambda from PPS amplitude A_s.

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

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_shifted (float or np.ndarray) – Inflaton field phi shifted by phi0 such that left potential minimum is at zero.

  • phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

N – Number of e-folds N until end of inflation.

Return type:

float or np.ndarray

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Keyword Arguments:

phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Double-Well potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

V(phi)

V(phi) = Lambda**4 * (1 - (phi/phi0)**p)**2.

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

d2V(phi)

V’’(phi) = 2p*Lambda**4 * (1-p+(2*p-1)*(phi/phi0)**p) * phi**(p-2) / phi0**p.

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

d3V(phi)

V’’’(phi) = 2p(p-1)Lambda**4 * (2-p+(4*p-2)*(phi/phi0)**p) * phi**(p-3) / phi0**p.

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

dV(phi)

V’(phi) = 2p*Lambda**4 * (-1 + (phi / phi0)**p) * phi**(p - 1) / phi0**p.

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

inv_V(V)

phi(V) = phi0 * (1 - sqrt(V) / Lambda**2)**(1/p).

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

static phi_end_squared(phi0)[source]

Get inflaton at end of inflation using slow-roll.

Parameters:

phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

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

Return type:

float

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_shifted (float or np.ndarray) – Inflaton field phi shifted by phi0 such that left potential minimum is at zero.

  • phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

N – Number of e-folds N until end of inflation.

Return type:

float or np.ndarray

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

Get potential amplitude Lambda from PPS amplitude A_s.

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.

Parameters:
  • A_s (float or np.ndarray) – Amplitude A_s of the primordial power spectrum.

  • phi_star (float or None) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float or None) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

Keyword Arguments:

phi0 (float) – Inflaton distance between local maximum and minima.

Returns:

  • Lambda (float or np.ndarray) – Amplitude parameter Lambda for the Double-Well potential.

  • phi_star (float) – Inflaton value at horizon crossing of the pivot scale.

  • N_star (float) – Number of observable e-folds of inflation N_star from horizon crossing till the end of inflation.

V(phi)

V(phi) = Lambda**4 * (1 - (phi/phi0)**p)**2.

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

d2V(phi)

V’’(phi) = 2p*Lambda**4 * (1-p+(2*p-1)*(phi/phi0)**p) * phi**(p-2) / phi0**p.

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

d3V(phi)

V’’’(phi) = 2p(p-1)Lambda**4 * (2-p+(4*p-2)*(phi/phi0)**p) * phi**(p-3) / phi0**p.

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

dV(phi)

V’(phi) = 2p*Lambda**4 * (-1 + (phi / phi0)**p) * phi**(p - 1) / phi0**p.

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

inv_V(V)

phi(V) = phi0 * (1 - sqrt(V) / Lambda**2)**(1/p).

primpy.events module

primpy.events: 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 [-1, 0, +1], optional, default 0) – The direction of the root finding (if any)

  • terminal (bool, optional, default False) – Whether to stop at this root

  • value (float, optional, default 0) – Offset to root

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

Bases: Event

Stop after a given amount of time t has passed.

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

Bases: Event

Stop after a given number of e-folds N.

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

Bases: Event

Track inflation start/end.

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

Bases: Event

Go a bit past the end of inflation.

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.

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

Bases: Event

Track zero crossings of inflaton phi.

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

Bases: Event

Track when mode exits the horizon aH.

primpy.initialconditions module

primpy.initialconditions: initial conditions for inflation.

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

Bases: object

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.

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

Bases: object

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.

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.

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.

primpy.equations module

primpy.equations: 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.

idx

dictionary mapping variable names to indices in the solution vector y

Type:

dict

independent_variable

name of independent variable

Type:

string

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:

*args (str) – Name of the dependent variables

primpy.inflation module

primpy.inflation: general setup for equations for cosmic inflation.

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

Bases: Equations, ABC

Base class for inflation equations.

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.solver module

primpy.solver: general setup for running scipy.integrate.solve_ivp().

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:

ic (primordial.initialconditions.InitialConditions) – Initial conditions specifying relevant equations, variables, and initial numerical values.

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

Returns:

sol – Solution to the inverse value problem. Monkey-patched version of the Bunch type usually returned by scipy.integrate.solve_ivp().

Return type:

Bunch object

(c) modified from “primordial” by Will Handley.

primpy.perturbations module

primpy.time.perturbations: comoving curvature perturbations w.r.t. time t.

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

Bases: object

Primordial Power spectrum of curvature perturbations.

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

Bases: ABC

Perturbation for wavenumber k.

oscode_postprocessing(oscode_sol, **kwargs)[source]

Post-processing for pyoscode.solve() solution.

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

Parameters:

oscode_sol (list) – 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.

abstract mukhanov_sasaki_frequency_damping()[source]

Frequency and damping term of the Mukhanov-Sasaki equations.

abstract get_vacuum_ic_RST()[source]

Get initial conditions for scalar modes for RST vacuum.

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

Bases: Mode, ABC

Template for scalar modes.

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

Bases: Mode, ABC

Template for tensor modes.

primpy.oscode_solver module

primpy.oscode_solver: setup for running pyoscode.solve().

primpy.oscode_solver.solve_oscode(background, k, **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:
  • background (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.

  • k (int, float, np.ndarray) – Comoving wavenumber used to evolve the Mukhanov-Sasaki equation.

Keyword Arguments:
  • 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

  • rtol (float) – Tolerance passed to pyoscode. default : 5e-5

  • fac_beg (int, 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_end (int, 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

  • even_grid (bool) – Set this to True if the grid of the independent variable is equally spaced. default : False

  • vacuum (tuple) – Set of vacuum initial conditions to be computed. Choose any of (‘RST’, ). default : (‘RST’, )

  • drop_closed_large_scales (bool) – 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

Returns:

sol – 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().

Return type:

Bunch object

primpy.bigbang module

primpy.bigbang: 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_reh(N1, N2, log_cHH1, log_cHH2)[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)[source]

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

Parameters:
  • N (float, 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_m0 (float) – matter density parameter today

  • Omega_K0 (float) – curvature density parameter today

  • h (float) – dimensionless Hubble parameter today, “little h”

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.

Returns:

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

Return type:

float

primpy.bigbang.no_Big_Bang_line(Omega_m0)[source]

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

Parameters:

Omega_m0 (float) – matter density parameter today

Returns:

Omega_L0 – 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).

Return type:

float

primpy.bigbang.expand_recollapse_line(Omega_m0)[source]

Return Omega_L0 for dividing line between expanding/recollapsing universes.

Parameters:

Omega_m0 (float) – matter density parameter today

Returns:

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

Return type:

float

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:
  • N (float, 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_m0 (float) – matter density parameter today

  • Omega_K0 (float) – curvature density parameter today

  • h (float) – dimensionless Hubble parameter today, “little h”

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

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.

Returns:

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

Return type:

float

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_start (float) – 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).

  • N (float, 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_m0 (float) – matter density parameter today

  • Omega_K0 (float) – curvature density parameter today

  • h (float) – dimensionless Hubble parameter today, “little h”

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.

Returns:

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

Return type:

float, np.ndarray

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_m0 (float) – matter density parameter today

  • Omega_K0 (float) – curvature density parameter today

  • h (float) – dimensionless Hubble parameter today, “little h”

  • b_forward (Bunch 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_backward (Bunch object same as returned by scipy.integrate.solve_ivp()) – Additional solution returned by primpy.solver.solve(). This second solution is assumed to be an integration from inflation start backwards in time. optional, default : None

Returns:

ratio – 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.

Return type:

float