analysis Package

analysis Package

analysis package - Provides modules to analyse results from cosmomodels runs.

adiabatic Module

pyflation.analysis.adiabatic - Tools to calculate the spectra of adiabatic perturbations

pyflation.analysis.adiabatic.Pphi_matrix(modes, axis)[source]

Return the cross correlation of scalar perturbations P_{I,J}

For multifield systems the full crossterm matrix is returned which has shape nfields*nfields.

Parameters :

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pphi_matrix : array_like, dtype: float64

array of Pphi values with the same shape as input variable modes.

pyflation.analysis.adiabatic.Pphi_modes(m)[source]

Return the modes of the scalar perturbations P_{I,J}.

This is a helper function which wraps the full calculation in Pphi_matrix and requires only the model as an argument. Provided for compatibility with previous versions.

Parameters :

m : Cosmomodels instance

Model class instance from which the yresult variable will be used to calculate P_{I,J}.

Returns :

Pphi_modes : array_like, dtype: float64

array of Pphi values flattened to have the same shape as m.yresult[:,m.dps_ix].

pyflation.analysis.adiabatic.Pr(m, tix=None, kix=None)[source]

Return the spectrum of (first order) curvature perturbations P_R1 for a model m.

For a multifield model this is given by:

P_\mathcal{R} = (\sum_K \dot{\phi_K}^2 )^{-2} 
    \sum_{I,J} \dot{\phi_I} \dot{\phi_J} P_{IJ}

where P_{IJ} = \sum_K \chi_{IK} \chi_{JK} and \chi are the mode matrix elements.

This is the unscaled version P_R which is related to the scaled version by \mathcal{P}_\mathcal{R} = k^3/(2\pi^2) P_\mathcal{R}.

This function is a wrapper function which only requires the model as a parameter. The full calculation is done in the Pr_spectrum function in the pyflation.analysis.adiabatic module.

Parameters :

m : Cosmomodels instance

Model class instance from which the yresult variable will be used to calculate P_R.

tix : integer, optional

index of timestep at which to calculate, defaults to full range of steps.

kix : integer, optional

integer of k value at which to calculate, defaults to full range of ks.

Returns :

Pr : array_like, dtype: float64

Array of Pr values for requested timesteps and kmodes

pyflation.analysis.adiabatic.Pr_spectrum(phidot, modes, axis)[source]

Return the spectrum of (first order) curvature perturbations P_R1 for each k, given the first order perturbation modes.

For a multifield model this is given by:

P_\mathcal{R} = (\sum_K \dot{\phi_K}^2 )^{-2} \sum_{I,J} \dot{\phi_I} \dot{\phi_J} P_{IJ}

where P_{IJ} = \sum_K \chi_{IK} \chi_{JK} and \chi are the mode matrix elements.

This is the unscaled version P_R which is related to the scaled version by \mathcal{P}_\mathcal{R} = k^3/(2\pi^2) P_\mathcal{R}.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pr : array_like, dtype: float64

Array of Pr values for all timesteps and k modes

pyflation.analysis.adiabatic.Pr_spectrum_from_Pphimodes(phidot, Pphi_modes, axis)[source]

Return the spectrum of (first order) curvature perturbations P_R1 for each k.

For a multifield model this is given by:

\rm{Pr} = (\sum_K \dot{\phi_K}^2 )^{-2} 
    \sum_{I,J} \dot{\phi_I} \dot{\phi_J} P_{IJ}

where P_{IJ} = \sum_K \chi_{IK} \chi_{JK} and \chi are the mode matrix elements.

This is the unscaled version P_\mathcal{R} which is related to the scaled version by \mathcal{P}_\mathcal{R} = k^3/(2pi^2) P_\mathcal{R}.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

Pphi_modes : array_like

Mode matrix of first order perturbations spectra. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pr : array_like, dtype: float64

Array of Pr values for all timesteps and k modes

