libaspic

NAME
SYNOPSIS
DESCRIPTION
EXAMPLES
LIST OF MODULES AND ROUTINES
LIST OF MODELS
ADDING A MODEL
NOTES
AUTHORS
REPORTING BUGS
COPYRIGHT
SEE ALSO

NAME

libaspic - a scientific fortran library dedicated to (a)ccurate (s)low-roll (p)redictions for (i)nflationary (c)osmology.

SYNOPSIS

libaspic is a collection of fortran modules, distributed as a shared and static library, which provides various computational functions associated with various models of inflation. The list of models, modules and available functions, as well as their usage, are detailed below.

In order to access the library, it has to be linked to your code by using the -l flag during compilation, i.e. by appending the command -laspic to the others flags. Inside your source code, all modules and their associated routines can be imported from a global include file that is accessible in fortran by:

include ’aspic.h’

or using the C pre-processor with

#include <aspic.h>

Both the library and the global include file are installed in the directory you specified during the installation. Compiling your own source file with, for instance, gfortran should only require:

gfortran -c mysrcfile.f90
gfortran mysrcfile.o -o myprog -laspic

In the situation in which the library and modules are installed in a non-standard location, you can use the -I command line flag for specifying the path to the module and include file. The same remark holds for the path to the library itself, which can be specified using the -L command line flag. Compiling with gfortran would now read:

gfortran -I/usr/include/aspic -c mysrcfile.f90
gfortran -L/usr/lib64 mysrcfile.o -o myprog -laspic

Alternatively, you can decide to import only one (or more) specific module with a subset of its associated routines with the fortran instruction ’use’. As an example, accessing the first and second Hubble flow functions for the model ’foo’ reads:

use foosr, only : foo_epsilon_one, foo_epsilon_two

DESCRIPTION

libaspic is a library for computing various observable quantities used in Cosmology from definite single field inflationary models. It aims at providing an efficient, extendable and accurate way of comparing theoretical predictions with cosmological data. As the list of inflationary models is always increasing, you are encouraged to add support for any model that would not yet be implemented; see section ADDING A MODEL.

By observable quantities, we currently refer to as the Hubble flow functions, up to second order in the slow-roll approximation, which are in direct correspondence with the spectral index, tensor-to-scalar ratio and the running of the primordial power spectrum. The library also provides the field potential, its first and second derivatives, the energy density at the end of inflation, the energy density at the end of reheating, and the field value (or e-fold value) at which the pivot scale crossed the Hubble radius during inflation. All these quantities are computed in a way which is consistent with the existence of a reheating phase. Some functions may therefore requires as an input both the amplitude of the scalar perturbations and the reheating parameters.

EXAMPLES

For each model implemented, assuming its acronym to be ’foo’, the source tarball contains a worked out example, ’foomain.f90’, which calls the library to extract the spectral index and tensor-to-scalar ratio scanning all possible reheating energy scales. These example codes are not compiled by default as they are for debugging and pedagogical purposes only. Their compilation can be forced by using

make check

when building libaspic from the tarball source file. All these programs are automatically compiled and executed by a

make test

that will return an error message if one of them aborts unexpectedly.

LIST OF MODULES AND ROUTINES

For a given inflationary model, let’s say ’foo’ inflation, libaspic provides two modules foosr and fooreheat respectively dedicated to slow-roll and reheating related functions. These modules make use of common modules that you should not need for standard usage but that you must call if you decide to add a new model.

SLOW-ROLL MODULE
The foosr module encapsulates the following public functions of kind real(kp)

foo_norm_potential(x,...)
foo_norm_deriv_potential
(x,...)
foo_norm_deriv_second_potential
(x,...)
foo_epsilon_one
(x,...)
foo_epsilon_two
(x,...)
foo_epsilon_three
(x,...)
foo_x_endinf
(...)
foo_efold_primitive
(x,...)
foo_x_trajectory
(N-Nend,xend,...)

All these functions take as an input the dimensionless field value ’x’ and/or a relative e-fold number ’N’ and some extra model parameter values, all of the kind real(kp). The field values ’x’ are usually assumed to be in reduced Planck mass units. For some specific models, the unit may however be in another more relevant energy scale. Please check out the explicit list of models below, their respective man pages (man aspic_foo(3)) or the header of the module file ’foosr.f90’ for more details. Each of these functions usually requires some additional input, referred to as ’...’ above, which are the (dimensionless) potential parameter values. Again, their number, units and kind can be either found in the man pages or in the header of the module file ’foosr.f90’.

