REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
TRestAxionField Class Reference

Detailed Description

A basic class to define analytical axion-photon conversion calculations for axion helioscopes.

TRestAxionField is a class used to calculate the axion-photon mixing and determine the probability of the particle being in a photon state after propagating inside a magnetic field.

The class allows to assign a buffer gas using a TRestAxionBufferGas and a magnetic field map loaded through TRestAxionMagneticField. If these objects have been assigned, the methods implemented in this class may consider those - if needed - inside the calculation.

In practice this class provides different versions of the method GammaTransmissionProbability that allows to calculate the axion-photon probability using different strategies.

Calculating axion-photon probability in a constant field.

1. In vacuum

For calculations inside a constant magnetic field one may simply invoke the following code which will launch the calculation in vacuum, for a field of 2T, a coherence length of 10000mm, and an axion energy and mass of 4.2keV and 0.1eV, respectively.

field->GammaTransmissionProbability( 2, 10000, 4.2, 0.1 );
A basic class to define analytical axion-photon conversion calculations for axion helioscopes.
Double_t GammaTransmissionProbability(Double_t ma, Double_t mg=0, Double_t absLength=0)
Performs the calculation of axion-photon conversion probability using directly equation (11) from van...

It is possible to reduce the number of arguments we are giving to this function by assigning few data member values present inside TRestAxionField as follows:

field.SetMagneticField(2);
field.SetCoherenceLength(10000);
field.SetAxionEnergy(4.2);

Indeed, the complete version with all arguments will also update the data member values, so that next calls to this method will use the same magnetic field, coherence length and energy. In the following example the second call to GammaTransmissionProbability will use again the same values for magnetic field, coherence length and energy from the preceding call.

field.GammaTransmissionProbability( 2, 10000, 4.2, 0.1 );

2. In a buffer gas medium

The axion-photon probability can also be calculated in the sin of a gaseous medium. For that we need to assign a buffer gas instance containing the relevant gas properties in the form of a TRestAxionBufferGas.

The following piece of code will assign the gas properties to the field calculation as defined inside the bufferGases.rml official file. The gas provides the absorption and equivalent photon mass as a funtion of the energy.

TRestAxionBufferGas gas("bufferGases.rml", "helium");
field.AssignBufferGas( &gas );
field.SetMagneticField(2);
field.SetCoherenceLength(10000);
field.SetAxionEnergy(4.2);
A metadata class to define the gas properties used in axion search calculations.
void AssignBufferGas(TRestAxionBufferGas *buffGas)
It assigns a gas buffer medium to the calculation.

Once we have assigned a buffer gas to the class we may return back to the vacuum state by just assigning a nullptr to the gas.

field.AssignBufferGas(nullptr);

Calculating axion-photon probability in an unhomogeneous field.

There are two main strategies presently implemented inside this class to integrate the axion field along an unhomogeneous magnetic field. The first one uses a simple trapezoidal integration method with the advantage that the computing time is controlled by length of the integration step, but no error calculation is given. The second method advantage is that in principle it allows to fix the accuracy of our calculation and an error of the integration is returned by the GSL integration method used. The main disadvantage is that by fixing the accuracy we lose control on the computation time required to reach such accuracy. Here we describe both methods.

1. Using external field profile

Using a external field profile requires to prepare the field profile and pass it along to TRestAxionField::GammaTransmissionProbability method performing the integration. The profile will be placed in a std::vector<Double_t> that describes the magnetic field profile with a fixed length step. The method TRestAxionMagneticField::GetTransversalComponentAlongPath can be used to generate such profile and calculate the field using the first method, as in the following example:

Double_t deltaL = 100; // mm
Double_t Ea = 4.2; // keV
Double_t ma = 0.001; // eV
// Recovering the B-field profile along the axis
TRestAxionMagneticField mField("fields.rml", "babyIAXO_2024");
std::vector<Double_t> bField = mField.GetTransversalComponentAlongPath( TVector3(0,0,-6000),
TVector3(0,0,6000), deltaL );
aField.GammaTransmissionProbability( bField, deltaL, Ea, ma );
A class to load magnetic field maps and evaluate the field on those maps including interpolation.

