mimes is hosted by Hepforge, IPPP Durham

Classes in python

All python classes can be used similarly to C++. Obviously, since there are no tempates and pointers in python, the syntax is a bit different. Regarding the template arguments, they are chosen at compile-time using some macros, which can be adjusted from the Definitions.mk file inside the root directory of MiMeS. Furthermore, in the corresponding Axion class, instead of pointer to AxionMass instance we just pass the instance by value.

In order to be able to load the various modules, the following has to be included at the top of the script

    from sys import path as sysPath
    sysPath.append('path_to_src')

where 'path_to_src' is the relative or absolute path to MiMeS/src.

The Cosmo class

The Cosmo class is used to interpolate various quantities of the plasma.Planck Collaboration. The class can be imported by running

    from interfacePy.Cosmo import Cosmo
Its constructor is

    Cosmo(cosmo_PATH, minT=0, maxT=mP)
The argument cosmo_PATH is the path (a string) of a data file that contains TT (in GeV), heffh_{eff}, geffg_{eff}, with accenting TT. The second and third arguments, minT and maxT, are minimum and maximum interpolation temperatures, with the interpolation being between the closest temperatures in the data file. Moreover, beyond these limits, both heffh_{eff} and geffg_{eff} are assumed to be constants. It is important to note that the class creates a void pointer that gets recasted to mimes::Cosmo< LD > in order to call the various member functions. This means that once an instance of Cosmo is no longer needed, it must be deleted, in order to free the memory that it occupies. An instance, say cosmo, is deleted using

    del cosmo

Interpolation of the RDOF, allows us to define various quantities related to the plasma. These quantities are given as the member functions:

  1. Cosmo.heff(T): heffh_{eff} as a function of TT.
  2. Cosmo.geff(T): geffg_{eff} as a function of TT.
  3. Cosmo.dheffdT(T): dheff/dTdh_{eff}/dT as a function of TT.
  4. Cosmo.dgeffdT(T): dgeff/dTdg_{eff}/dT as a function of TT.
  5. Cosmo.dh(T): δh=1+13dlogheffdlogT\delta_h = 1 + \frac{1}{3} \frac{d\log h_{eff}}{d\log T} as a function of TT.
  6. Cosmo.s(T): The entropy density of the plasma as a function of TT.
  7. Cosmo.rhoR(T): The energy density of the plasma as a function of TT.
  8. Cosmo.Hubble(T): The Hubble parameter assuming radiation dominated expansion as a function of TT.

Moreover, there are several cosmological quantities are given as members variables:1

  1. Cosmo.T0: CMB temperature today in GeV.
  2. Cosmo.h_hub: Dimensionless Hubble constant.
  3. Cosmo.rho_crit: Critical density in GeV3^3.
  4. Cosmo.relicDM_obs: Central value of the measured DM relic abundance.
  5. Cosmo.mP: Planck mass in GeV.
Note that these values can be directly imported from the module, without declaring an instance of the class, as

    from interfacePy.Cosmo import T0, h_hub, rho_crit, relicDM_obs, mP

The AxionMass class

The AxionMass class is defined in the module with the same name that can be found in the directory MiMeS/src/interfacePy/AxionMass. This class is responsible for the definition of the axion mass. This module loads the corresponding shared library from MiMeS/lib/libma.so, which is created by compiling MiMeS/src/AxionMass/AxionMass.cpp using make lib/libma.so. Its usage is described in the Examples. Moreover, this class is used in the same way as mimes::AxionMass< LD >. However, we should append in this section the definition of its member functions in python .

The constructor is

    AxionMass(*args)

and can be used in two different ways.

First, one can pass three arguments, as

    AxionMass(chi_PAT, minT, maxT)
The first argument is the path to a data file that contains two columns; TT (in GeV) and chichi (in GeV4^4), with increasing TT. The arguments minT and maxT are the interpolation limits. These limits are used in order to stop the interpolation in the closest temperatures in chi_PATH. That is the actual interpolation limits are TminT_{min}\geqminT and TmaxT_{max}\leqmaxT. Beyond these limits, by default, the axion mass is assumed to be constant. However, this can be changed by using the member functions

    AxionMass.set_ma2_MIN(ma2_MIN)
    AxionMass.set_ma2_MAX(ma2_MAX)

Here, ma2_MIN(T,fa) and ma2_MAX(T,fa), are functions (not any callable object), should take as arguments T and fa, and return the axion mass squared beyond the interpolation limits. In order to ensure that the axion mass is continuous, usually we need TminT_{min}, TmaxT_{max}, χ(Tmin)\chi(T_{min}), and χ(Tmax)\chi(T_{max}). These values can be obtained using the member functions

  1. AxionMass.getTMin(): This function returns the minimum interpolation temperature, TminT_{min}.
  2. AxionMass.getTMax(): This function returns the maximum interpolation temperature, TmaxT_{max}.
  3. AxionMass.getChiMin(): This function returns χ(Tmin)\chi(T_{min}).
  4. AxionMass.getChiMax(): This function returns χ(Tmax)\chi(T_{max}).