The first three functions return the normalized field potential and its derivatives with respect to the input variable ’x’. They do not include any multiplicative factor (such as M^4) because it is always uniquely determined by the amplitude of the cosmological perturbations. As such, the normalization is never counted as a model parameter in the classification below. The next three functions, the ones containing the wording "epsilon" in their name, return the value of the three first Hubble flow functions at leading order (the so-called "epsilonV" functions)

Finally, the last three functions allows to determine uniquely the slow-roll trajectory at leading order. Calling foo_x_endinf() returns the field value ’xend’ at the end of inflation. For most models, this happens when the first Hubble flow function equals unity but it can also be given by more specific mechanisms. Notice that, for some models, the end of inflation may be an extra model parameter. For those, this function is obviously not provided. The function foo_efold_primitive() returns, up to a constant, the value of the integral SUM{V(x)/V’(x)dx} from which the slow-roll trajectory can be calculated. This is precisely the output of foo_x_trajectory() which returns the field value ’x’ from the input of ’xend’ and the number of e-folds ’N-Nend’ before the end of inflation at which you want this value.

REHEATING MODULES
The fooreheat module encapsulates the following public functions of kind real(kp)

foo_lnrhoreh_max(...,xend,Pstar)
foo_x_star
(...,xend,w,lnRhoReh,Pstar)
foo_x_star
(...,xend,w,lnRhoReh,Pstar, deltaNstar)
foo_x_rrad
(...,xend,lnRrad,Pstar)
foo_x_rrad
(...,xend,lnRrad,Pstar, deltaNstar)
foo_x_rreh
(...,xend,lnR)
foo_x_rreh
(...,xend,lnR, deltaNstar)

These functions take as first arguments the dimensionless potential parameters ’...’ as specified in the man pages aspic_foo(3) and in the header of the module files ’foosr.f90’ and ’fooreheat.f90’. They also require to specify the field value at the end of inflation ’xend’.

The function foo_lnrhoreh_max() returns the natural logarithm of the total energy density of the universe at the end of reheating when it occurs instantaneously at the end of inflation. There, this is also the energy density at the end of inflation when reheating is instantaneous, or radiation dominated. This number depending on the physical normalization of the potential, you need to input ’Pstar’ the measured amplitude of the scalar power spectrum evaluated at ’kstar’ the pivot wave-number. Its default value has been set to 0.05 Mpc^-1 (see below for specifying another value).

The function foo_x_star() returns the field value ’xstar’ at which the pivot wave-number ’kstar’ crossed the Hubble radius during inflation. Plugging this field value into the Hubble flow functions immediately gives the observable slow-roll parameters, spectral index, running, tensor-to-scalar ratio. As an input, this function requires some assumptions on how the reheating proceeded. It needs the mean equation of state parameter ’w’ during (pre)reheating, together with the logarithm of total energy density ’lnRhoReh’ of the universe when the reheating ends. Finally, in order to determine the correct normalization of the inflationary potential, you have to input ’Pstar’ again. The same routine can be called with an additional real(kp), optional, intent(out) argument ’deltaNstar’ which contains on return the value of ’Nstar-Nend’, the number of e-folds before the end of inflation at which the pivot wave-number crossed the Hubble radius (negative).

The functions foo_x_rrad() and foo_x_rreh() are in all points similar to the previous one, i.e. they return the field value ’xstar’ at which the pivot wave-number ’kstar’ crossed the Hubble radius during inflation. They take as input the reheating parameter ’lnRrad’, or the rescaled reheating parameter ’lnR’, respectively. These parametrizations are most generic as they are the combination of reheating parameters the CMB is sensitive to. For more details, see the references below.

The srreheat module is not model specific and its source files are located under the directory ’src/common/’. Unless otherwise specified, this module encapsulates functions of kind real(kp) which are called by all the above-described modules. As such their usage should be necessary only if you decide to add a new model:

