Berkeley Nuclear Data Software
GammaEfficiencyAnalysis.h
Go to the documentation of this file.
1 #include "SpectrumAnalysis.h"
2 #include "TFile.h"
3 #include "TTree.h"
4 #include "TChain.h"
5 #include "TH2.h"
6 #include "TCanvas.h"
7 #include <vector>
8 #include <string>
9 #include <map>
10 #include <random>
11 using std::vector;
12 using std::string;
13 using std::map;
14 using std::pair;
15 #include <Eigen/Dense>
20 {
21 public:
23  void setVerbosity(int a_verbosity);
25  void loadExperimentalSpectra(string a_fileName,
26  string a_isotope,
27  double a_countTime = 1,
28  double a_activity = 1,
29  int a_numLines = 5);
30  void useTFIntDiffCorrection(bool a_useCorr);
31  void loadBackgroundSpectra(string a_fileName,
32  double a_countTime);
35  int addSimulatedData(string a_fileName, int a_numSimEvents);
37  void setSimHistNumBins(int a_bins);
40  void setSimHistMaxE(double a_maxEnergy);
44  void processSimulatedData(string a_configFileName,
45  bool a_resampleEffCorr=false,
46  double a_acceptanceWindow=2.5e-3);
48  TH1F* getSimHist(DetType a_detType,
49  int a_detID,
50  int a_elemID=0,
51  bool a_addBack=false);
52  double getSimEff(double a_gammaEnergy,
53  DetType a_detType,
54  int a_detID,
55  int a_elemID=0,
56  bool a_addBack=false);
59  void fitTriggerProb(string a_fileNameHighTF,
60  string a_fileNameLowTF,
61  double a_countTimeHigh,
62  double a_countTimeLow,
63  string a_fileNameBKG="",
64  double a_countTimeBKG=-1,
65  string a_isotope="Eu152",
66  int a_numLines = 16);
68  double calcTriggerProb(DetType a_detType,
69  int a_detID,
70  int a_elemID,
71  double a_energy);
77  void computeAllExpEffs(bool a_drawFits=false);
78  void drawEffPeakFits(DetType a_detType,
79  int a_detID,
80  int a_elemID=0,
81  bool a_addBack=false);
82  vector<vector<double>> getExpEffs(DetType a_detType,
83  int a_detID,
84  int a_elemID=0,
85  bool a_addBack=false);
88  map<int, vector<vector<double>>> calcAddBackFactors();
90  void plotAddBackSimExpComp();
92  void plotSimExpComp(DetType a_detType,
93  int a_detID,
94  int a_elemID=0,
95  bool a_addBack=false);
96  map<double, PeakInfo> fitIsotopePeaks(TH1* a_hist,
97  string a_isotope,
98  double a_activity,
99  double a_countTime,
100  int a_numLines,
101  bool a_draw=false
102  );
103  map<double, std::pair<double,double>> calcEfficiency(TH1* a_hist,
104  string a_isotope,
105  double a_activity,
106  double a_countTime,
107  int a_numLines,
108  vector<PeakInfo>& a_peakFits
109  );
110  map<double, std::pair<double,double>> calcEfficiency(
111  string a_isotope,
112  double a_activity,
113  double a_countTime,
114  int a_numLines,
115  vector<PeakInfo>& a_peakFits
116  );
117  map<double, std::pair<double,double>> calcEfficiency(map<double, PeakInfo> a_peakFits,
118  string a_isotope,
119  double a_activity,
120  double a_countTime,
121  int a_numLines
122  );
123  void setNewPeakFit(PeakInfo a_peakInfo,
124  DetType a_detType,
125  int a_detID,
126  int a_elemID=0,
127  bool a_addBack=false);
128  vector<PeakInfo> getPeakFits(DetType a_detType,
129  int a_detID,
130  int a_elemID=0,
131  bool a_addBack=false);
132  map<double, PeakInfo> m_lastPeakFits;
133  vector<double> m_gammaEnergies;
134  map<DetType, map<int, map<int, Eigen::VectorXf>>> m_corrParamsD0;
135  map<DetType, map<int, map<int, Eigen::MatrixXf>>> m_corrParamsDCov;
136  map<DetType, map<int, map<int, Eigen::VectorXf>>> m_corrParamsD;
137  map<DetType, map<int, map<int, TH2F*>>> m_angEffHists;
138 private:
139  SpectrumAnalysis m_specAna;
140  int m_verbosity;
141  TChain* m_simTree;
142  bool m_hasSimData;
143  int m_numSimEvents;
144  int m_numSimEnergyBins;
145  double m_maxSimEnergy;
146  int m_numResamples;
148  map<string, map<string, vector<double>>> m_expData;
150  // map<DetType, map<int, map<int, Eigen::VectorXf>>> m_corrParamsD;
151  // map<DetType, map<int, map<int, Eigen::MatrixXf>>> m_corrParamsDCov;
152  // map<DetType, map<int, map<int, Eigen::VectorXf>>> m_corrParamsD0;
154  map<DetType, map<int, map<int, vector<double>>>> m_simHistsM1;
155  map<DetType, map<int, map<int, vector<double>>>> m_simHistsM2;
157  map<DetType, map<int, map<int, vector<vector<double>>>>> m_allExpEffs;
158  map<DetType, map<int, map<int, vector<PeakInfo>>>> m_allExpPeakFits;
159 
160  std::mt19937_64 m_randGen;
161 };
DetType
Definition: ConfigClasses.h:20
Definition: GammaEfficiencyAnalysis.h:20
vector< PeakInfo > getPeakFits(DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
Definition: GammaEfficiencyAnalysis.cpp:1343
void plotAddBackSimExpComp()
plots simulated and experimental add back factors
Definition: GammaEfficiencyAnalysis.cpp:865
map< int, vector< vector< double > > > calcAddBackFactors()
Definition: GammaEfficiencyAnalysis.cpp:828
map< double, std::pair< double, double > > calcEfficiency(TH1 *a_hist, string a_isotope, double a_activity, double a_countTime, int a_numLines, vector< PeakInfo > &a_peakFits)
Definition: GammaEfficiencyAnalysis.cpp:1190
int addSimulatedData(string a_fileName, int a_numSimEvents)
Definition: GammaEfficiencyAnalysis.cpp:52
void fitTriggerProb(string a_fileNameHighTF, string a_fileNameLowTF, double a_countTimeHigh, double a_countTimeLow, string a_fileNameBKG="", double a_countTimeBKG=-1, string a_isotope="Eu152", int a_numLines=16)
Definition: GammaEfficiencyAnalysis.cpp:416
void loadExperimentalSpectra(string a_fileName, string a_isotope, double a_countTime=1, double a_activity=1, int a_numLines=5)
load calibrated experimental spectra
Definition: GammaEfficiencyAnalysis.cpp:38
map< double, PeakInfo > fitIsotopePeaks(TH1 *a_hist, string a_isotope, double a_activity, double a_countTime, int a_numLines, bool a_draw=false)
Definition: GammaEfficiencyAnalysis.cpp:1130
void setSimHistNumBins(int a_bins)
set number of bins for simulated hists
Definition: GammaEfficiencyAnalysis.cpp:92
void processSimulatedData(string a_configFileName, bool a_resampleEffCorr=false, double a_acceptanceWindow=2.5e-3)
Definition: GammaEfficiencyAnalysis.cpp:100
void plotSimExpComp(DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
plot simulated and experimental efficiency + percent difference
Definition: GammaEfficiencyAnalysis.cpp:1036
TH1F * getSimHist(DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
returns simulated histogram (if it exists)
Definition: GammaEfficiencyAnalysis.cpp:321
void drawEffPeakFits(DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
double getSimEff(double a_gammaEnergy, DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
Definition: GammaEfficiencyAnalysis.cpp:381
map< DetType, map< int, map< int, Eigen::MatrixXf > > > m_corrParamsDCov
Definition: GammaEfficiencyAnalysis.h:135
void setSimHistMaxE(double a_maxEnergy)
Definition: GammaEfficiencyAnalysis.cpp:96
void generateNewTriggerProbs()
Definition: GammaEfficiencyAnalysis.cpp:635
void computeAllExpEffs(bool a_drawFits=false)
Definition: GammaEfficiencyAnalysis.cpp:658
void setNewPeakFit(PeakInfo a_peakInfo, DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
Definition: GammaEfficiencyAnalysis.cpp:1306
void loadBackgroundSpectra(string a_fileName, double a_countTime)
Definition: GammaEfficiencyAnalysis.cpp:46
void useTFIntDiffCorrection(bool a_useCorr)
Definition: GammaEfficiencyAnalysis.cpp:38
map< DetType, map< int, map< int, TH2F * > > > m_angEffHists
Definition: GammaEfficiencyAnalysis.h:137
map< DetType, map< int, map< int, Eigen::VectorXf > > > m_corrParamsD0
Definition: GammaEfficiencyAnalysis.h:134
vector< vector< double > > getExpEffs(DetType a_detType, int a_detID, int a_elemID=0, bool a_addBack=false)
Definition: GammaEfficiencyAnalysis.cpp:780
double calcTriggerProb(DetType a_detType, int a_detID, int a_elemID, double a_energy)
calculate the trigger probability for a given detector at an energy
Definition: GammaEfficiencyAnalysis.cpp:598
void setVerbosity(int a_verbosity)
Definition: GammaEfficiencyAnalysis.cpp:33
GammaEfficiencyAnalysis()
Definition: GammaEfficiencyAnalysis.cpp:22
map< DetType, map< int, map< int, Eigen::VectorXf > > > m_corrParamsD
Definition: GammaEfficiencyAnalysis.h:136
map< double, PeakInfo > m_lastPeakFits
Definition: GammaEfficiencyAnalysis.h:132
vector< double > m_gammaEnergies
Definition: GammaEfficiencyAnalysis.h:133
Definition: SpectrumAnalysis.h:21
Definition: SpectrumAnalysis.h:40