pyflation Package

pyflation Package

Pyflation - Cosmological simulations in Python

Pyflation is a python package to simulate cosmological perturbations in the early universe. Using the Klein-Gordon equations for both first and second order perturbations, the evolution and behaviour of these perturbations can be studied.

The main entry point to this package is the cosmomodels module, which contains the main simulation classes. Configuration of the simulation runs can be changed in the configuration.py and run_config.py files.

For more information please visit http://pyflation.ianhuston.net.

cmpotentials Module

cmpotentials.py - Cosmological potentials for cosmomodels.py

Provides functions which can be used with cosmomodels.py. Default parameter values are included but can also be specified as a dictionary.

pyflation.cmpotentials.bump_nothirdderiv(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 ( 1 + c*sech((phi-phi_b) / d) where m is the mass of the inflaton field and c, d and phi_b are provided. Form is taken from Chen etal. arxiv:0801.3295.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.bump_potential(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 ( 1 + c*sech((phi-phi_b) / d) where m is the mass of the inflaton field and c, d and phi_b are provided. Form is taken from Chen etal. arxiv:0801.3295.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.hilltopaxion(y, params=None)[source]

Return the potential and its first three derivatives for a hilltop axion model.

The potential is given by

V = 0.5 m^2 \varphi^2 + \Lambda^4 (1 - \cos(2\pi\chi/f))

where the parameters are Lambda, m, and f . Needs nfields=2.

Parameters :

y: array :

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params: dict :

Dictionary of parameter values labelled “Lambda” , “m”, “f”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3: tuple of arrays :

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.hybrid2and4(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for hybrid potential V = -m^2/2 phi^2 +lambda/4 phi^4

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameters “mass” and “lambda” which specifies the variables.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

lambda can be specified in the dictionary params or otherwise it defaults to the value as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

mass can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.hybridquadratic(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V = 1/2 m1^2 phi^2 + 1/2 m2^2 chi^2 where m1 and m2 are the masses of the fields. Needs nfields=2.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameters “m1” and “m2” specified above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.hybridquartic(y, params=None)[source]

Return the potential and its first three derivatives for the hybrid quartic model.

The potential is given by

V = \Lambda^4 [ (1-\chi^2/v^2)^2 + \phi^2/\mu^2 
            + 2\phi^2\chi^2/(\phi_c^2 v^2) ]

where the parameter are \Lambda, v, \mu and \phi_c. Needs nfields=2.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values labelled “lambda” , “v”, “mu”, “phi_c”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.inflection(y, params=None)[source]

Return the potential and its first three derivatives for an inflection point model.

The potential is given by

V = V_0 + 0.5 m^2 \phi^2 + g \chi + 1/6 \lambda \chi^3 + \lambda/(8 r) \chi^4

where V_0 = 0.75gr + \lambda/24 r^3 + g/(4r^3) and the parameters are \lambda, m, g and r. Needs nfields=2.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values labelled “lambda” , “m”, “g”, “r”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.lambdaphi4(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/4 lambda phi^4 for a specified lambda.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “lambda” which specifies lambda above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

lambda can be specified in the dictionary params or otherwise it defaults to the value as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.linde(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for Linde potential V = -m^2/2 phi^2 +lambda/4 phi^4 + m^4/4lambda

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameters “mass” and “lambda” which specifies the variables.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

lambda can be specified in the dictionary params or otherwise it defaults to the value as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

mass can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.msqphisq(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 where m is the mass of the inflaton field.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above.

Returns :

——- :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.msqphisq_withV0(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 + V0 where m is the mass of the inflaton field.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.nflation(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V = sum_lpha 1/2 m^2 phi_lpha^2 where m is the mass of each of the fields.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above. The number of fields is specified through “nfields”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.phi2over3(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V= sigma phi^(2/3) for a specified sigma.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “sigma” which specifies lambda above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

sigma can be specified in the dictionary params or otherwise it defaults to the value as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.productexponential(y, params=None)[source]

Return the potential and its first three derivatives for a product exponential potential.

The potential is given by

V = V_0 \phi^2 \exp(-\lambda \chi^2)

where the parameters are V_0, \lambda. Needs nfields=2.

Parameters :

y: array :

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params: dict :

Dictionary of parameter values labelled “lambda” , “V_0”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3: tuple of arrays :

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.quartictwofield(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V= 1/2(m1^2 phi^2 + 1/2 l1 phi^4 + m2^2 chi^2 + 1/2 l2 chi^4) where m1, m2, l1, l2 are parameters. Needs nfields=2.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameters “m1”, “m2”, “l1”, “l2”, as specified above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.resonance(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 ( 1 + c*sin(phi / d) ) where m is the mass of the inflaton field and c, d and phi_b are provided. Form is taken from Chen etal. arxiv:0801.3295.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above, and the parameters “c” and “d” which tune the oscillation.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

pyflation.cmpotentials.ridge_twofield(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=V0 - g phi - 1/2 m^2 chi^2 where g is a parameter and m is the mass of the chi field. Needs nfields=2.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameters “V0”, “g”, “m”.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

pyflation.cmpotentials.step_potential(y, params=None)[source]

Return (V, dV/dphi, d2V/dphi2, d3V/dphi3) for V=1/2 m^2 phi^2 ( 1 + c*tanh((phi-phi_s) / d) where m is the mass of the inflaton field and c, d and phi_s are provided. Form is taken from Chen etal. arxiv:0801.3295.

Parameters :

y : array

Array of variables with background phi as y[0] If you want to specify a vector of phi values, make sure that the first index still runs over the different variables, using newaxis if necessary.

params : dict

Dictionary of parameter values in this case should hold the parameter “mass” which specifies m above.

Returns :

U, dUdphi, d2Udphi2, d3Udphi3 : tuple of arrays

Tuple of the potential and its first three derivatives.

Notes

m can be specified in the dictionary params or otherwise it defaults to the mass as normalized with the WMAP spectrum Pr = 2.457e-9 at the WMAP pivot scale of 0.002 Mpc^-1.

configuration Module

Configuration file

The main configuration options are for logging. By changing _debug to 1 (default is 0) much more debugging information will be added to the log files. The overall logging level can also be set using the LOGLEVEL variable. This level can be overridden using command line options to the scripts.

cosmographs Module

cosmographs.py - Graphing functions for pyflation package

This module provides helper functions for graphing the results of pyflation simulations using the Matplotlib package (http://matplotlib.sf.net).

Especially useful is the multi_format_save function which saves the specified figure to different formats as requested.

exception pyflation.cosmographs.CosmoGraphError[source]

Bases: exceptions.StandardError

Generic error for graphing facilities.

class pyflation.cosmographs.LogFormatterTeXExponent(*args, **kwargs)[source]

Bases: matplotlib.ticker.LogFormatter, object

Extends pyplot.LogFormatter to use tex notation for tick labels.

pyflation.cosmographs.calN_figure(ts, ys, fig=None, plot_fn=None, models_legends=None, ylabel=None, size='large', ls=['r-', 'g--', 'b:'], halfticks=False)[source]

Create a figure using mathcal{N} on the x-axis.

Parameters :

ts: list of arrays :

x-axis values

ys: list of arrays :

y-axis values, each element should be the same length as corresponding ts element.

fig: Pylab figure instance, optional :

Figure to draw on, default is to create a new figure.

plot_fn: function, optional :

The Pylab plotting function to use, defaults to P.plot.

models_legends: list, optional :

List of raw strings (including LaTeX) to use in legend of plot. Defaults to not plotting legend.

ylabel: string, optional :

Raw string (including LaTeX) for y-axis, defaults to empty.

size: string, optional :

One of “small”, “large”, “half”, which specifies the size of the figure by using pyflation.cosmographs.set_size. Defaults to “large”.

ls: list, optional :

List of linestyle strings to use with each successive line. Defaults to [“r-”, “g–”, “b:”]

halfticks: boolean, optional :

If True only show half the ticks on the y-axis

Returns :

——- :

fig: the figure instance :

pyflation.cosmographs.generic_figure(xs, ys, fig=None, plot_fn=None, models_legends=None, xlabel=None, ylabel=None, size='large', ls=['r-', 'g--', 'b:'], halfticks=False)[source]

Create a generic figure with standard x-axis.

Parameters :

xs: list of arrays :

x-axis values

ys: list of arrays :

y-axis values, each element should be the same length as corresponding xs element.

fig: Pylab figure instance, optional :

Figure to draw on, default creates a new figure.

plot_fn: function, optional :

The Pylab plotting function to use, defaults to P.plot.

models_legends: list, optional :

List of raw strings (including LaTeX) to use in legend of plot. Defaults to not plotting legend.

xlabel: string, optional :

Raw string (including LaTeX) for x-axis, defaults to empty.

ylabel: string, optional :

Raw string (including LaTeX) for y-axis, defaults to empty.

size: string, optional :

One of “small”, “large”, “half”, which specifies the size of the figure by using pyflation.cosmographs.set_size. Defaults to “large”.

ls: list, optional :

List of linestyle strings to use with each successive line. Defaults to [“r-”, “g–”, “b:”]

halfticks: boolean, optional :

If True only show half the ticks on the y-axis

Returns :

fig: the figure instance :

pyflation.cosmographs.makeklegend(fig, k)[source]

Attach a legend to the figure specified outlining the k modes used.

pyflation.cosmographs.multi_format_save(filenamestub, fig=None, formats=None, overwrite=False, **kwargs)[source]

Save figure in multiple formats at once.

Parameters :

filenamestub: String :

Filename with path to which will be added the appropriate extension for each different file saved.

fig: Matplotlib figure object, optional :

Figure to save to disk. Defaults to current figure.

formats: list-like, optional :

List of format specifiers to save to, default is [“pdf”, “eps”, “png”, “svg”] Must be of types supported by current installation.

Other kwargs: These are provided to matplotlib.pylab.savefig. :

See there for documentation.

Returns :

filenames: list :

list of saved file names

pyflation.cosmographs.reversexaxis(a=None)[source]

Reverse the direction of the x axis for the axes object a.

pyflation.cosmographs.save(fname, dir=None, fig=None)[source]
pyflation.cosmographs.save_with_prompt(fname)[source]
pyflation.cosmographs.set_size(fig, size='large')[source]
pyflation.cosmographs.set_size_half(fig)[source]
pyflation.cosmographs.set_size_large(fig)[source]
pyflation.cosmographs.set_size_small(fig)[source]

cosmomodels Module

cosmomodels.py - Cosmological Model simulations

Provides classes for modelling cosmological inflationary scenarios. Especially important classes are:

The make_wrapper_model() function takes a filename and returns a model instance corresponding to the one stored in the file. This allows easier access to the results than through a direct inspection of the HDF5 results file.

class pyflation.cosmomodels.BasicBgModel(ystart=array([ 0.1, 0.1, 0.1]), tstart=0.0, tend=120.0, tstep_wanted=0.02, solver='rkdriver_tsix')[source]

Bases: pyflation.cosmomodels.CosmologicalModel

Basic model with background equations
Array of dependent variables y is given by:

y[0] - phi_0 : Background inflaton y[1] - dphi_0/deta : First deriv of phi y[2] - a : Scale Factor

derivs(y, t, **kwargs)[source]

Basic background equations of motion. dydx[0] = dy[0]/deta etc

potentials(y, pot_params=None)[source]

Return value of potential at y, along with first and second derivs.

class pyflation.cosmomodels.CanonicalBackground(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.PhiModels

Basic model with background equations in terms of n
Array of dependent variables y is given by:

y[0] - phi_0 : Background inflaton y[1] - dphi_0/dn : First deriv of phi y[2] - H: Hubble parameter

derivs(y, t, **kwargs)[source]

Basic background equations of motion. dydx[0] = dy[0]/dn etc

class pyflation.cosmomodels.CanonicalFirstOrder(k=None, ainit=None, *args, **kwargs)[source]

Bases: pyflation.cosmomodels.PhiModels

First order model using efold as time variable with multiple fields.

nfields holds the number of fields and the yresult variable is then laid out as follows:

yresult[0:nfields*2] : background fields and derivatives yresult[nfields*2] : Hubble variable H yresult[nfields*2 + 1:] : perturbation fields and derivatives

derivs(y, t, **kwargs)[source]

Return derivatives of fields in y at time t.

setfieldindices()[source]

Set field indices. These can be used to select only certain parts of the y variable, e.g. y[self.bg_ix] is the array of background values.

class pyflation.cosmomodels.CanonicalHomogeneousSecondOrder(k=None, ainit=None, *args, **kwargs)[source]

Bases: pyflation.cosmomodels.PhiModels

Second order homogeneous model using efold as time variable. y[0] - delta arphi_2 : Second order perturbation [Real Part] y[1] - delta arphi_2^prime : Derivative of second order perturbation [Real Part] y[2] - delta arphi_2 : Second order perturbation [Imag Part] y[3] - delta arphi_2^prime : Derivative of second order perturbation [Imag Part]

derivs(y, t, **kwargs)[source]

Equation of motion for second order perturbations including source term

class pyflation.cosmomodels.CanonicalRampedSecondOrder(k=None, ainit=None, *args, **kwargs)[source]

Bases: pyflation.cosmomodels.PhiModels

Second order model using efold as time variable. y[0] - delta arphi_2 : Second order perturbation [Real Part] y[1] - delta arphi_2^prime : Derivative of second order perturbation [Real Part] y[2] - delta arphi_2 : Second order perturbation [Imag Part] y[3] - delta arphi_2^prime : Derivative of second order perturbation [Imag Part]

derivs(y, t, **kwargs)[source]

Equation of motion for second order perturbations including source term

class pyflation.cosmomodels.CanonicalSecondOrder(k=None, ainit=None, *args, **kwargs)[source]

Bases: pyflation.cosmomodels.PhiModels

Second order model using efold as time variable. y[0] - delta arphi_2 : Second order perturbation [Real Part] y[1] - delta arphi_2^prime : Derivative of second order perturbation [Real Part] y[2] - delta arphi_2 : Second order perturbation [Imag Part] y[3] - delta arphi_2^prime : Derivative of second order perturbation [Imag Part]

derivs(y, t, **kwargs)[source]

Equation of motion for second order perturbations including source term

class pyflation.cosmomodels.CombinedCanonicalFromFile(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.MultiStageDriver

Model class for combined first and second order data, assumed to be used with a file wrapper.

deltaphi[source]

Return the calculated values of $deltaphi$ for all times and modes.

The result is stored as the instance variable self.deltaphi but will be recomputed if recompute is True.

Parameters :

recompute : boolean, optional

Should the values be recomputed? Default is False.

Returns :

deltaphi : array_like

Array of $deltaphi$ values for all timesteps and k modes.

dp1[source]

Return (and save) the first order perturbation.

dp2[source]

Return (and save) the first order perturbation.

class pyflation.cosmomodels.CosmologicalModel(ystart=None, simtstart=0.0, tstart=0.0, tstartindex=None, tend=83.0, tstep_wanted=0.01, solver='rkdriver_tsix', potential_func=None, pot_params=None, nfields=1, **kwargs)[source]

Bases: object

Generic class for cosmological model simulations. Contains run() method which chooses a solver and runs the simulation.

All cosmological model classes are subclassed from this one.

callingparams()[source]

Returns list of parameters to save with results.

createhdf5structure(filename, grpname='results', yresultshape=None, hdf5complevel=2, hdf5complib='blosc')[source]

Create a new hdf5 file with the structure capable of holding results.

Parameters :

filename : string

Path including filename of file to create

grpname : string, optional

Name of the HDF5 group to create, default is “results”

yresultshape : tuple, optional

Shape of yresult variable to store

hdf5complevel : integer, optional

Compression level to use with PyTables, default 2.

hdf5complib : string, optional

Compression library to use with PyTables, default “blosc”.

Returns :

rf : file handle

Handle of file created

derivs(yarray, t)[source]

Return an array of derivatives of the dependent variables yarray at timestep t

findH(potential, y)[source]

Return value of comoving Hubble variable given potential and y.

gethf5paramsdict()[source]

Describes the fields required to save the calling parameters.

potentials(y, pot_params=None)[source]

Return a 4-tuple of potential, 1st, 2nd and 3rd derivs given y.

run(saveresults=True)[source]

Execute a simulation run using the parameters already provided.

saveallresults(filename=None, filetype='hf5', yresultshape=None, **kwargs)[source]

Saves results already calculated into a file.

saveresultsinhdf5(rf, grpname='results')[source]

Save simulation results in a HDF5 format file with filename.

Parameters :

rf : filelike

File to save results in

grpname : string, optional

Name of the HDF5 group to create in the file

solverlist = ['rkdriver_tsix']
class pyflation.cosmomodels.FOCanonicalTwoStage(bgystart=None, tstart=0.0, tstartindex=None, tend=83.0, tstep_wanted=0.01, k=None, ainit=None, solver='rkdriver_tsix', bgclass=None, foclass=None, potential_func=None, pot_params=None, simtstart=0, nfields=1, **kwargs)[source]

Bases: pyflation.cosmomodels.MultiStageDriver

Uses a background and firstorder class to run a full (first-order) simulation. Main additional functionality is in determining initial conditions. Variables finally stored are as in first order class.

H[source]

Hubble parameter

a[source]

Scale factor of the universe

deltaphi(recompute=False)[source]

Return the calculated values of $deltaphi$ for all times, fields and modes. For multifield systems this is the quantum matrix of solutions:

hat{deltaphi} = Sum_{lpha, I} xi_{lpha I} hat{a}_I

The result is stored as the instance variable m.deltaphi but will be recomputed if recompute is True.

Parameters :

recompute : boolean, optional

Should the values be recomputed? Default is False.

Returns :

deltaphi : array_like, dtype: complex128

Array of $deltaphi$ values for all timesteps, fields and k modes.

dpdotmodes[source]

Quantum modes of derivatives of first order perturbations

dpmodes[source]

Quantum modes of first order perturbations

getdeltaphi()[source]
getfoystart(ts=None, tsix=None)[source]

Model dependent setting of ystart

phidots[source]

Derivatives of background fields w.r.t N phi_i^dagger

phis[source]

Background fields phi_i

run(saveresults=True)[source]

Run the full model.

The background model is first run to establish the end time of inflation and the start times for the k modes. Then the initial conditions are set for the first order variables. Finally the first order model is run and the results are saved if required.

Parameters :

saveresults : boolean, optional

Should results be saved at the end of the run. Default is False.

Returns :

filename : string

name of the results file if any

runbg()[source]

Run bg model after setting initial conditions.

runfo()[source]

Run first order model after setting initial conditions.

setfieldindices()[source]

Set field indices. These can be used to select only certain parts of the y variable, e.g. y[self.bg_ix] is the array of background values.

setfoics()[source]

After a bg run has completed, set the initial conditions for the first order run.

class pyflation.cosmomodels.FONoPhase(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.FOCanonicalTwoStage

First order two stage class which does not include a phase in the initial conditions for the first order field.

getfoystart(ts=None, tsix=None)[source]

Model dependent setting of ystart

class pyflation.cosmomodels.FOSuppressOneField(suppress_ix=0, *args, **kwargs)[source]

Bases: pyflation.cosmomodels.FOCanonicalTwoStage

First order two stage class which does not include a phase in the initial conditions for the first order field.

getfoystart(ts=None, tsix=None)[source]

Model dependent setting of ystart

class pyflation.cosmomodels.FixedainitTwoStage(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.FOCanonicalTwoStage

First order driver class with ainit fixed to a specified value. This is useful for comparing models which have different H behaviour and therefore different ainit values.

It should be remembered that these models with fixed ainit are equivalent to changing the number of efoldings from the end of inflation until today.

setfoics()[source]

After a bg run has completed, set the initial conditions for the first order run.

exception pyflation.cosmomodels.ModelError[source]

Bases: exceptions.StandardError

Generic error for model simulating. Attributes include current results stack.

class pyflation.cosmomodels.MultiStageDriver(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.CosmologicalModel

Parent of all multi (2 or 3) stage models. Contains methods to determine ns, k crossing and outlines methods to find Pr that are implemented in children.

Pr[source]

The power spectrum of comoving curvature perturbation. This is the unscaled spectrum P_R calculated for all timesteps and ks. Calculated using the pyflation.analysis package.

Pzeta[source]

The power spectrum of the curvature perturbation on uniform energy density hypersurfaces.

Calculated using the pyflation.analysis package.

callingparams()[source]

Returns list of parameters to save with results.

deltaphi[source]

Return the value of deltaphi for this model, recomputing if necessary.

findHorizoncrossings(factor=1)[source]

Find horizon crossing time indices and efolds for all ks

Parameters :

factor : float

Value of coefficient to calculate crossing time, k=a*H*factor

Returns :

hcrossings : array

Array of (kcrindex, kcrefold) pairs of time index and efold number pairs

find_efolds_after_inflation(Hend, Hreh=None)[source]

Calculate the number of efolds after inflation given the reheating temperature and assuming standard calculation of radiation and matter phases.

Parameters :

Hend : scalar, value of Hubble parameter at end of inflation

Hreh : scalar (default=Hend), value of Hubble parameter at end of reheating

Returns :

N : scalar, number of efolds after the end of inflation until today.

N = ln (a_today/a_end) where a_end is scale factor at end of inflation.

References

See Huston, arXiv: 1006.5321, Liddle and Lyth, Cambridge University Press 2000, or Peiris and Easther, JCAP 0807 (2008) 024, arXiv:0805.2154, for more details on calculation of post-inflation expansion.

finda_0(Hend, Hreh=None, a_end=None)[source]

Given the Hubble parameter at the end of inflation and at the end of reheating, and the scale factor at the end of inflation, calculate the scale factor today.

Parameters :

Hend : scalar

value of Hubble parameter at end of inflation

Hreh : scalar, optional

value of Hubble parameter at end of reheating, default is Hend.

a_end : scalar, optional

value of scale factor at the end of inflation, default is calculated from results.

Returns :

a_0 : scalar

scale factor today

finda_end(Hend, Hreh=None, a_0=1)[source]

Given the Hubble parameter at the end of inflation and at the end of reheating calculate the scale factor at the end of inflation.

This function assumes by default that the scale factor = 1 today and should be used with caution. A more correct approach is to call find_efolds_after_inflation directly and to use the result as required.

Parameters :

Hend : scalar

value of Hubble parameter at end of inflation.

Hreh : scalar, optional

value of Hubble parameter at end of reheating, default is Hend.

a_0 : scalar, optional

value of scale factor today, default is 1.

Returns :

a_end : scalar

scale factor at the end of inflation.

findallkcrossings(t, H)[source]

Iterate over findkcrossing to get full list

Parameters :

t : array

Array of t values to calculate over

H : array

Array of Hubble parameter values, should be the same shape as t

Returns :

kcrossings : array

Array of (kcrindex, kcrefold) pairs of index (in to t) and efold number at which each k in self.k crosses the horizon (k=a*H).

findkcrossing(k, t, H, factor=None)[source]

Given k, time variable and Hubble parameter, find when mode k crosses the horizon.

Parameters :

k : float

Single k value to compute crossing time with

t : array

Array of time values

H : array

Array of values of the Hubble parameter

factor : float, optional

coefficient of crossing k = a*H*factor

Returns :

kcrindex, kcrefold : tuple

Tuple containing k cross index (in t variable) and the efold number e.g. t[kcrindex]

getfoystart()[source]

Return model dependent setting of ystart

gethf5paramsdict()[source]

Describes the fields required to save the calling parameters.

scaled_Pr[source]

The power spectrum of comoving curvature perturbation.

Calculated using the pyflation.analysis package.

scaled_Pzeta[source]

The power spectrum of comoving curvature perturbation.

Calculated using the pyflation.analysis package.

class pyflation.cosmomodels.PhiModels(*args, **kwargs)[source]

Bases: pyflation.cosmomodels.CosmologicalModel

Parent class for models implementing the scheme in Malik 06[astro-ph/0610864]

findH(U, y)[source]

Return value of Hubble variable, H at y for given potential.

findinflend()[source]

Find the efold time where inflation ends, i.e. the hubble flow parameter epsilon >1. Returns tuple of endefold and endindex (in tresult).

getepsilon()[source]

Return an array of epsilon = -dot{H}/H values for each timestep.

potentials(y, pot_params=None)[source]

Return value of potential at y, along with first and second derivs.

class pyflation.cosmomodels.SOCanonicalThreeStage(second_stage, soclass=None, ystart=None, **soclassargs)[source]

Bases: pyflation.cosmomodels.MultiStageDriver

Runs third stage calculation (typically second order perturbations) using a two stage model instance which could be wrapped from a file.

H[source]

Hubble parameter

a[source]

Scale factor of the universe

deltaphi[source]

Return the calculated values of $deltaphi$ for all times and modes.

The result is stored as the instance variable self.deltaphi but will be recomputed if recompute is True.

Parameters :

recompute : boolean, optional

Should the values be recomputed? Default is False.

Returns :

deltaphi : array_like

Array of $deltaphi$ values for all timesteps and k modes.

dp2dotmodes[source]

Quantum modes of derivatives of second order perturbations

dp2modes[source]

Quantum modes of second order perturbations

dpdotmodes[source]

Quantum modes of derivatives of first order perturbations

dpmodes[source]

Quantum modes of first order perturbations

phidots[source]

Derivatives of background fields w.r.t N phi_i^dagger

phis[source]

Background fields phi_i

run(saveresults=True)[source]

Run simulation and save results.

runso()[source]

Run second order model.

setfieldindices()[source]

Set field indices. These can be used to select only certain parts of the y variable, e.g. y[self.bg_ix] is the array of background values.

setup_soclass()[source]

Initialize the second order class that will be used to run simulation.

class pyflation.cosmomodels.SOHorizonStart(second_stage, soclass=None, ystart=None, **soclassargs)[source]

Bases: pyflation.cosmomodels.SOCanonicalThreeStage

Runs third stage calculation (typically second order perturbations) using a two stage model instance which could be wrapped from a file.

Second order calculation starts at horizon crossing.

class pyflation.cosmomodels.TestModel(ystart=array([ 1., 1.]), tstart=0.0, tend=1.0, tstep_wanted=0.01)[source]

Bases: pyflation.cosmomodels.CosmologicalModel

Test class defining a very simple function

derivs(y, t, **kwargs)[source]

Very simple set of ODEs

pyflation.cosmomodels.make_wrapper_model(modelfile, *args, **kwargs)[source]

Return a wrapper class that provides the given model class from a file.

pyflation.cosmomodels.profile(f)[source]

helpers Module

helpers.py - Helper functions

Provides helper functions for use by the package elsewhere.

pyflation.helpers.cartesian(arrays, out=None)[source]

Generate a cartesian product of input arrays.

Parameters :

arrays : list of array-like

1-D arrays to form the cartesian product of.

out : ndarray

Array to place the cartesian product in.

Returns :

out : ndarray

2-D array of shape (M, len(arrays)) containing cartesian products formed of input arrays.

Examples

>>> cartesian(([1, 2, 3], [4, 5], [6, 7]))
array([[1, 4, 6],
       [1, 4, 7],
       [1, 5, 6],
       [1, 5, 7],
       [2, 4, 6],
       [2, 4, 7],
       [2, 5, 6],
       [2, 5, 7],
       [3, 4, 6],
       [3, 4, 7],
       [3, 5, 6],
       [3, 5, 7]])
pyflation.helpers.cartesian_product(lists, previous_elements=[])[source]

Generator of cartesian products of lists.

pyflation.helpers.ensurepath(path)[source]

Check that the path for given directory exists and create it if not.

pyflation.helpers.eto10(number)[source]

Convert scientific notation e.g. 1e-5 to 1x10^{-5} for use in LaTeX, converting to string.

pyflation.helpers.find_nearest(array, value)[source]

Find the index of the number in array which is nearest to value.

pyflation.helpers.find_nearest_ix(array, value)[source]

Find the index of the number in array which is nearest to value.

pyflation.helpers.getintfunc(x)[source]

Return the correct function to integrate with.

Checks the given set of values and returns either scipy.integrate.romb or scipy.integrate.simps. This depends on whether the number of values is a power of 2 + 1 as required by romb.

Parameters :

x: array_like :

Array of x values to check

Returns :

intfunc: function object :

Correct integration function depending on length of x.

fnargs: dictionary :

Dictionary of arguments to integration function.

pyflation.helpers.getkend(kinit, deltak, numsoks)[source]

Correct kend value given the values of kinit, deltak and numsoks.

Parameters :

kinit: float :

Initial k value in range.

deltak: float :

Difference between two k values in range.

numsoks: integer :

The number of second order k modes required.

Returns :

kend: float :

The value of the end of the first order k range required so that the second order k range has numsoks modes in it.

pyflation.helpers.invmpc2mpl(x=1)[source]

Convert from Mpc^-1 to Mpl (reduced Planck Mass)

pyflation.helpers.ispower2(n)[source]

Returns the log base 2 of n if n is a power of 2, zero otherwise.

Note the potential ambiguity if n==1: 2**0==1, interpret accordingly.

pyflation.helpers.klegend(ks, mpc=False)[source]

Return list of string representations of k modes for legend.

pyflation.helpers.mpl2invmpc(x=1)[source]

Convert from Mpl (reduced Planck Mass) to Mpc^-1

pyflation.helpers.nanfillstart(a, l)[source]

Return an array of length l by appending array a to end of block of NaNs along axis 0.

Parameters :

a: numpy array :

array to append to end of block of NaNs.

l: integer :

length of the new array to return.

Returns :

c: numpy_array :

array of length l with initial values filled with NaNs and array a at the end.

pyflation.helpers.removedups(l)[source]

Return an array with duplicates removed but order retained.

An array is returned no matter what the input type. The first of each duplicate is retained.

Parameters :

l: array_like :

Array (or list etc.) of values with duplicates that are to be removed.

Returns :

retlist: ndarray :

Array of values with duplicates removed but order intact.

pyflation.helpers.seq(min=0.0, max=None, inc=1.0, type=<type 'float'>, return_type='NumPyArray')[source]

Generate numbers from min to (and including!) max, with increment of inc. Safe alternative to arange. The return_type string governs the type of the returned sequence of numbers (‘NumPyArray’, ‘list’, or ‘tuple’).

pyflation.helpers.startlogging(log, logfile, loglevel=20, consolelevel=None)[source]

Start the logging system to store rotational file based log.

rk4 Module

rk4.py - Runge-Kutta ODE solver

Provides Runge-Kutta based ODE solvers for use with pyflation models.

exception pyflation.rk4.SimRunError[source]

Bases: exceptions.StandardError

Generic error for model simulating run. Attributes include current results stack.

pyflation.rk4.rk4stepks(x, y, h, dydx, dargs, derivs)[source]

Do one step of the classical 4th order Runge Kutta method, starting from y at x with time step h and derivatives given by derivs

pyflation.rk4.rkdriver_tsix(ystart, simtstart, tsix, tend, allks, h, derivs)[source]

Driver function for classical Runge Kutta 4th Order method. Uses indexes of starting time values instead of actual times. Indexes are number of steps of size h away from initial time simtstart.

sohelpers Module

sohelpers.py - Second order helper functions

Provides helper functions for second order data from cosmomodels.py.

pyflation.sohelpers.combine_source_and_fofile(sourcefile, fofile, newfile=None)[source]

Copy source term to first order file in preparation for second order run.

Parameters :

sourcefile: String :

Full path and filename of source file.

fofile: String :

Full path and filename of first order results file.

newfile: String, optional :

Full path and filename of combined file to be created. Default is to place the new file in the same directory as the source file with the filename having “src” replaced by “foandsrc”.

Returns :

newfile: String :

Full path and filename of saved combined results file.