In a similar way to the constant field case, it would suffice to assign a buffer gas to the TRestAxionField instance to have the calculation performed in a buffer gas.

2. Using a field map instance

A new method has been added in which the integration method is given a pointer to the function that will allow evaluation of the magnetic field along the B-profile parameterized track. Thus, it is the integration method the one requesting the points on the parametric curve that will be necessary to perform the integration. The required field evaluation is encoded inside TRestAxionField::Integrand which accesses the magnetic field instance TRestAxionField::fMagneticField in order to recover a parametric track along the path.

The method name using this technique requires thus direct access to a TRestAxionMagneticField instance that must be assigned beforehand. For computing the integral following this strategy one needs to write the following code where we consider the most generic case including a gas buffer medium.

TRestAxionBufferGas gas("bufferGases.rml", "helium");
aField.AssignBufferGas( &gas );
TRestAxionMagneticField mField("fields.rml", "babyIAXO_2024");
aField.AssignMagneticField( &mField );
Double_t ma = 0.001;
Double_t Ea = 4.2;
Double_t accuracy = 0.1;
// The track that will be integrated must be predefined before inside the magnetic field class
mField.SetTrack( TVector3(0,0,-6000), TVector3(0,0,1) );
std::pair<Double_t, Double_t> prob = aField.GammaTransmissionFieldMapProbability(Ea, ma, accuracy);
std::pair< Double_t, Double_t > GammaTransmissionFieldMapProbability(Double_t Ea, Double_t ma, Double_t accuracy=1.e-1, Int_t num_intervals=100, Int_t qawo_levels=20)
Performs the calculation of axion-photon conversion probability using directly equation (28) from J....
void AssignMagneticField(TRestAxionMagneticField *mField)
It assigns a magnetic field to the calculation.

Where the later call will return a std::pair with the probability and its corresponding error.

Determining density steps for continuous scanning

A method may be used to determine the masses and gas densities to achieve a continuous scanning TRestAxionField::GetMassDensityScanning. This method will place a new mass or gas density setting at FWHM/2 till it reaches the maximum axion mass specified in the argument.

The following code recovers the axion masses and density settings required for a continuous scan till 0.2eV.

std::pair <Double_t,Double_t> aField.GetMassDensityScanning( "He", 0.2 );
std::vector< std::pair< Double_t, Double_t > > GetMassDensityScanning(std::string gasName="He", Double_t maMax=0.15, Double_t rampDown=5.0)
This method determines the proper densities to be used in an axion helioscope experiment in order to ...

RESTsoft - Software for Rare Event Searches with TPCs

History of developments:

2019-March: First concept and implementation of TRestAxionField class. Javier Galan 2023-December: Implementing methods to recover the axion mass scanning nodes Fran Candon

Author
Javier Galan
Fran Candon

Definition at line 30 of file TRestAxionField.h.

#include <TRestAxionField.h>

Inheritance diagram for TRestAxionField:

Public Member Functions

void AssignBufferGas (TRestAxionBufferGas *buffGas)
 It assigns a gas buffer medium to the calculation. More...
 
void AssignMagneticField (TRestAxionMagneticField *mField)
 It assigns a magnetic field to the calculation. More...
 
Double_t AxionAbsorptionProbability (Double_t Bmag, Double_t Lcoh, Double_t Ea, Double_t ma, Double_t mg=0, Double_t absLength=0)
 On top of calculating the axion absorption probability it will assign new values for the magnetic field (Bmag/T), coherence length (Lcoh/mm) and axion energy (Ea/keV). More...
 
Double_t AxionAbsorptionProbability (Double_t ma, Double_t mg=0, Double_t absLength=0)
 Performs the calculation of axion-photon absorption probability using directly equation (18) from van Bibber, Phys Rev D Part Fields. 1989. More...
 
Double_t BL (Double_t Bmag, Double_t Lcoh)
 Performs the calculation of (BL) factor in natural units. More...
 
