primpy.time package
Calculations for the primordial Universe w.r.t. time t.
primpy.time.inflation module
Differential equations for inflation w.r.t. time t.
- class primpy.time.inflation.InflationEquationsT(K, potential, track_eta=False, verbose=False)[source]
Bases:
InflationEquationsBackground equations during inflation w.r.t. time t.
Solves background variables in cosmic time for curved and flat universes using the Klein-Gordon and Friedmann equations.
- Independent variable:
t: cosmic time- Dependent variables:
_N: number of e-foldsphi: inflaton fielddphidt: d(phi)/dteta: conformal time (optional)
Methods
H2(x, y)Hubble parameter squared.
__call__(x, y)System of coupled ODEs for underlying variables.
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(N, 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.
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.
- static get_dH_H(N, 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
- sol(sol, **kwargs)[source]
Post-processing of
scipy.integrate.solve_ivp()solution.
primpy.time.perturbations module
Curvature perturbations with respect to time t.
- class primpy.time.perturbations.PerturbationT(background, k, **kwargs)[source]
Bases:
PerturbationCurvature perturbation for wavenumber k with respect to time t.
Solves the Mukhanov–Sasaki equations w.r.t. cosmic time for curved universes.
- Parameters:
- backgroundBunch object same as returned by
scipy.integrate.solve_ivp() Background solution as returned by
primpy.time.inflation.InflationEquationsT.sol().- kfloat
wavenumber
- backgroundBunch object same as returned by
- class primpy.time.perturbations.ScalarModeT(background, k, **kwargs)[source]
Bases:
ScalarModeTemplate for scalar modes.
Methods
__call__(x, y)Vector of derivatives.
Get initial conditions for scalar modes for HD vacuum w.r.t.
Get initial conditions for scalar modes for RST vacuum w.r.t.
Get initial conditions for scalar modes for HD approximation w.r.t.
Frequency and damping term of the Mukhanov-Sasaki equations for scalar modes.
- mukhanov_sasaki_frequency_damping()[source]
Frequency and damping term of the Mukhanov-Sasaki equations for scalar modes.
Frequency and damping term of the Mukhanov-Sasaki equations for the comoving curvature perturbations R w.r.t. time t, where the e.o.m. is written as ddR + 2 * damping * dR + frequency**2 * R = 0.
- get_vacuum_ic_k()[source]
Get initial conditions for scalar modes for HD approximation w.r.t. cosmic time t.
- class primpy.time.perturbations.TensorModeT(background, k, **kwargs)[source]
Bases:
TensorModeTemplate for tensor modes.
Methods
__call__(x, y)Vector of derivatives.
Get initial conditions for tensor modes for HD vacuum w.r.t.
Get initial conditions for tensor modes for RST vacuum w.r.t.
Get initial conditions for tensor modes for HD approximation w.r.t.
Frequency and damping term of the Mukhanov-Sasaki equations for tensor modes.
- mukhanov_sasaki_frequency_damping()[source]
Frequency and damping term of the Mukhanov-Sasaki equations for tensor modes.
Frequency and damping term of the Mukhanov-Sasaki equations for the tensor perturbations h w.r.t. time t, where the e.o.m. is written as ddh + 2 * damping * dh + frequency**2 * h = 0.
- get_vacuum_ic_k()[source]
Get initial conditions for tensor modes for HD approximation w.r.t. cosmic time t.