An alternative way to define the axion mass is via the constructor

    AxionMass(ma2)

Here, ma2(T,fa) is a function (not any callable object) that takes TT (in GeV) and faf_a, and returns ma2(T)m_a^2(T) (in GeV). As in the other python classes, once the instances of this class are no longer needed, they must be deleted using the destructor, del.

The member function that returns ma2(T)m^2_a(T) is

    AxionMassma2(T,fa)
We should note that another ma2 can be changed using the following member function:

    AxionMass.set_ma2(ma2)

Again, ma2(T,fa) is a function (not any callable object) that takes TT (in GeV) and faf_a, and returns ma2(T)m_a^2(T) (in GeV).

The AxionMass class

The class Axion, solves the axion EOM. This class is defined in MiMeS/interfacePy/Axion/Axion.py, and imports the corresponding shared library. This library is compiled by running make lib/Axion_py.so, and its source file is MiMeS/src/Axion/Axion-py.cpp. As in the previous classes, its usage is similar to the C++ version.

Its constructor is

    Axion(theta_i, fa, umax, TSTOP, 
            ratio_ini, unsigned int N_convergence_max, convergence_lim, 
            inputFile, axionMass, initial_step_size=1e-2, 
            minimum_step_size=1e-8, maximum_step_size=1e-2, absolute_tolerance=1e-8, 
            relative_tolerance=1e-8, beta=0.9, fac_max=1.2, fac_min=0.8, 
            unsigned int maximum_No_steps=10000000)

Again, the various arguments are discussed in the C++ case. Notice one important difference between this and the C++ verison of this class; the instance of AxionMass, axionMass, is passed by value. However, internally, the constructor converts this instance to a pointer, which is then passed to the underlying C function responsible for creating the relevant instance.

The member function responsible for solving the EOM is

    Axion.solveAxion()

Once this function is finished, the following member variables are available

  1. Axion.T_osc: the oscillation temperature, ToscT_{osc}, in GeV.
  2. Axion.a_osc: aai\dfrac{a}{a_i} at the oscillation temperature.
  3. Axion.theta_osc: thetaosctheta_{osc}, i.e. θ\theta at ToscT_{osc}.
  4. Axion.gamma: the entropy injection between the last peak (T=TpeakT=T_{peak}) and today (T=T0T=T_0), γ\gamma (entropy injection factor from TiT_i).
  5. Axion.relic The relic abundance of the axion.
The evolution of a/ai,T,θ,ζ,ρaa/a_i, \ T, \ \theta, \ \zeta, \rho_a at the integration steps, is not automatically accessible to user, but they can be made so using

    Axion.getPoints()
Then, the following member variables are filled
  1. Axion.a: The scale factor over its initial value, aai\dfrac{a}{a_i}.
  2. Axion.T: The temperature in GeV.
  3. Axion.theta: The axion angle, θ\theta.
  4. Axion.zeta: The derivative of θ\theta, ζdθdlog(a/ai)\zeta \equiv \dfrac{d \theta}{d \log (a/a_i)}.
  5. Axion.rho_axion: The axion energy density in GeV4^4.
Moreover, the function

    Axion.getPeaks()

fills the (numpy) arrays Axion.a_peak, Axion.T_peak, Axion.theta_peak, Axion.zeta_peak, Axion.rho_axion_peak, and Axion.adiabatic_invariant with the quantities a/ai,T,θ,ζ,ρa,Ja/a_i, \ T, \ \theta, \ \zeta, \rho_a, \ J, at the peaks of the oscillation. These points are computed using linear interpolation between two integration points with a change in the sign of ζ\zeta.

The local integration errors for θ\theta and ζ\zeta are stored in Axion.dtheta and Axion.dzeta, after the following function is run

    Axion.getErrors()
Another initial condition, θi\theta_i, can be used without declaring a new instance using

    Axion.setTheta_i(theta_i)

We should note that running this function all variables are cleared.

As in the previous py classes, once an instance of the Axion class is no longer needed, it needs to be deleted, by calling the destructor, del.

Important difference between the C++ version Since the axion mass is passed by value in the constructor, a change of the AxionMass instance has no effect on the Axion instance that uses it. Therefore, if the definition of the axion mass changes, one has to declare a new instance of the Axion class. The new instance can be named using the name of the previous one, if the latter is deleted by running its destructor.


  1. Taken from PRD and the Planck Collaboration.