Double_t BLHalfSquared (Double_t Bmag, Double_t Lcoh)
 Performs the calculation of (BL/2)^2 factor in natural units. More...
 
std::pair< Double_t, Double_t > GammaTransmissionFieldMapProbability (Double_t Ea, Double_t ma, Double_t accuracy=1.e-1, Int_t num_intervals=100, Int_t qawo_levels=20)
 Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf. More...
 
Double_t GammaTransmissionFWHM (Double_t step=0.00001)
 Performs the calculation of the FWHM for the axion-photon conversion probability computed in TRestAxionField::GammaTransmissionProbability. More...
 
Double_t GammaTransmissionProbability (Double_t Bmag, Double_t Lcoh, Double_t Ea, Double_t ma, Double_t mg=0, Double_t absLength=0)
 On top of calculating the gamma transmission probability it will assign new values for the magnetic field (Bmag/T), coherence length (Lcoh/mm) and axion energy (Ea/keV). More...
 
Double_t GammaTransmissionProbability (Double_t ma, Double_t mg=0, Double_t absLength=0)
 Performs the calculation of axion-photon conversion probability using directly equation (11) from van Bibber, Phys Rev D Part Fields. 1989. More...
 
Double_t GammaTransmissionProbability (std::vector< Double_t > Bmag, Double_t deltaL, Double_t Ea, Double_t ma, Double_t mg=0, Double_t absLength=0)
 Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf. More...
 
Double_t GetAxionEnergy () const
 
Double_t GetCoherenceLength () const
 
Double_t GetMagneticField () const
 
std::vector< std::pair< Double_t, Double_t > > GetMassDensityScanning (std::string gasName="He", Double_t maMax=0.15, Double_t rampDown=5.0)
 This method determines the proper densities to be used in an axion helioscope experiment in order to achieve a continuous axion mass scan. More...
 
void SetAxionEnergy (Double_t e)
 
void SetBufferGas (TRestAxionBufferGas *buffGas)
 It assigns a gas buffer medium to the calculation. More...
 
void SetCoherenceLength (Double_t l)
 
void SetDebug (Bool_t v)
 It enables/disables debug mode. More...
 
void SetMagneticField (Double_t b)
 
 TRestAxionField ()
 Default constructor. More...
 
 ~TRestAxionField ()
 Default destructor. More...
 

Static Public Member Functions

static double Integrand (double x, void *params)
 Integrand used for axion-photon probability integration. More...
 

Private Member Functions

std::pair< Double_t, Double_t > ComputeOffResonanceIntegral (Double_t q, Double_t Gamma, Double_t accuracy, Int_t num_intervals, Int_t qawo_levels)
 Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf. More...
 
std::pair< Double_t, Double_t > ComputeResonanceIntegral (Double_t Gamma, Double_t accuracy, Int_t num_intervals)
 Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf. More...
 
void Initialize ()
 Initialization of TRestAxionField class. More...
 

Private Attributes

Double_t fBmag = 2.5
 The magnetic field in Teslas (used for constant field formulas) More...
 
TRestAxionBufferGasfBufferGas = nullptr
 A pointer to the buffer gas definition. More...
 
Bool_t fDebug = false
 
Double_t fEa = 4.2
 The energy of the axion in keV. More...
 
Double_t fLcoh = 10000
 The coherence lenght (in mm) where the magnetic field is defined (for constant field) More...
 
TRestAxionMagneticFieldfMagneticField = nullptr
 A pointer to the magnetic field definition. More...
 

Constructor & Destructor Documentation

◆ TRestAxionField()

TRestAxionField::TRestAxionField ( )

Default constructor.

Definition at line 227 of file TRestAxionField.cxx.

◆ ~TRestAxionField()

TRestAxionField::~TRestAxionField ( )

Default destructor.

Definition at line 232 of file TRestAxionField.cxx.

Member Function Documentation

◆ AssignBufferGas()

void TRestAxionField::AssignBufferGas ( TRestAxionBufferGas buffGas)
inline