pyflation.analysis.adiabatic.Pzeta(m, tix=None, kix=None)[source]

Return the spectrum of (first order) curvature perturbations P_zeta for each k.

For a multifield model this is given by:

\rm{Pzeta} = \mathcal{P}_{\delta \rho} / (\rho^\dagger)^2

This is the unscaled version P_zeta which is related to the scaled version by \mathcal{P}_\zeta = k^3/(2\pi^2) P_\zeta.

Parameters :

m : Cosmomodels instance

model containing yresult with which to calculate spectrum

tix : integer, optional

index of timestep at which to calculate, defaults to full range of steps.

kix : integer, optional

integer of k value at which to calculate, defaults to full range of ks.

Returns :

Pzeta : array_like, dtype: float64

Array of Pzeta values for all timesteps and k modes

pyflation.analysis.adiabatic.Pzeta_spectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Return the spectrum of (first order) curvature perturbations P_zeta for each k.

For a multifield model this is given by:

\rm{Pzeta} = \mathcal{P}_{\delta \rho} / (\rho^\dagger)^2

This is the unscaled version P_zeta which is related to the scaled version by \mathcal{P}_\zeta = k^3/(2\pi^2) P_\zeta.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pzeta_spectrum : array_like, dtype: float64

Array of Pzeta values

pyflation.analysis.adiabatic.findHorizoncrossings(m, factor=1)[source]

Find horizon crossing for all ks

Parameters :

m : Cosmomodels instance

Model for which to calculate horizon crossings

factor : float

Coefficient for which to calculate horizon crossings e.g. k=a*H*factor

Returns :

hcross : array

horizon crossing values for all ks.

pyflation.analysis.adiabatic.findns(sPr, k, kix, running=False)[source]

Return the value of n_s

Parameters :

sPr : array_like

Power spectrum of scalar curvature perturbations at a specific time This should be the scaled power spectrum i.e.

\rm{sPr} = k^3/(2*\pi)^2 P_\mathcal{R},

<\mathcal{R}(k)\mathcal{R}(k')> = (2\pi^3) \delta(k+k') P_\mathcal{R}

The array should be one-dimensional indexed by the k value.

k : array_like

Array of k values for which sPr has been calculated.

kix : integer

Index value of k for which to return n_s.

running : boolean, optional

Whether running should be allowed or not. If true, a quadratic polynomial fit is made instead of linear and the value of the running alpha_s is returned along with n_s. Defaults to False.

Returns :

n_s : float

The value of the spectral index at the requested k value and timestep

\rm{n_s} = 1 - d \ln(\rm{sPr}) / d \ln(k)

evaluated at k[kix]

This is calculated using a polynomial least squares fit with numpy.polyfit. If running is True then a quadratic polynomial is fit, otherwise only a linear fit is made.

alpha_s : float, present only if running = True

If running=True the alpha_s value at k[kix] is returned in a tuple along with n_s.

pyflation.analysis.adiabatic.scaled_Pr(m, tix=None, kix=None)[source]

Return the spectrum of (first order) curvature perturbations \mathcal{P}_\mathcal{R} for each timestep and k mode.

This is the scaled power spectrum which is related to the unscaled version by \mathcal{P}_\mathcal{R} = k^3/(2\pi^2) P_\mathcal{R}.

Parameters :

m : Cosmomodels instance

Model class instance from which the yresult variable will be used to calculate P_R.

tix : integer, optional

index of timestep at which to calculate, defaults to full range of steps.

kix : integer, optional

integer of k value at which to calculate, defaults to full range of ks.

Returns :

calPr : array_like

Array of Pr values for all timesteps and k modes

pyflation.analysis.adiabatic.scaled_Pzeta(m, tix=None, kix=None)[source]

Return the spectrum of scaled (first order) curvature perturbations \mathcal{P}_\zeta for each timestep and k mode.

