Berkeley Nuclear Data Software
GenEfficiencyAnalysis.h
Go to the documentation of this file.
1 #include "GENESISConfig.h"
2 
3 #include "TFile.h"
4 #include "TTree.h"
5 #include "TH2.h"
6 #include "TCanvas.h"
7 #include "TChain.h"
8 #include <vector>
9 #include <string>
10 #include <map>
11 #include <random>
12 using std::vector;
13 using std::string;
14 using std::map;
15 #include <Eigen/Dense>
21 {
22  public:
24 
25  GenEfficiencyAnalysis(string a_configFileName);
28  void setRNGSeed(int a_seed);
30  void addConfigFile(string a_configFileName);
32  void loadSimulationEvents(string a_simFileName, int a_numEvents);
36  void addLYCovarianceData(string a_fileName);
46  void processSimData(bool a_useExpTimeRes=false,
47  double a_timeResolution=1.,
48  bool a_useLYResolution=true,
49  bool a_resampleLYCalib=false,
50  double a_minLight=-1);
52  void loadSinglesData(string a_singlesFileName,
53  double a_countingTime,
54  double a_minTime=0);
56  void loadCoincidenceData(string a_coinFileName,
57  double a_countingTime,
58  double a_minTime=0);
60  void loadExpGammaGammaData(string a_ggFileName);
62  int loadCoincidenceTreeData(string a_coinFileName,
63  double a_countingTime);
65  void prepareCoinTreeData(bool a_useADC=false);
67  TH2F* getExpTOFLightHist(int a_detID);
69  TH2F* getTotalCoinHist(int a_detID, string a_histNameRoot);
71  double calcAverageGammaEff(int a_detID,
72  double a_minEnergyBin=1,
73  double a_maxEnergyBin=-1);
75  vector<double> plotTOFCf252Comparison(int a_detNum,
76  bool a_normalize=false);
77 
80  vector<TH2F*> plotEDepCf252Comparison(int a_detNum);
84  void calcMeanTOF(int a_meanSigma=0);
87  map<double, double> respKDE1DInt(double a_mean,
88  double a_sig,
89  double a_nbins,
90  double a_max);
93  void calcCompEff();
94 
96  vector<double> plotSliceEDep(vector<TH2F*> a_hists, int a_bin, bool a_normalize=false);
99  void plotSinglesComparison();
102  TCanvas* plotCumulativeEDepSlice(vector<TH2F*> a_hists, int a_detNum);
104  TCanvas* plotCumulativeEDepSlice(int a_detNum);
107  void plotCumulativeEDepSliceAllDet(int a_simExp);
110  vector<double> plotPISingles(int a_detNum, bool a_draw=false);
116  map<int, std::pair<double,double>> plotSinglesComparison(int a_xAxisType);
119  void plotChi2SliceEDep(vector<TH2F*> a_hists, int a_detNum);
120 
121  void plotExpSinglesEff(int a_xAxisType);
122 
125  map<double,std::pair<double,double>> plotRelativeEnergyEff(int a_detNum,
126  bool a_drawSim=true,
127  std::pair<double,double> a_singlesEffs={1.,5.});
129  TH1F* interpolateCf252Spectrum(vector<double> a_binEdges,
130  vector<double> a_flux,
131  double a_numNeutrons,
132  bool a_draw=false);
137  double calcWattIntegral(double a_lowerE, double a_upperE);
139  TH2F* getSimHist(int a_detNum);
143  map<double,std::pair<double,double>> calcUniformEfficiency(int a_detNum,
144  double a_lowEnergy,
145  double a_upperEnergy);
147  TH2F* buildResponseFunction(int a_detNum);
148  map<int, Eigen::MatrixXf> m_respFuncs;
149 private:
150  GenesisConfig m_conf;
151  // TTree* m_simTree;
152  TChain* m_simTree;
153  double m_simTime;
154  double m_simLY;
155  int m_simChNum;
156  map<int, TH2F*> m_simHists;
157  map<int, TH2F*> m_simHistsEnerLight;
158  map<int, map<int, TH1*>> m_expTimeResHists;
159  vector<string> m_singlesFilesNames;
160  vector<string> m_coinFilesNames;
161  vector<string> m_coinTreeFilesNames;
162  TChain* m_expCoinTree;
163  double m_expCADC;
164  double m_expCEgam;
165  double m_expCTOF;
166  int m_expCDet;
167  bool m_hasCoinTreeFile=false;
168  map<int, TH2F*> m_expCoinTreeHists_tofLight;
169  map<int, map<int,TH2F*>> m_expCoinTreeHists_egamLight;
170  map<int, TH2F*> m_expCoinTreeBKGHists;
171 
172  double m_sourceStrengthFrac = .8303;//.8133;//.8551;
173  double m_PSDThreshold=.135;
174  long int m_numSimEvents;
175 
176  vector<double> m_singlesTimes;
177  vector<double> m_coinsTimes;
178  vector<double> m_minimumTimesCoin;
179  vector<double> m_minimumTimesSingles;
180 
181  double m_totalSinglesTime=0.;
182  double m_totalCoinsTime=0.;
183  double m_totalCoinsTreeTime=0.;
184 
185  int m_nBins_tof = 5000;
186  // int m_nBins_tof = 5;
187  int m_tofMinMax = 500;
188  int m_nBins_egam = 100;
189  int m_egamMax = 65536;//10;
190  int m_nBins_light = 4096;
191  // int m_nBins_light = 8;
192  int m_lightMax = 10;
193 
194  map<int, Eigen::MatrixXf> m_lyCovariances;
195  map<int, Eigen::VectorXf> m_lyParams0;
196 
197  std::mt19937_64 m_randGen;
198 // map<int,int> m_detConv = {{0,8},{1,7},{2,6},{3,5},{4,4},{5,3},{6,2},{7,1},{8,0},{9,9},{10,10},
199 // {11,18},{13,16},{14,15},{15,14},{16,13},{17,17},{18,11},
200 // // {19,19},{20,25},{21,24},{22,23},{23,22},{24,21},{25,20},{26,26}};
201 // {19,19},{20,20},{21,21},{22,22},{23,23},{24,24},{25,25}};//,{26,26}};
202 
203 };
Definition: GenEfficiencyAnalysis.h:21
void plotChi2SliceEDep(vector< TH2F * > a_hists, int a_detNum)
Definition: GenEfficiencyAnalysis.cpp:1951
void calcCompEff()
Definition: GenEfficiencyAnalysis.cpp:1017
map< double, double > respKDE1DInt(double a_mean, double a_sig, double a_nbins, double a_max)
Definition: GenEfficiencyAnalysis.cpp:971
TH1F * interpolateCf252Spectrum(vector< double > a_binEdges, vector< double > a_flux, double a_numNeutrons, bool a_draw=false)
get Cf252 spectrum interpolated for correct TOF bins
Definition: GenEfficiencyAnalysis.cpp:2076
double calcAverageGammaEff(int a_detID, double a_minEnergyBin=1, double a_maxEnergyBin=-1)
correction for coin data by spectrum-average LaBr efficiency
Definition: GenEfficiencyAnalysis.cpp:615
TCanvas * plotCumulativeEDepSlice(vector< TH2F * > a_hists, int a_detNum)
Definition: GenEfficiencyAnalysis.cpp:1141
vector< double > plotSliceEDep(vector< TH2F * > a_hists, int a_bin, bool a_normalize=false)
takes output of plotEDepCf252Comparison {sim, exp} and plots time-slices
Definition: GenEfficiencyAnalysis.cpp:1072
void generateResampledConfig()
modifies the stored GenesisConfig using the LY covariances
Definition: GenEfficiencyAnalysis.cpp:129
void addLYCovarianceData(string a_fileName)
Definition: GenEfficiencyAnalysis.cpp:89
void prepareCoinTreeData(bool a_useADC=false)
make the experimental tof vs light histograms
Definition: GenEfficiencyAnalysis.cpp:509
void loadCoincidenceData(string a_coinFileName, double a_countingTime, double a_minTime=0)
loads hists from GenCoinAnalysis::buildOutgoingTOFHist()
Definition: GenEfficiencyAnalysis.cpp:332
vector< double > plotTOFCf252Comparison(int a_detNum, bool a_normalize=false)
plot sim and exp TOF plots
Definition: GenEfficiencyAnalysis.cpp:670
vector< TH2F * > plotEDepCf252Comparison(int a_detNum)
Definition: GenEfficiencyAnalysis.cpp:810
map< double, std::pair< double, double > > calcUniformEfficiency(int a_detNum, double a_lowEnergy, double a_upperEnergy)
Definition: GenEfficiencyAnalysis.cpp:2710
void addConfigFile(string a_configFileName)
loads in a config file
Definition: GenEfficiencyAnalysis.cpp:54
void plotSinglesComparison()
Definition: GenEfficiencyAnalysis.cpp:1234
double calcWattIntegral(double a_lowerE, double a_upperE)
Definition: GenEfficiencyAnalysis.cpp:2102
vector< double > plotPISingles(int a_detNum, bool a_draw=false)
Definition: GenEfficiencyAnalysis.cpp:1433
void loadSinglesData(string a_singlesFileName, double a_countingTime, double a_minTime=0)
loads hists from GENESISPostProcessing::accumulateSinglesQuantities()
Definition: GenEfficiencyAnalysis.cpp:322
TH2F * getTotalCoinHist(int a_detID, string a_histNameRoot)
sums up histograms from each loaded coincidence data file
Definition: GenEfficiencyAnalysis.cpp:583
void setRNGSeed(int a_seed)
Definition: GenEfficiencyAnalysis.cpp:49
void loadSimulationEvents(string a_simFileName, int a_numEvents)
loads simulation, needs total number of events simulated
Definition: GenEfficiencyAnalysis.cpp:68
void calcMeanTOF(int a_meanSigma=0)
Definition: GenEfficiencyAnalysis.cpp:898
int loadCoincidenceTreeData(string a_coinFileName, double a_countingTime)
load experimental n/g scint-scint coincidence tree
Definition: GenEfficiencyAnalysis.cpp:442
map< int, Eigen::MatrixXf > m_respFuncs
Definition: GenEfficiencyAnalysis.h:148
void processSimData(bool a_useExpTimeRes=false, double a_timeResolution=1., bool a_useLYResolution=true, bool a_resampleLYCalib=false, double a_minLight=-1)
Definition: GenEfficiencyAnalysis.cpp:176
void loadExpGammaGammaData(string a_ggFileName)
loads gamma-gamma experimental data and prepares slices of TOF in light
Definition: GenEfficiencyAnalysis.cpp:341
map< double, std::pair< double, double > > plotRelativeEnergyEff(int a_detNum, bool a_drawSim=true, std::pair< double, double > a_singlesEffs={1., 5.})
Definition: GenEfficiencyAnalysis.cpp:2133
TH2F * getSimHist(int a_detNum)
returns simulation histogram: tof vs light
Definition: GenEfficiencyAnalysis.cpp:2637
TH2F * buildResponseFunction(int a_detNum)
return neutron energy vs light 2d histogram – built in processSimData
Definition: GenEfficiencyAnalysis.cpp:2759
GenEfficiencyAnalysis()
Definition: GenEfficiencyAnalysis.cpp:21
void plotExpSinglesEff(int a_xAxisType)
Definition: GenEfficiencyAnalysis.cpp:1991
void plotCumulativeEDepSliceAllDet(int a_simExp)
Definition: GenEfficiencyAnalysis.cpp:1838
TH2F * getExpTOFLightHist(int a_detID)
get tof vs light hist for given detector
Definition: GenEfficiencyAnalysis.cpp:573
Definition: GENESISConfig.h:15