It assigns a gas buffer medium to the calculation.

Definition at line 73 of file TRestAxionField.h.

◆ AssignMagneticField()

void TRestAxionField::AssignMagneticField ( TRestAxionMagneticField mField)
inline

It assigns a magnetic field to the calculation.

Definition at line 76 of file TRestAxionField.h.

◆ AxionAbsorptionProbability() [1/2]

Double_t TRestAxionField::AxionAbsorptionProbability ( Double_t  Bmag,
Double_t  Lcoh,
Double_t  Ea,
Double_t  ma,
Double_t  mg = 0,
Double_t  absLength = 0 
)

On top of calculating the axion absorption probability it will assign new values for the magnetic field (Bmag/T), coherence length (Lcoh/mm) and axion energy (Ea/keV).

Definition at line 743 of file TRestAxionField.cxx.

◆ AxionAbsorptionProbability() [2/2]

Double_t TRestAxionField::AxionAbsorptionProbability ( Double_t  ma,
Double_t  mg = 0,
Double_t  absLength = 0 
)

Performs the calculation of axion-photon absorption probability using directly equation (18) from van Bibber, Phys Rev D Part Fields. 1989.

If not provided as argument, m_gamma it will be attempted to obtain the value from the buffer gas definition. If no buffer gas has been assigned the medium will be assumed to be vacuum.

Ea in keV, ma in eV, mgamma in eV, Lcoh in mm, Bmag in T

mg in eV, absLength in cm-1

The returned value is given for g_ag = 10^-10 GeV-1

Definition at line 682 of file TRestAxionField.cxx.

◆ BL()

double TRestAxionField::BL ( Double_t  Bmag,
Double_t  Lcoh 
)

Performs the calculation of (BL) factor in natural units.

Lcoh should be expressed in mm, and Bmag in T. The result will be given for an axion-photon coupling of 10^{-10} GeV^{-1}

Definition at line 251 of file TRestAxionField.cxx.

◆ BLHalfSquared()

double TRestAxionField::BLHalfSquared ( Double_t  Bmag,
Double_t  Lcoh 
)

Performs the calculation of (BL/2)^2 factor in natural units.

Lcoh should be expressed in mm, and Bmag in T. The result will be given for an axion-photon coupling of 10^{-10} GeV^{-1}

Definition at line 268 of file TRestAxionField.cxx.

◆ ComputeOffResonanceIntegral()

std::pair< Double_t, Double_t > TRestAxionField::ComputeOffResonanceIntegral ( Double_t  q,
Double_t  Gamma,
Double_t  accuracy,
Int_t  num_intervals,
Int_t  qawo_levels 
)
private

Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf.

It integrates the Integrand function defined on the header of TRestAxionField.

This method uses the GSL QAWO method for oscillatory functions. That is the case when the q is not zero, in the off-resonance case.

See https://www.gnu.org/software/gsl/doc/html/integration.html for more details on the integration parameters.

The Gamma function should be expressed in mm-1 since the field map accessed in the integrand is evaluated using mm.

Definition at line 609 of file TRestAxionField.cxx.

◆ ComputeResonanceIntegral()

std::pair< Double_t, Double_t > TRestAxionField::ComputeResonanceIntegral ( Double_t  Gamma,
Double_t  accuracy,
Int_t  num_intervals 
)
private

Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf.

It integrates the Integrand function defined on the header of TRestAxionField.

This method uses the GSL QAG integration method. We use this method when the cosine function vanishes when q = 0.

See https://www.gnu.org/software/gsl/doc/html/integration.html for more details on the integration parameters.

The Gamma function should be expressed in mm-1 since the field map accessed in the integrand is evaluated using mm.

Definition at line 549 of file TRestAxionField.cxx.

◆ GammaTransmissionFieldMapProbability()

std::pair< Double_t, Double_t > TRestAxionField::GammaTransmissionFieldMapProbability ( Double_t  Ea,
Double_t  ma,
Double_t  accuracy = 1.e-1,
Int_t  num_intervals = 100,
Int_t  qawo_levels = 20 
)

Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf.

