REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionSolarFlux.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 _TRestAxionSolarFlux
24 #define _TRestAxionSolarFlux
25 
26 #include <TCanvas.h>
27 #include <TH1F.h>
28 #include <TH2F.h>
29 #include <TRandom3.h>
30 #include <TRestMetadata.h>
31 
34  private:
36  std::string fCouplingType = ""; //<
37 
39  Double_t fCouplingStrength = 0; //<
40 
42  Int_t fSeed = 0; //<
43 
45  Bool_t fTablesLoaded = false;
46 
47  protected:
49  TCanvas* fCanvas = nullptr;
50 
52  TRandom3* fRandom = nullptr;
53 
55  TRestAxionSolarFlux(const char* cfgFileName, std::string name = "");
56 
58  virtual Bool_t LoadTables() = 0;
59 
60  public:
62  void Initialize();
63 
65  virtual Double_t IntegrateFluxInRange(TVector2 eRange = TVector2(-1, -1), Double_t mass = 0) = 0;
66 
68  virtual Double_t GetTotalFlux(Double_t mass = 0) = 0;
69 
71  virtual std::pair<Double_t, Double_t> GetRandomEnergyAndRadius(TVector2 eRange = TVector2(-1, -1),
72  Double_t mass = 0) = 0;
73 
75  virtual TH1F* GetEnergySpectrum(Double_t m = 0) = 0;
76 
77  virtual TCanvas* DrawSolarFlux();
78 
79  virtual void ExportTables(Bool_t ascii = false) {
80  RESTWarning << "TRestAxionSolarFlux::ExportTables must be re-implemented in the inherited class"
81  << RESTendl;
82  }
83 
84  Bool_t AreTablesLoaded() { return fTablesLoaded; }
85 
86  TH1F* GetFluxHistogram(std::string fname, Double_t binSize = 0.01);
87  TCanvas* DrawFluxFile(std::string fname, Double_t binSize = 0.01);
88 
89  void PrintMetadata();
90 
92 
93  ClassDef(TRestAxionSolarFlux, 2);
94 };
95 #endif
A metadata class to load tabulated solar axion fluxes.
void Initialize()
It is required in order to load solar flux tables into memory.
TRestAxionSolarFlux()
Default constructor.
TCanvas * fCanvas
A canvas pointer for drawing.
~TRestAxionSolarFlux()
Default destructor.
virtual Bool_t LoadTables()=0
It defines how to read the solar tables at the inhereted class.
Double_t fCouplingStrength
Axion coupling strength.
virtual std::pair< Double_t, Double_t > GetRandomEnergyAndRadius(TVector2 eRange=TVector2(-1, -1), Double_t mass=0)=0
It defines how to generate Monte Carlo energy and radius values to reproduce the flux.
virtual Double_t GetTotalFlux(Double_t mass=0)=0
It returns the total integrated flux at earth in cm-2 s-1.
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 ....
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.
std::string fCouplingType
Axion coupling. Defines coupling type and strength.
Int_t fSeed
Seed used in random generator.
TRandom3 * fRandom
Random number generator.
virtual TH1F * GetEnergySpectrum(Double_t m=0)=0
It returns an energy integrated spectrum in cm-2 s-1 keV-1.
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
virtual Double_t IntegrateFluxInRange(TVector2 eRange=TVector2(-1, -1), Double_t mass=0)=0
It returns the integrated flux at earth in cm-2 s-1 for the given energy range.
Bool_t fTablesLoaded
A metadata member to control if this class has been initialized.
virtual TCanvas * DrawSolarFlux()
It draws the contents of a .flux file. This method just receives the.
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
endl_t RESTendl
Termination flag object for TRestStringOutput.