This is the scaled power spectrum which is related to the unscaled version by \mathcal{P}_\zeta = k^3/(2\pi^2) P_\zeta.

Parameters :

m : Cosmomodels instance

model containing yresult with which to calculate spectrum

tix : integer, optional

index of timestep at which to calculate, defaults to full range of steps.

kix : integer, optional

integer of k value at which to calculate, defaults to full range of ks.

Returns :

scaled_Pzeta : array_like

Array of Pzeta values for all timesteps and k modes

nonadiabatic Module

pyflation.analysis.nonadiabatic - Module to calculate relative pressure perturbations.

pyflation.analysis.nonadiabatic.Pdots(Vphi, phidot, H)[source]

Derivative of pressure of the background fields

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

Returns :

Pdotalpha : array

Derivative of the pressure for the background fields

Notes

All the arguments should have the same number of dimensions. Vphi and phidot should be arrays of the same size, but H should have a dimension of size 1 corresponding to the “field” dimension of the other variables.

pyflation.analysis.nonadiabatic.S_alternate(phidot, Pphi_modes, axis)[source]

Return the alternate spectrum of (first order) isocurvature perturbations P_S for each k. This is only available for a two field model and is given by:

\rm{P_S} = (H/\dot{\sigma})^2 \langle\delta s \delta s*\rangle

where

\dot{\sigma} = \sqrt{\dot{\phi}^2 + \dot{\chi}^2}

\delta s = - \dot{\chi}/\dot{\sigma} \delta \phi
             + \dot{\phi}/\dot{\sigma} \delta \chi

This is the unscaled version P_S which is related to the scaled version by \mathcal{P}_S = k^3/(2\pi^2) P_S.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

Pphi_modes : array_like

Mode matrix of first order perturbation power spectrum given by adiabatic.Pphi_matrix. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pr : array_like, dtype: float64

Array of Pr values for all timesteps and k modes

pyflation.analysis.nonadiabatic.Smodes(Vphi, phidot, H, modes, modesdot, axis)[source]

Isocurvature perturbation S of the fields given as quantum mode functions.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

result : array

Isocurvature perturbation S of the fields

pyflation.analysis.nonadiabatic.Sspectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Power spectrum of the isocurvature perturbation S.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Sspectrum : array

Spectrum of the isocurvature perturbation S.

pyflation.analysis.nonadiabatic.deltaPmatrix(Vphi, phidot, H, modes, modesdot, axis)[source]

Matrix of the first order perturbed pressure of the field components.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

result : array_like

The matrix of the first order perturbed pressure.

pyflation.analysis.nonadiabatic.deltaPnadmodes(Vphi, phidot, H, modes, modesdot, axis)[source]

Perturbed non-adiabatic pressure of the fields given as quantum mode functions.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

pyflation.analysis.nonadiabatic.deltaPnadspectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Power spectrum of the full perturbed non-adiabatic pressure.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

deltaPnadspectrum : array

Spectrum of the non-adiabatic pressure perturbation

pyflation.analysis.nonadiabatic.deltaPrelmodes(Vphi, phidot, H, modes, modesdot, axis)[source]

Perturbed relative pressure of the fields given as quantum mode functions.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

result : array

Perturbed relative pressure of the fields.

pyflation.analysis.nonadiabatic.deltaPrelspectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Power spectrum of the full perturbed relative pressure.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

deltaPrelspectrum : array

Spectrum of the perturbed relative pressure

pyflation.analysis.nonadiabatic.deltaPspectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Power spectrum of the full perturbed relative pressure.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

deltaPspectrum : array

Spectrum of the perturbed pressure

pyflation.analysis.nonadiabatic.deltarhosmatrix(Vphi, phidot, H, modes, modesdot, axis)[source]

Matrix of the first order perturbed energy densities of the field components.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

result : array_like

The matrix of the first order perturbed energy densities.

pyflation.analysis.nonadiabatic.deltarhospectrum(Vphi, phidot, H, modes, modesdot, axis)[source]