logical :: slowroll_validity(epsOne,epsTwo)
potential_normalization
(Pstar,epsOneStar,Vstar)
primscalar
(M,epsOneStar,Vstar)
log_energy_reheat_ingev
(lnRhoReh)
ln_rho_endinf
(Pstar,epsOneStar,epsOneEnd,VendOverVstar)
ln_rho_reheat
(w,Pstar,epsOneStar,epsOneEnd,deltaNstar,VendOverVstar)
find_reheat
(nuStar,calFplusNuEnd,w,epsStar,Vstar)
get_calfconst
(lnRhoReh,Pstar,w,epsEnd,potEnd, lnOmega4End )
find_reheat_rrad
(nuStar,calFplusNuEnd,epsStar,Vstar)
get_calfconst_rrad
(lnRrad,Pstar,epsEnd,potEnd)
find_reheat_rreh
(nuStar,calFplusNuEnd,Vstar)
get_calfconst_rreh
(lnR,epsEnd,potEnd, lnOmega4End )
get_lnrrad_rhow
(lnRhoReh,w,lnRhoEnd)
get_lnrreh_rhow
(lnRhoReh,w,lnRhoEnd)
get_lnrrad_rreh
(lnR,lnRhoEnd)
get_lnrreh_rrad
(lnRrad,lnRhoEnd)

All of these functions take as input real(kp) kind arguments. The first function slowroll_validity() returns .true. or .false. according to the values of the first and second Hubble flow functions to assess the validity of the slow-roll approximation and numerical precision. The second function potential_normalization() returns the potential normalization factor required to get the correct amplitude of the CMB anisotropies. This factor is commonly denoted as ’M^4’ and this function returns the mass scale ’M’ in Planck units. Conversely, the function primscalar() returns the amplitude of the primordial scalar perturbations at the pivot scale from the input of the potential mass scale ’M’. The next function log_energy_reheat_ingev() is for convenience and simply returns the logarithm in base 10 of the energy density at the end of reheating from the its natural logarithmic value in Planck units (used elsewhere). The next functions are at the root of the reheating related calculations and are fully model independent. The function ln_rho_endinf() returns the logarithm of the energy density at the end of inflation, ln_rho_reheat() returns the logarithm of the energy density at the end of reheating, while find_reheat() and get_calfconst() solve algebraic equations necessary to get the reheating parameter assuming slow-roll. Some function takes an optional argument, ’lnOmega4End’, which is the logarithm of the power fourth of the conformal factor for scalar-tensor theories of gravity. For more details on what are these quantities, see the references at the end of this section. The next four functions equally solve the reheating equations but take as input either the reheating parameter ’lnRrad’, or the rescaled one ’lnR’. Finally, the last four functions allow to pass from one reheating variable to the others. For instance, get_lnrrad_rhow() gives the reheating parameter ’lnRrad’ from the value of ’lnRhoReh’ and ’w’. Notice that the energy scale at which inflation ends, ’lnRhoEnd’, is a required input for all the conversion functions but can be computed with ln_rho_endinf().

All these routines are valid for any slow-roll inflationary models. The quantity ’Pstar’ stands for the primordial power spectrum amplitude at the pivot, ’w’ the mean equation of state during (pre)reheating, ’epsOneStar’ and ’epsOneEnd’ are the first Hubble flow function respectively evaluated at the time the pivot mode crossed the Hubble radius during inflation, and at the end of inflation. The argument ’VendOverVstar’ is the ratio between the field potential, evaluated at those two times. All those arguments are of real(kp) kind.

The srflow module provides some potentially useful functions to get other cosmological observables as well as higher order corrections on the Hubble flow functions. Its source file is located under the directory ’src/common/’. In particular, the module has the public functions of kind real(kp)

scalar_spectral_index(epsH)
tensor_to_scalar_ratio
(epsH)
scalar_running
(epsH)

All of these functions take as input a

real(kp), dimension (:) :: epsH

vector assumed to contain the value of the successive Hubble flow parameters ’epsilonH_i’, with ’i’ increasing. The calculations are consistent with the size of the input vector. For instance, calling scalar_spectral_index() with a dimension two vector containing the values of the first and second Hubble flow parameters returns the spectral index computed at first order in a Hubble flow expansion. If you input a dimension three vector, the calculations are performed at second order. The same holds for the tensor-to-scalar ratio and the running of the spectral index (which is non-zero at second order only).