m_gamma will be obtainned from the buffer gas definition. If no buffer gas has been assigned then the medium will be assumed to be vacuum.

Ea in keV, ma in eV

The mgamma and absorption lengths are the ones defined by the gas.

The returned value is given for g_ag = 10^-10 GeV-1

Note
The density is for the moment homogeneous. We would need to implemnent a double integral to solve the problem with a density profile.

Definition at line 482 of file TRestAxionField.cxx.

◆ GammaTransmissionFWHM()

Double_t TRestAxionField::GammaTransmissionFWHM ( Double_t  step = 0.00001)

Performs the calculation of the FWHM for the axion-photon conversion probability computed in TRestAxionField::GammaTransmissionProbability.

Ea in keV, ma in eV, mgamma in eV, deltaL in mm, Bmag in T.

Default value for the parameters are: Bmag = 2.5 T, Lcoh = 10000 mm, Ea = 4 keV eV,

IMPORTANT: In the case that the buffer gas is not defined, this method will return the mass at which the probability reaches half of the maximum vacuum probability.

Scanning towards the right (valid also for vacuum)

Definition at line 764 of file TRestAxionField.cxx.

◆ GammaTransmissionProbability() [1/3]

Double_t TRestAxionField::GammaTransmissionProbability ( Double_t  Bmag,
Double_t  Lcoh,
Double_t  Ea,
Double_t  ma,
Double_t  mg = 0,
Double_t  absLength = 0 
)

On top of calculating the gamma transmission probability it will assign new values for the magnetic field (Bmag/T), coherence length (Lcoh/mm) and axion energy (Ea/keV).

Definition at line 353 of file TRestAxionField.cxx.

◆ GammaTransmissionProbability() [2/3]

Double_t TRestAxionField::GammaTransmissionProbability ( Double_t  ma,
Double_t  mg = 0,
Double_t  absLength = 0 
)

Performs the calculation of axion-photon conversion probability using directly equation (11) from van Bibber, Phys Rev D Part Fields. 1989.

If m_gamma (mg) is not given as an argument, i.e. it is equal to zero, then m_gamma will be obtainned from the buffer gas definition. If no buffer gas has been assigned then the medium will be assumed to be vacuum.

Ea in keV, ma in eV, mgamma in eV, Lcoh in mm, Bmag in T

mg in eV, absLength in cm-1

The returned value is given for g_ag = 10^-10 GeV-1

Definition at line 294 of file TRestAxionField.cxx.

◆ GammaTransmissionProbability() [3/3]

Double_t TRestAxionField::GammaTransmissionProbability ( std::vector< Double_t >  Bmag,
Double_t  deltaL,
Double_t  Ea,
Double_t  ma,
Double_t  mg = 0,
Double_t  absLength = 0 
)

Performs the calculation of axion-photon conversion probability using directly equation (28) from J. Redondo and A. Ringwald, Light shinning through walls. https://arxiv.org/pdf/1011.3741.pdf.

If m_gamma (mg) is not given as an argument, i.e. it is equal to zero, then m_gamma will be obtainned from the buffer gas definition. If no buffer gas has been assigned then the medium will be assumed to be vacuum.

Ea in keV, ma in eV, mgamma in eV, deltaL in mm, Bmag in T

mg in eV, absLength in cm-1

The returned value is given for g_ag = 10^-10 GeV-1

Note
The density is for the moment homogeneous. We would need to implemnent a double integral to solve the problem with a density profile. TOBE implemented in a new method if needed, where Gamma is not constant and \integral{q(z)} is integrated at each step.

We integrate following the Midpoint rule method. (Other potential options : Trapezoidal, Simpsons)

Definition at line 381 of file TRestAxionField.cxx.

◆ GetAxionEnergy()

Double_t TRestAxionField::GetAxionEnergy ( ) const
inline

Definition at line 64 of file TRestAxionField.h.

◆ GetCoherenceLength()

Double_t TRestAxionField::GetCoherenceLength ( ) const
inline