Power spectrum of the perturbed energy density.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

deltarhospectrum : array

Spectrum of the perturbed energy density

pyflation.analysis.nonadiabatic.dprel_from_model(m, tix=None, kix=None)[source]

Get the spectrum of delta Prel from a model instance.

Parameters :

m : Cosmomodels model instance

The model instance with which to perform the calculation

tix : integer

Index for timestep at which to perform calculation. Default is to calculate over all timesteps.

kix : integer

Index for k mode for which to perform the calculation. Default is to calculate over all k modes.

Returns :

spectrum : array

Array of values of the power spectrum of the relativistic pressure perturbation

pyflation.analysis.nonadiabatic.fullPdot(Vphi, phidot, H, axis=-1)[source]

Combined derivative in e-fold time of the pressure of the fields.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

axis : integer, optional

Specifies which axis is the field dimension, default is the last one.

Returns :

fullPdot : array

The deriative of the pressure of the fields summed over the fields.

pyflation.analysis.nonadiabatic.fullrhodot(phidot, H, axis=-1)[source]

Combined derivative in e-fold time of the energy density of the field.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

axis : integer, optional

Specifies which axis is the field dimension, default is the last one.

Returns :

fullrhodot : array

The derivative of the energy density summed over the fields.

pyflation.analysis.nonadiabatic.rhodots(phidot, H)[source]

Derivative in e-fold time of the energy densities of the individual fields.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

Both arrays should have the same number of dimensions, but H should have a :

dimension of size 1 corresponding to the field dimension of phidot. :

Returns :

rhodots : array

The derivative of the energy densities of the fields

pyflation.analysis.nonadiabatic.scaled_S_alternate_spectrum(phidot, Pphi_modes, axis, k)[source]

Return the scaled alternate spectrum of (first order) isocurvature perturbations for each k. This is only available for a two field model and is given by:

P_{\bar{S}} = (H/\dot{\sigma})^2 <\delta s \delta s*>

\dot{\sigma} = \sqrt{\dot{\phi}^2 + \dot{\chi}^2} 

\delta s = - \dot{\chi}/\dot{\sigma} \delta \phi
             + \dot{\phi}/\dot{\sigma} \delta \chi

This is the scaled version \mathcal{P}_{\bar{S}} which is related to the unscaled version by \mathcal{P}_{\bar{S}} = k^3/(2\pi^2) P_{\bar{S}}.

Parameters :

phidot : array_like

First derivative of the field values with respect to efold number N.

Pphi_modes : array_like

Mode matrix of first order perturbation power spectrum given by adiabatic.Pphi_matrix. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

Pr : array_like, dtype: float64

Array of Pr values for all timesteps and k modes

pyflation.analysis.nonadiabatic.scaled_S_from_model(m, tix=None, kix=None)[source]

Return the spectrum of isocurvature perturbations \mathcal{P}_\mathcal{S} for each timestep and k mode.

This is the scaled power spectrum which is related to the unscaled version by \mathcal{P}_\mathcal{S} = k^3/(2\pi^2) P_\mathcal{S}.

Parameters :

m : Cosmomodels instance

Model class instance from which the yresult variable will be used to calculate P_R.

tix : integer, optional

index of timestep at which to calculate, defaults to full range of steps.

kix : integer, optional

integer of k value at which to calculate, defaults to full range of ks.

Returns :

scaled_S : array_like

Array of spectrum values for all timesteps and k modes

pyflation.analysis.nonadiabatic.scaled_S_spectrum(Vphi, phidot, H, modes, modesdot, axis, k)[source]

Power spectrum of S scaled with k^3/(2*pi^2)

Assumes that k dimension is last.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

k : array

The values of k to scale the result with.

Returns :

scaled_S_spectrum : array

Scaled spectrum of the isocurvature perturbation S

pyflation.analysis.nonadiabatic.scaled_dP_spectrum(Vphi, phidot, H, modes, modesdot, axis, k)[source]