References:

arXiv:1303.3787 (section 2.2)
arXiv:1302.6013 (section 2.2)
arXiv:1301.1778 (section IIA)
arXiv:1202.3022 (section 2)
arXiv:1009.4157 (section IIB)
arXiv:1004.5525 (whole paper)
arXiv:0711.4307 (section 2.4)
astro-ph/0703486 (section 4.1)
astro-ph/0605367 (section 4.1)

COSMOPAR MODULE
The cosmopar module encapsulates some public parameters of the kind real(kp) which encodes some measured cosmological parameters today, observational choices such as the pivot scale, and some particle physics constant. More explicitly, they are

HubbleSquareRootOf3OmegaRad
HubbleSquareRootOf2OmegaRad
RelatDofRatio
lnRhoNuc
lnMpcToKappa
lnMpinGeV
QrmsOverT
kpivot
PowerAmpScalar
HiggsVeV
HiggsMass
HiggsCoupling

The first two are the Hubble parameter today times the square root of the double (or triple) density parameter of radiation today, the second is the ratio between the number of entropic relativistic species at the end of reheating and today (gives only small corrections). The constant lnRhoNuc stands for the natural logarithm of the energy density of the universe just before Big-Bang Nucleosynthesis. Next lnMpcToKappa is the logarithm of the Einstein equation coupling (8piG/c^4) expressed in mega-parsecs. The parameter lnMpinGev is the reduced Planck mass in GeV, QrmsOverT stands for the effective quadrupole moment, kpivot is the pivot scale at which the amplitude of the scalar primordial power spectrum is measured. A default amplitude is stored in the parameter PowerAmpScalar (mean value from PLANCK 2013), that very same quantity has been referred to as ’Pstar’ in some functional arguments above. The effective quadrupole moment ’QrmsOverT’ is such that the amplitude of the power spectrum matches ’Pstar’. As such it may not correspond to the real quadrupole moment, which is still slightly lower :-). HiggsVeV is the vacuum expectation value of the Higgs field, in reduced Planck units. HiggsMass is the mass of the Higgs boson, again in reduced Planck mass. Finally HiggsCoupling is the quartic coupling constant of the Higgs field.

Notice that changing any of these constants requires edition of the source file ’src/common/cosmopar.f90’ and a recompilation of the whole library.

UTILITY MODULES
Finally, libaspic comes with some utility modules that you may find useful in performing some specific computations.

The inftools module encapsulates some public subroutines which are various modified Runge-Kutta numerical integrators based on the subroutine dverk(). The specialinf module encapsulates some special functions arising by analytically integrating some slow-roll trajectories. The hyp_2f1_module module encapsulates various functions and subroutines dedicated to the computation of the Gauss hyper-geometric function. All source files are located under the ’src/common/’ directory.

LIST OF MODELS

At the time of this writing, libaspic deals with the inflationary models listed below. Their respective potential parameters, conventions for field units and so on, are described in their man pages aspic_foo(3).

ZERO PARAMETER MODELS

Acronym

Model name

si

Starobinsky Inflation

hi

Higgs inflation

ONE PARAMETER MODELS

Acronym

Model name

rchi

radiatively corrected Higgs inflation

lfi

large field inflation

mlfi

mixed large field inflation

rcmi

radiatively corrected massive inflation

rcqi

radiatively corrected quartic inflation

ni

natural inflation

esi

exponential SUSY inflation

pli

power law inflation

kmii

Kahler moduli inflation I

hf1i

horizon flow inflation at first order

cwi

Coleman-Weinberg inflation

li

global SUSY with loop inflation

rpi1

R + R^2p inflation I

dwi

double well inflation

mhi

mutated hilltop inflation

rgi

radion gauge inflation

mssmi

minimal supersymmetric model inflation

ripi

renormalizable inflection point inflation

ai

arctan inflation

cnai

constant spectral index inflation A

cnbi

constant spectral index inflation B

osti

open string tachyonic inflation

wri

Witten-O’Raifeartaigh inflation

ccsi1

cublicly corrected Starobinski inflation I

ccsi3

cublicly corrected Starobinski inflation III

di

dual inflation

sbki

symmetry breaking Kahler inflation