Definition at line 63 of file TRestAxionField.h.

◆ GetMagneticField()

Double_t TRestAxionField::GetMagneticField ( ) const
inline

Definition at line 62 of file TRestAxionField.h.

◆ GetMassDensityScanning()

std::vector< std::pair< Double_t, Double_t > > TRestAxionField::GetMassDensityScanning ( std::string  gasName = "He",
Double_t  maMax = 0.15,
Double_t  rampDown = 5.0 
)

This method determines the proper densities to be used in an axion helioscope experiment in order to achieve a continuous axion mass scan.

The first scanning density is placed where the axion-photon vacuum probability reaches half the value, P_ag(max)/2. Once the first density, or step, has been obtained, the method calculates the FWHM resonance probability for each density/mass and moves the next scanning axion mass by a step of amplitude FWHM/factor, where the factor is a value that moves from 2 to 1 as the mass increases, and it falls down with a velocity related with the argument rampDown, where the factor follows this formula:

factor = TMath::Exp( - ma*rampDown ) + 1;

The method stops when the axion mass is bigger than the maximum axion mass given as an argument, maMax.

Default arguments: gasName="He", ma_max=0.15 eV, rampDown=5

Returns
It returns a vector of pair with the values for the scan, the first one is the axion mass and the second one is the density.

For additional info see PR: https://github.com/rest-for-physics/axionlib/pull/78

Setting mass-density pair for the first step

Definition at line 814 of file TRestAxionField.cxx.

◆ Initialize()

void TRestAxionField::Initialize ( )
private

Initialization of TRestAxionField class.

MOVED TO TRestAxionFieldPropagationProcess class faxion = SetComplexReal(1, 0); fAem = SetComplexReal(0, 0);

Definition at line 237 of file TRestAxionField.cxx.

◆ Integrand()

static double TRestAxionField::Integrand ( double  x,
void *  params 
)
inlinestatic

Integrand used for axion-photon probability integration.

Definition at line 100 of file TRestAxionField.h.

◆ SetAxionEnergy()

void TRestAxionField::SetAxionEnergy ( Double_t  e)
inline

Definition at line 60 of file TRestAxionField.h.

◆ SetBufferGas()

void TRestAxionField::SetBufferGas ( TRestAxionBufferGas buffGas)
inline

It assigns a gas buffer medium to the calculation.

Definition at line 79 of file TRestAxionField.h.

◆ SetCoherenceLength()

void TRestAxionField::SetCoherenceLength ( Double_t  l)
inline

Definition at line 59 of file TRestAxionField.h.

◆ SetDebug()

void TRestAxionField::SetDebug ( Bool_t  v)
inline

It enables/disables debug mode.

Definition at line 70 of file TRestAxionField.h.

◆ SetMagneticField()

void TRestAxionField::SetMagneticField ( Double_t  b)
inline

Definition at line 58 of file TRestAxionField.h.

Field Documentation

◆ fBmag

Double_t TRestAxionField::fBmag = 2.5
private

The magnetic field in Teslas (used for constant field formulas)

Definition at line 35 of file TRestAxionField.h.

◆ fBufferGas

TRestAxionBufferGas* TRestAxionField::fBufferGas = nullptr
private

A pointer to the buffer gas definition.

Definition at line 46 of file TRestAxionField.h.

◆ fDebug

Bool_t TRestAxionField::fDebug = false
private

Definition at line 32 of file TRestAxionField.h.

◆ fEa

Double_t TRestAxionField::fEa = 4.2
private

The energy of the axion in keV.

Definition at line 41 of file TRestAxionField.h.

◆ fLcoh

Double_t TRestAxionField::fLcoh = 10000
private

The coherence lenght (in mm) where the magnetic field is defined (for constant field)

Definition at line 38 of file TRestAxionField.h.

◆ fMagneticField

TRestAxionMagneticField* TRestAxionField::fMagneticField = nullptr
private

A pointer to the magnetic field definition.

Definition at line 49 of file TRestAxionField.h.


The documentation for this class was generated from the following files: