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

Detailed Description

A metadata class to load tabulated solar axion fluxes. Mass independent.

TRestAxionSolarQCDFlux will use a file in ASCII or binary format to initialize a solar flux table that will describe the solar flux spectrum as a function of the solar radius.

This class may serve to load any generic flux definition that is independent from the axion mass. However, since the design of the class was motivated to reproduce the standard QCD axion flux, therefore the name of the class.

Another scenario arises when the axion or an axion-like particle production mechanism depends on its mass, and we may need to introduce a flux description as the given in the class TRestAxionSolarHiddenPhotonFlux. Both classes are prototyped by a pure base class TRestAxionSolarFlux that defines common methods used to evaluate the flux, and generate Monte-Carlo events inside TRestAxionGeneratorProcess.

Basic use

Once the class has been initialized, the main use of this class will be provided by the method TRestAxionSolarQCDFlux::GetRandomEnergyAndRadius. This method will return a random axion energy and position inside the solar radius following the distributions given by the solar flux tables.

Description of the specific parameters accepted by this metadata class.

Additionally this class will be able to read .flux files that are the original files produced in 3-columns format (inner radius [solar units] / energy [keV] / flux [cm-2 s-1 keV-1]). The .flux files may contain the full information, continuum and spectral components. Those components will be splited into two independent contributions by TRestAxionSolarQCDFlux::ReadFluxFile to be managed internally. Two additional parameters will be required to translate the .flux files into the tables that are understood by this class.

Optionally, if we want to consider a different binning on the monochromatic/continuum histogram used internally for the calculation we may specify optionally a new parameter. In that case, fBinSize will be the binning of the internal histogram, while the new parameter will be the binning given inside the .flux file.

Pre-generated solar axion flux tables will be available at the axionlib-data repository. The different RML flux definitions used to load those tables will be found at the fluxes.rml file found at the axionlib-data repository.

Inside a local REST installation, the fluxes.rml file will be found at the REST installation directory, and it will be located automatically by the TRestMetadata::SearchFile method.

A basic RML definition

The following definition integrates an axion-photon component with a continuum spectrum using a Primakoff production model, and a dummy spectrum file that includes two monocrhomatic lines at different solar disk radius positions.

<TRestAxionSolarQCDFlux name="sunPrimakoff" verboseLevel="debug" >
<!-- Common parameters that belong to TRestAxionSolarFlux -->
<parameter name="couplingType" value="g_ag"/>
<parameter name="couplingStrength" value="1.e-10"/>
<!-- Specific parameters that belong to TRestAxionSolarQCDFlux -->
<parameter name="fluxDataFile" value="Primakoff_Gianotti_201904.dat"/>
<parameter name="fluxSptFile" value="Dummy_Galan_202202.spt"/>
A metadata class to load tabulated solar axion fluxes.
A metadata class to load tabulated solar axion fluxes. Mass independent.
Warning
When the flux is loaded manually inside the restRoot interactive shell, or inside a macro or script, after metadata initialization, it is necessary to call the method TRestAxionSolarQCDFlux::LoadTables to trigger the tables initialization.

Performing MonteCarlo tests using pre-loaded tables

In order to test the response of different solar flux definitions we may use the script solarPlot.py found at pipeline/metadata/solarFlux/. This script will generate a number of particles and it will assign to each particle an energy and solar disk location with the help of the method TRestAxionSolarQCDFlux::GetRandomEnergyAndRadius.

python3 solarPlotQCD.py --fluxname LennertHoofABC --N 1000000

By default, it will load the flux definition found at fluxes.rml from the axionlib-data repository, and generate a png image with the resuts from the Monte Carlo execution.

Reading solar flux tables from .flux files

In a similar way we may initialize the class using a .flux file. The .flux described previously will contain a high definition flux table measured in cm-2 s-1 keV-1 as a function of the energy (from 0 to 20keV) and the inner solar radius (from 0 to 1). The binning of this table will be typically between 1eV and 10eV.

The class TRestAxionSolarQCDFlux will be initialized directly using the .flux file provided under the fluxDataFile parameter. During the initialization, the flux will be splitted into two independent flux components. One smooth continuum component integrated in 100 eV steps, and a monochromatic peak components.

In order to help with the identification of peaks we need to define the binSize used in the .flux table and the peakSigma defining the number of sigmas over the average for a bin to be considered a peak.

<TRestAxionSolarQCDFlux name="LennertHoofABC_Flux" verboseLevel="warning" >
<parameter name="couplingType" value="g_ae"/>
<parameter name="couplingStrength" value="1.e-13"/>
<parameter name="fluxDataFile" value="ABC_LennertHoof_202203.flux"/>
<parameter name="binSize" value="10eV" />
<parameter name="peakSigma" value="10" />
<parameter name="seed" value="137" />

We will be able to load this file as usual, using the following recipe inside restRoot,

TRestAxionSolarQCDFlux *sFlux = new TRestAxionSolarQCDFlux("fluxes.rml", "LennertHoofABC")
sFlux->Initialize()
TCanvas *c = sFlux->DrawSolarFluxes()
c->Print("ABC_FluxTable.png" )
void Initialize()
It is required in order to load solar flux tables into memory.
TRestAxionSolarQCDFlux()
Default constructor.
void Print()
Implementing TObject::Print() method.

will generate the following figure.

Exporting the solar flux tables

On top of that, we will be able to export those tables to the TRestAxionSolarQCDFlux standard format to be used in later occasions.

TRestAxionSolarQCDFlux *sFlux = new TRestAxionSolarQCDFlux("fluxes.rml", "LennertHoofABC")
sFlux->Initialize()
sFlux->ExportTables()
void ExportTables(Bool_t ascii=false) override
It will create files with the continuum and spectral flux components to be used in a later ocasion.

which will produce two files, a binary table .N200f with the continuum flux, and an ASCII table containning the .spt monochromatic lines. The filename root will be extracted from the original .flux file. Optionally we may export the continuum flux to an ASCII file by indicating it at the TRestAxionSolarQCDFlux::ExportTables method call. The files will be placed at the REST user space, at $HOME/.rest/export/ directory.

TODO Implement the method TRestAxionSolarQCDFlux::InitializeSolarTable using a solar model description by TRestAxionSolarModel.

TODO Perhaps it would be interesting to replace fFluxTable for a TH2D


RESTsoft - Software for Rare Event Searches with TPCs

History of developments:

2023-May: Specific methods extracted from TRestAxionSolarFlux Javier Galan

Author
Javier Galan

Definition at line 30 of file TRestAxionSolarQCDFlux.h.

#include <TRestAxionSolarQCDFlux.h>

Inheritance diagram for TRestAxionSolarQCDFlux:
TRestAxionSolarFlux TRestMetadata

Public Member Functions

 ClassDefOverride (TRestAxionSolarQCDFlux, 1)
 
virtual TCanvas * DrawSolarFlux () override
 It draws the contents of a .flux file. This method just receives the.
 
void ExportTables (Bool_t ascii=false) override
 It will create files with the continuum and spectral flux components to be used in a later ocasion.
 
TH1F * GetContinuumSpectrum ()
 It builds a histogram with the continuum spectrum component. The flux will be expressed in cm-2 s-1 keV-1. Binned in 100eV steps.
 
TH1F * GetEnergySpectrum (Double_t m=0) override
 It returns an energy integrated spectrum in cm-2 s-1 keV-1.
 
TH1F * GetMonochromaticSpectrum ()
 It builds a histogram with the monochromatic spectrum component. The flux will be expressed in cm-2 s-1 eV-1. Binned in 1eV steps.
 
std::pair< Double_t, Double_t > GetRandomEnergyAndRadius (TVector2 eRange=TVector2(-1, -1), Double_t mass=0) override
 It defines how to generate Monte Carlo energy and radius values to reproduce the flux. More...
 
Double_t GetTotalFlux (Double_t mass=0) override
 It returns the total integrated flux at earth in cm-2 s-1.
 
TH1F * GetTotalSpectrum ()
 It builds a histogram adding the continuum and the monochromatic spectrum component. The flux will be expressed in cm-2 s-1 keV-1. Binned in 1eV steps.
 
void InitializeSolarTable (TRestAxionSolarModel *model)
 Tables might be loaded using a solar model description by TRestAxionSolarModel.
 
Double_t IntegrateFluxInRange (TVector2 eRange=TVector2(-1, -1), Double_t mass=0) override
 It returns the integrated flux at earth in cm-2 s-1 for the given energy range.
 
Bool_t isSolarSpectrumLoaded ()
 It returns true if monochromatic flux spectra was loaded.
 
Bool_t isSolarTableLoaded ()
 It returns true if continuum flux spectra was loaded.
 
Bool_t LoadTables () override
 It defines how to read the solar tables at the inhereted class. More...
 
void PrintContinuumSolarTable ()
 It prints on screen the table that has been loaded in memory.
 
void PrintIntegratedRingFlux ()
 It prints on screen the integrated solar flux per solar ring.
 
void PrintMetadata () override
 Prints on screen the information about the metadata members of TRestAxionSolarQCDFlux.
 
void PrintMonoChromaticFlux ()
 It prints on screen the spectral lines loaded in memory.
 
 TRestAxionSolarQCDFlux ()
 Default constructor.
 
 TRestAxionSolarQCDFlux (const char *cfgFileName, std::string name="")
 Constructor loading data from a config file. More...
 
 ~TRestAxionSolarQCDFlux ()
 Default destructor.
 
- Public Member Functions inherited from TRestAxionSolarFlux
Bool_t AreTablesLoaded ()
 
TCanvas * DrawFluxFile (std::string fname, Double_t binSize=0.01)
 It draws the contents of a .flux file. This method just receives the name of the .flux file and it works stand-alone.
 
TH1F * GetFluxHistogram (std::string fname, Double_t binSize=0.01)
 It builds a histogram using the contents of the .flux file given in the argument.
 
void Initialize ()
 It is required in order to load solar flux tables into memory. More...
 
 ~TRestAxionSolarFlux ()
 Default destructor.
 

Private Member Functions

void IntegrateSolarFluxes ()
 A helper method to initialize the internal class data members with the integrated flux for each solar ring. It will be called by TRestAxionSolarQCDFlux::Initialize.
 
void LoadContinuumFluxTable ()
 A helper method to load the data file containning continuum spectra as a function of the solar radius. It will be called by TRestAxionSolarQCDFlux::Initialize.
 
void LoadMonoChromaticFluxTable ()
 A helper method to load the data file containning monochromatic spectral lines as a function of the solar radius. It will be called by TRestAxionSolarQCDFlux::Initialize.
 
void ReadFluxFile ()
 It loads a .flux file. It will split continuum and monochromatic peaks, loading both internal flux tables.
 

Private Attributes

Double_t fBinSize = 0
 It will be used when loading .flux files to define the input file energy binsize in eV.
 
TH1F * fContinuumHist = nullptr
 A pointer to the continuum spectrum histogram.
 
std::string fFluxDataFile = ""
 The filename containning the solar flux table with continuum spectrum.
 
std::vector< Double_t > fFluxLineIntegrals
 Accumulative integrated solar flux for each monochromatic energy (renormalized to unity)
 
std::map< Double_t, TH1F * > fFluxLines
 The tabulated solar flux in cm-2 s-1 for a number of monochromatic energies versus solar radius.
 
Double_t fFluxRatio = 0
 The ratio between monochromatic and total flux.
 
std::string fFluxSptFile = ""
 The filename containning the solar flux spectra for monochromatic spectrum.
 
std::vector< TH1F * > fFluxTable
 The tabulated solar flux continuum spectra TH1F(200,0,20)keV in cm-2 s-1 keV-1 versus solar radius.
 
std::vector< Double_t > fFluxTableIntegrals
 Accumulative integrated solar flux for each solar ring for continuum spectrum (renormalized to unity)
 
TH1F * fMonoHist = nullptr
 A pointer to the monochromatic spectrum histogram.
 
Double_t fPeakSigma = 0
 It will be used when loading .flux files to define the threshold for peak identification.
 
Double_t fTotalContinuumFlux = 0
 Total solar flux for monochromatic contributions.
 
TH1F * fTotalHist = nullptr
 A pointer to the superposed monochromatic and continuum spectrum histogram.
 
Double_t fTotalMonochromaticFlux = 0
 Total solar flux for monochromatic contributions.
 

Additional Inherited Members

- Protected Member Functions inherited from TRestAxionSolarFlux
 TRestAxionSolarFlux ()
 Default constructor.
 
 TRestAxionSolarFlux (const char *cfgFileName, std::string name="")
 Constructor loading data from a config file. More...
 
- Protected Attributes inherited from TRestAxionSolarFlux
TCanvas * fCanvas = nullptr
 A canvas pointer for drawing.
 
TRandom3 * fRandom = nullptr
 Random number generator.
 

Constructor & Destructor Documentation

◆ TRestAxionSolarQCDFlux()

TRestAxionSolarQCDFlux::TRestAxionSolarQCDFlux ( const char *  cfgFileName,
std::string  name = "" 
)

Constructor loading data from a config file.

If no configuration path is defined using TRestMetadata::SetConfigFilePath the path to the config file must be specified using full path, absolute or relative.

The default behaviour is that the config file must be specified with full path, absolute or relative.

Parameters
cfgFileNameA const char* giving the path to an RML file.
nameThe name of the specific metadata. It will be used to find the corresponding TRestAxionMagneticField section inside the RML.

Definition at line 236 of file TRestAxionSolarQCDFlux.cxx.

Member Function Documentation

◆ GetRandomEnergyAndRadius()

std::pair< Double_t, Double_t > TRestAxionSolarQCDFlux::GetRandomEnergyAndRadius ( TVector2  eRange = TVector2(-1, -1),
Double_t  mass = 0 
)
overridevirtual

It defines how to generate Monte Carlo energy and radius values to reproduce the flux.

It returns a random solar radius position and energy according to the flux distributions defined inside the solar tables loaded in the class.

Implements TRestAxionSolarFlux.

Definition at line 661 of file TRestAxionSolarQCDFlux.cxx.

◆ LoadTables()

Bool_t TRestAxionSolarQCDFlux::LoadTables ( )
overridevirtual

It defines how to read the solar tables at the inhereted class.

It will load the tables in memory by using the filename information provided inside the metadata members.

Implements TRestAxionSolarFlux.

Definition at line 252 of file TRestAxionSolarQCDFlux.cxx.


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