ahi

axion hilltop inflation

pai

pure arctangent inflation

saai

superconformal alpha attractor A inflation

sati (alpha=1)

T-model inflation

TWO PARAMETERS MODELS

Acronym

Model name

sfi

small field inflation

ii

intermediate inflation

kmiii

Kahler moduli inflation II

lmi1

logamediate inflation I

rpi

R + R^2p inflation

twi

twisted inflation

hf2i

horizon flow inflation at second order

gmssmi

generalized minimal supersymmetric model inflation

gripi

generalized renormalizable point inflation

bsusybi

brane SUSY breaking inflation

ti

tip inflation

bei

beta exponential inflation

psni

pseudo natural inflation

ncki

non-canonical Kahler inflation

csi

constant spectrum inflation

oi

orientifold inflation

cnci

constant spectral index inflation C

sbi

supergravity brane inflation

ssbi

spontaneous symmetry breaking inflation

imi

inverse monomial inflation

bi

brane inflation

kklti

Kachru Kallosh Linde Trivedi inflation

nfi1

N-formalism inflation I

nfi3

N-formalism inflation III

ccsi2

cublicly corrected Starobinski inflation II

hni

hybrid natural inflation

vfmi

Viatcheslav Fyodorovich Mukhanov inflation

fi

fibre inflation

hbi

hyberbolic inflation

shi

smeared Higgs inflation

dei

double exponential inflation

sdi

s-dual inflation gwdi generalized double well inflation

nmlfi

non-minimal large field inflation

sabi

superconformal alpha attractor B inflation

sati

superconformal alpha attractor T inflation

THREE PARAMETERS MODELS

Acronym

Model name

lmi2

logamediate inflation II

rmi

running mass inflation

vhi

valley hybrid inflation

dsi

dynamical supersymmetric inflation

gmlfi

generalized mixed large field inflation

lpi

logarithmic potential inflation

cndi

constant spectral index inflation D

saii3

string axion inflation II

rclfi

radiatively corrected large field inflation

ncli

non-renormalizable corrected loop inflation

hni

hybrid natural inflation

nfi

N-formalism inflation

rcipi

radiatively corrected inflection point inflation

ADDING A MODEL

Before deciding to add a model, you should first check that its potential is not already encoded within the existing modules. From our experience, it is frequent in the literature that different theoretical motivations lead to exactly the same effective potential. As a result, identical models often share different names. If you encounter such a situation, please let us know, or even better, send us an updated man page for the relevant module by adding the alternative names under which this potential is known.

In the opposite situation, importing a new model, let’s say ’convoluted wow loop inflation’, of acronym wowi is equivalent to write the source codes of the two modules wowisr and wowireheat as well as updating various autoconf files, namely ’Makefile.am’ and ’configure.ac’, and finally writing a very short documentation.

This can be done step by step along the following lines:

Create the sub-directory ’src/wooi’ containing five new files, ’wooimain.f90’, ’wooisr.f90’, ’wooireheat.f90’, ’aspic_wooi.3’ and ’Makefile.am’.

Edit the file ’Makefile.am’ such as it now reads

SRC = wooisr.f90 wooireheat.f90
MOD = wooisr.$(FC_MODEXT) wooireheat.$(FC_MODEXT)

check_PROGRAMS = wooimain
wooimain_SOURCES = $(SRC) wooimain.f90
wooimain_FCFLAGS = -I../$(SRCOMMDIR)
wooimain_LDADD = ../$(SRCOMMDIR)/libsrcommon.a

noinst_LTLIBRARIES = libwooi.la
libwooi_la_SOURCES = $(SRC)
libwooi_la_FCFLAGS = -I../$(SRCOMMDIR) $(AM_FCFLAGS)
libwooi_la_includedir = $(includedir)/$(SRINCDIR)
libwooi_la_include_HEADERS = $(MOD)

man_MANS = aspic_wooi.3

clean-local: clean-modules clean-outfiles
clean-modules:
test -z "$(FC_MODEXT)" || $(RM) *.$(FC_MODEXT)
clean-outfiles:
test -z "$(DATEXT)" || $(RM) *.$(DATEXT)  
.NOTPARALLEL:

Edit the files ’wooisr.f90’ and ’wooireheat.f90’ such that they respectively provide the wooisr and wooireheat modules and their respective public functions starting with the wowi acronym. The best way to do this is to copy-paste the files of one of the existing model and modify them accordingly. You must use the already common routines for this, such as zbrent() is you need to solve algebraic equations or get_calfconst() and find_reheat() to solve for the reheating. You may also need some special functions that are already encoded in the specialinf module. In the unlikely situation in which you would need a special function or another solver, you should add it into the relevant modules (located in ’src/common’) and render public those new functions.

Write the test program ’wooimain.f90’ to check that your code is actually working and produce sensible results. Again you may be inspired by the already encoded models.

Document your model, i.e. write the mini man page in the file ’aspic_wooi.3’ summarizing the potential functional shape, the number and kind of the parameters, as well as the physical units used.

Add your model to the library by editing the parent Makefile ’src/Makefile.am’. Update the environment variable libaspic_la_LIBADD by adding the line ’wooi/libwooi.la’ and append to SUBDIRS the name of the new sub-directory ’wooi’.

Finally, edit the global ’configure.ac’ file and run the command autoreconf such that the autoconf tools can automatically generate the various makefiles.

And send us your code, we will be happy to add it, as your name, in the next release of libaspic

NOTES

Please help us to maintain this library readable. As such, we strongly encourage the use of modern fortran and will not accept routines written in f66 or f77. The only exception might be for the fantastic two-century old hyper fast routines, under the condition that you provide them enclosed into a module box with a maximal amount of private routines. If you are not (yet) familiar with fortran 90/95/03/08 and later revisions, check out the tutorials from the IDRIS (in french).

AUTHORS

libaspic has been written by:

Name

Affiliation

Jerome Martin

Institut d’Astrophysique de Paris (France)

Christophe Ringeval

Cosmology, Universe and Relativity at Louvain, Louvain University (Belgium)

Vincent Vennin

Laboratoire de Physique de l’Ecole Normale Superieure, ENS, Paris (France)

REPORTING BUGS

Please contact us in case of bugs.

COPYRIGHT

GNU GENERAL PUBLIC LICENSE Version 3

SEE ALSO

aspic_si(3), aspic_hi(3), aspic_lfi(3), aspic_rcmi(3), aspic_rcqi(3), aspic_ni(3), aspic_esi(3), aspic_pli(3), aspic_kmii(3), aspic_hf1i(3), aspic_cwi(3), aspic_li(3), aspic_rpi1(3), aspic_dwi(3), aspic_mhi(3), aspic_rgi(3), aspic_mssmi(3), aspic_ripi(3), aspic_ai(3), aspic_cnai(3), aspic_cnbi(3), aspic_osti(3), aspic_wri(3), aspic_ccsi1(3), aspic_ccsi3(3), aspic_di(3), aspic_sbki(3), aspic_ahi(3), aspic_pai(3), aspic_saai(3).

aspic_sfi(3), aspic_ii(3), aspic_kmiii(3), aspic_lmi1(3), aspic_rpi2(3), aspic_rpi3(3), aspic_twi(3), aspic_hf2i(3), aspic_gmssmi(3), aspic_gripi(3), aspic_bsusybi(3), aspic_ti(3), aspic_bei(3), aspic_psni(3), aspic_ncki(3), aspic_csi(3), aspic_oi(3), aspic_cnci(3), aspic_sbi(3), aspic_ssbi(3), aspic_imi(3), aspic_bi(3), aspic_kklti(3), aspic_nfi1(3), aspic_nfi3(3), aspic_ccsi2(3), aspic_hni(3), aspic_vfmi(3), aspic_sdi(3), aspic_fi(3), aspic_sati(3), aspic_saci(3), aspic_hbi(3), aspic_shi(3), aspic_dei(3), aspic_sdi(3), aspic_gdwi(3), aspic_nmlfi(3), aspic_sabi(3), aspic_sati(3).

aspic_lmi2(3), aspic_rmi(3), aspic_vhi(3), aspic_dsi(3), aspic_gmlfi(3), aspic_lpi(3), aspic_cndi(3), aspic_saii3(3), aspic_rclfi(3), aspic_ncli(3), aspic_hni(3), aspic_nfi(3), aspic_rcipi(3).