REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionSolarQCDFlux.h
1 /*************************************************************************
2  * This file is part of the REST software framework. *
3  * *
4  * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5  * For more information see http://gifna.unizar.es/trex *
6  * *
7  * REST is free software: you can redistribute it and/or modify *
8  * it under the terms of the GNU General Public License as published by *
9  * the Free Software Foundation, either version 3 of the License, or *
10  * (at your option) any later version. *
11  * *
12  * REST is distributed in the hope that it will be useful, *
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15  * GNU General Public License for more details. *
16  * *
17  * You should have a copy of the GNU General Public License along with *
18  * REST in $REST_PATH/LICENSE. *
19  * If not, see http://www.gnu.org/licenses/. *
20  * For the list of contributors see $REST_PATH/CREDITS. *
21  *************************************************************************/
22 
23 #ifndef _TRestAxionSolarQCDFlux
24 #define _TRestAxionSolarQCDFlux
25 
26 #include <TRestAxionSolarFlux.h>
27 #include <TRestAxionSolarModel.h>
28 
31  private:
33  std::string fFluxDataFile = ""; //<
34 
36  std::string fFluxSptFile = ""; //<
37 
39  Double_t fBinSize = 0; //<
40 
42  Double_t fPeakSigma = 0; //<
43 
45  std::vector<TH1F*> fFluxTable;
46 
48  std::map<Double_t, TH1F*> fFluxLines;
49 
51  std::vector<Double_t> fFluxTableIntegrals;
52 
54  std::vector<Double_t> fFluxLineIntegrals;
55 
57  Double_t fTotalMonochromaticFlux = 0;
58 
60  Double_t fTotalContinuumFlux = 0;
61 
63  Double_t fFluxRatio = 0;
64 
66  TH1F* fContinuumHist = nullptr;
67 
69  TH1F* fMonoHist = nullptr;
70 
72  TH1F* fTotalHist = nullptr;
73 
74  void ReadFluxFile();
77  void IntegrateSolarFluxes();
78 
79  public:
81  Bool_t isSolarTableLoaded() { return fFluxTable.size() > 0; }
82 
84  Bool_t isSolarSpectrumLoaded() { return fFluxLines.size() > 0; }
85 
87  Double_t IntegrateFluxInRange(TVector2 eRange = TVector2(-1, -1), Double_t mass = 0) override;
88 
90  std::pair<Double_t, Double_t> GetRandomEnergyAndRadius(TVector2 eRange = TVector2(-1, -1),
91  Double_t mass = 0) override;
92 
94  Bool_t LoadTables() override;
95 
97  Double_t GetTotalFlux(Double_t mass = 0) override {
99  }
100 
102  TH1F* GetEnergySpectrum(Double_t m = 0) override { return GetTotalSpectrum(); }
103 
104  TH1F* GetContinuumSpectrum();
105  TH1F* GetMonochromaticSpectrum();
106  TH1F* GetTotalSpectrum();
107 
108  virtual TCanvas* DrawSolarFlux() override;
109 
112  // TOBE implemented
113  // This method should initialize the tables fFluxTable and fFluxLines
114  }
115 
116  void ExportTables(Bool_t ascii = false) override;
117 
118  void PrintMetadata() override;
119 
122  void PrintMonoChromaticFlux();
123 
125  TRestAxionSolarQCDFlux(const char* cfgFileName, std::string name = "");
127 
128  ClassDefOverride(TRestAxionSolarQCDFlux, 1);
129 };
130 #endif
A metadata class to load tabulated solar axion fluxes.
A metadata class to define theoretical axion models and calculations related.
A metadata class to load tabulated solar axion fluxes. Mass independent.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarQCDFlux.
void InitializeSolarTable(TRestAxionSolarModel *model)
Tables might be loaded using a solar model description by TRestAxionSolarModel.
std::string fFluxSptFile
The filename containning the solar flux spectra for monochromatic spectrum.
Double_t fPeakSigma
It will be used when loading .flux files to define the threshold for peak identification.
TRestAxionSolarQCDFlux()
Default constructor.
~TRestAxionSolarQCDFlux()
Default destructor.
void ReadFluxFile()
It loads a .flux file. It will split continuum and monochromatic peaks, loading both internal flux ta...
TH1F * GetTotalSpectrum()
It builds a histogram adding the continuum and the monochromatic spectrum component....
TH1F * fMonoHist
A pointer to the monochromatic spectrum histogram.
TH1F * fTotalHist
A pointer to the superposed monochromatic and continuum spectrum histogram.
void LoadContinuumFluxTable()
A helper method to load the data file containning continuum spectra as a function of the solar radius...
void PrintIntegratedRingFlux()
It prints on screen the integrated solar flux per solar ring.
std::vector< Double_t > fFluxTableIntegrals
Accumulative integrated solar flux for each solar ring for continuum spectrum (renormalized to unity)
Double_t GetTotalFlux(Double_t mass=0) override
It returns the total integrated flux at earth in cm-2 s-1.
Bool_t isSolarTableLoaded()
It returns true if continuum flux spectra was loaded.
Bool_t isSolarSpectrumLoaded()
It returns true if monochromatic flux spectra was loaded.
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 LoadTables() override
It defines how to read the solar tables at the inhereted class.
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.
TH1F * GetContinuumSpectrum()
It builds a histogram with the continuum spectrum component. The flux will be expressed in cm-2 s-1 k...
void IntegrateSolarFluxes()
A helper method to initialize the internal class data members with the integrated flux for each solar...
TH1F * fContinuumHist
A pointer to the continuum spectrum histogram.
void LoadMonoChromaticFluxTable()
A helper method to load the data file containning monochromatic spectral lines as a function of the s...
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.
std::vector< Double_t > fFluxLineIntegrals
Accumulative integrated solar flux for each monochromatic energy (renormalized to unity)
TH1F * GetMonochromaticSpectrum()
It builds a histogram with the monochromatic spectrum component. The flux will be expressed in cm-2 s...
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.
std::map< Double_t, TH1F * > fFluxLines
The tabulated solar flux in cm-2 s-1 for a number of monochromatic energies versus solar radius.
virtual TCanvas * DrawSolarFlux() override
It draws the contents of a .flux file. This method just receives the.
void PrintMonoChromaticFlux()
It prints on screen the spectral lines loaded in memory.
std::string fFluxDataFile
The filename containning the solar flux table with continuum spectrum.
Double_t fFluxRatio
The ratio between monochromatic and total flux.
TH1F * GetEnergySpectrum(Double_t m=0) override
It returns an energy integrated spectrum in cm-2 s-1 keV-1.
void PrintContinuumSolarTable()
It prints on screen the table that has been loaded in memory.
Double_t fBinSize
It will be used when loading .flux files to define the input file energy binsize in eV.
Double_t fTotalContinuumFlux
Total solar flux for monochromatic contributions.
Double_t fTotalMonochromaticFlux
Total solar flux for monochromatic contributions.