Power spectrum of delta P scaled with k^3/(2*pi^2)

Assumes that k dimension is last.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

k : array

The values of k to scale the result with.

Returns :

scaled_dP_spectrum : array

Scaled spectrum of the perturbed pressure

pyflation.analysis.nonadiabatic.scaled_dPnad_spectrum(Vphi, phidot, H, modes, modesdot, axis, k)[source]

Power spectrum of delta Pnad scaled with k^3/(2*pi^2)

Assumes that k dimension is last.

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

modes : array_like

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot : array_like

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis : integer

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

k : array

The values of k to scale the result with.

Returns :

scaled_dPnad_spectrum : array

Scaled spectrum of the non-adiabatic pressure perturation.

pyflation.analysis.nonadiabatic.slope_of_S_spectrum(scaled_S, k, kix=None, running=False)[source]

Return the value of the slope of the k-scaled spectrum of S

Parameters :

scaled_S : array_like

Power spectrum of isocurvature perturbations at a specific time This should be the scaled power spectrum i.e.

\rm{scaled_S} = k^3/(2\pi)^2 P_S,

<S(k)S(k')> = (2\pi^3) \delta(k+k') P_S

The array should be one-dimensional indexed by the k value.

k : array_like

Array of k values for which sPr has been calculated.

kix : integer

Index value of k for which to return spec_S.

running : boolean, optional

Whether running should be allowed or not. If true, a quadratic polynomial fit is made instead of linear and the value of the running is returned along with the slope. Defaults to False.

Returns :

spec_S : float

The value of the spectral index of S at the requested k value and timestep. Normalised in the same way as n_s

\rm{spec_S} = 1 - d \log(\rm{scaled_S}) / d \log(k)

evaluated at k[kix].

This is calculated using a polynomial least squares fit with numpy.polyfit. If running is True then a quadratic polynomial is fit, otherwise only a linear fit is made.

running_S : float, present only if running = True

If running=True the value of the derivative of the slope at k[kix] is returned in a tuple along with spec_S.

pyflation.analysis.nonadiabatic.soundspeeds(Vphi, phidot, H)[source]

Sound speeds of the background fields

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

Returns :

calphasq : array

sound speed (squared) of the background fields

Notes

All the arguments should have the same number of dimensions. Vphi and phidot should be arrays of the same size, but H should have a dimension of size 1 corresponding to the “field” dimension of the other variables.

pyflation.analysis.nonadiabatic.totalsoundspeed(Vphi, phidot, H, axis)[source]

Total sound speed of the fluids

Parameters :

Vphi : array_like

First derivative of the potential with respect to the fields

phidot : array_like

First derivative of the field values with respect to efold number N.

H : array_like

The Hubble parameter

axis : integer

Index of dimension to sum over (field dimension).

All the arguments should have the same number of dimensions. Vphi and phidot :

should be arrays of the same size, but H should have a dimension of size 1 :

corresponding to the “field” dimension of the other variables. :

Returns :

csq : array_like

The total sound speed of the fluid, csq = P’/rho’

utilities Module

pyflation.analysis.utilities - Helper module for analysis package

pyflation.analysis.utilities.components_from_model(m, tix=None, kix=None)[source]

Get the component variables of delta Prel from a model instance.

Parameters :

m: Cosmomodels model instance :

The model instance with which to perform the calculation

tix: integer :

Index for timestep at which to perform calculation. Default is to calculate over all timesteps.

kix: integer :

Index for k mode for which to perform the calculation. Default is to calculate over all k modes.

Returns :

Tuple containing: :

Vphi: array_like :

First derivative of the potential with respect to the fields

phidot: array_like :

First derivative of the field values with respect to efold number N.

H: array_like :

The Hubble parameter

modes: array_like :

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot: array_like :

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis: integer :

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

pyflation.analysis.utilities.correct_shapes(Vphi, phidot, H, modes, modesdot, axis)[source]

Return variables with the correct shapes for calculation.

Parameters :

Vphi: array_like :

First derivative of the potential with respect to the fields

phidot: array_like :

First derivative of the field values with respect to efold number N.

H: array_like :

The Hubble parameter

modes: array_like :

Mode matrix of first order perturbations. Component array should have two dimensions of length nfields.

modesdot: array_like :

Mode matrix of N-derivative of first order perturbations. Component array should have two dimensions of length nfields.

axis: integer :

Specifies which axis is first in mode matrix, e.g. if modes has shape (100,3,3,10) with nfields=3, then axis=1. The two mode matrix axes are assumed to be beside each other so (100,3,10,3) would not be valid.

Returns :

result: tuple :

Tuple of the variables with correct shapes in the form (Vphi, phidot, H, modes, modesdot, axis)

pyflation.analysis.utilities.flattenmodematrix(modematrix, nfields, ix1=None, ix2=None)[source]

Flatten the mode matrix given into nfield^2 long vector.

Parameters :

modematrix: array :

Array of values with two nfields long dimension

nfields: integer :

Number of fields

ix1: integer, optional :

Index of first nfields long dimension

ix2: integer, optional :

Index of second nfields long dimension

Returns :

mreshaped: array :

Reshaped mode matrix array

pyflation.analysis.utilities.getmodematrix(y, nfields, ix=None, ixslice=None)[source]

Helper function to reshape flat nfield^2 long y variable into nfield*nfield mode matrix. Returns a view of the y array (changes will be reflected in underlying array).

Parameters :

y: array :

Array of y values in which is nfields^2 long in dimension specified by ix

nfields: integer :

Number of fields

ix: integer :

Index of dimension which is nfields^2 long

ixslice: index slice, optional :

The index slice of y to use, defaults to full extent of y.

Returns :

result: view of y array with shape nfield*nfield structure :

pyflation.analysis.utilities.kscaling(k, kix=None)[source]

Return the k scaling for power spectrum calculations.

Parameters :

k: array :

Array of k values to use

kix: integer, optional :

Index of k value in array to use, defaults to None, which uses full array.

Returns :

kscaling: array :

k scaling values, kscaling = k^3/(2*pi^2)

Note for a single variable (not a spectrum) you may want the square root :

of this result. :

pyflation.analysis.utilities.makespectrum(modes_I, axis)[source]

Calculate spectrum of input.

Spectrum = Sum_I modes_I * modes_I.conjugate, where sum over I is on the dimension denoted by axis

Parameters :

modes_I: array_like :

Array of (complex) values to calculate spectrum with.

axis: integer :

Dimension along which to sum the spectrum over.

Returns :

spectrum: array_like, float :

Real valued array of calculated spectrum.

pyflation.analysis.utilities.spectral_index(y, k, kix, running=False)[source]

Return the value of spectral index (and running) of y at k[kix]

Parameters :

y: array_like :

Array of values The array should be one-dimensional indexed by the x values.

k: array_like :

Array of k values for which y has been calculated.

kix: integer :

Index value of k for which to return the spectral index.

running: boolean, optional :

Whether running should be allowed or not. If true, a quadratic polynomial fit is made instead of linear and the value of the running is returned along with the spectral index. Defaults to False.

Returns :

spec_index: float :

The value of the spectral index at the requested x value. This is calculated as with the spectral index, so the value is 1 for a flat line.

spec_index = 1 - d ln(y) / d ln(k) evaluated at k[kix]

This is calculated using a polynomial least squares fit with numpy.polyfit. If running is True then a quadratic polynomial is fit, otherwise only a linear fit is made.

alpha: float, present only if running = True :

If running=True the alpha_s value at k[kix] is returned in a tuple along with spec_index.

Table Of Contents

Previous topic

sourceterm Package

Next topic

sourceterm Package

This Page