LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
FlashHypothesisComparison.cxx
Go to the documentation of this file.
1 
12 #include "TTree.h"
13 #include "TH1F.h"
14 
16  TH1F* h_h_p, TH1F* h_s_p, TH1F* h_c_p,
17  TH1F* h_h_l, TH1F* h_s_l, TH1F* h_c_l,
18  TH1F* h_h_t, TH1F* h_s_t, TH1F* h_c_t,
19  const unsigned int n_opdet,
20  bool fill)
21 {
22  fTree = tree;
23  fFillTree = fill;
24  fTree->SetName("ctree");
25 
26  fTree->Branch("run",&fRun,"run/i");
27  fTree->Branch("event",&fEvent,"event/i");
28 
29  fHypHist_p = h_h_p;
30  fSimHist_p = h_s_p;
31  fCompareHist_p = h_c_p;
32 
33  fHypHist_p->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
34  fSimHist_p->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
35  fCompareHist_p->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
36 
37  fHypHist_p->SetNameTitle("hHypHist_p","Hypothesis (Prompt);Opdet;PEs");
38  fSimHist_p->SetNameTitle("hSimHist_p","SimPhoton (Prompt);Opdet;PEs");
39  fCompareHist_p->SetNameTitle("hCompareHist_p","Comparison (Hyp - Sim) (Prompt);Opdet;PEs");
40 
41  fTree->Branch("hyp_PEs_p",&fHypPEs_p,"hyp_PEs_p/F");
42  fTree->Branch("hyp_PEsError_p",&fHypPEsError_p,"hyp_PEsError_p/F");
43  fTree->Branch("sim_PEs_p",&fSimPEs_p,"sim_PEs_p/F");
44 
45  fTree->Branch("hyp_Y_p",&fHypY_p,"hyp_Y_p/F");
46  fTree->Branch("sim_Y_p",&fSimY_p,"sim_Y_p/F");
47 
48  fTree->Branch("hyp_RMSY_p",&fHypRMSY_p,"hyp_RMSY_p/F");
49  fTree->Branch("sim_RMSY_p",&fSimRMSY_p,"sim_RMSY_p/F");
50 
51  fTree->Branch("hyp_Z_p",&fHypZ_p,"hyp_Z_p/F");
52  fTree->Branch("sim_Z_p",&fSimZ_p,"sim_Z_p/F");
53 
54  fTree->Branch("hyp_RMSZ_p",&fHypRMSZ_p,"hyp_RMSZ_p/F");
55  fTree->Branch("sim_RMSZ_p",&fSimRMSZ_p,"sim_RMSZ_p/F");
56 
57  fTree->Branch("comp_total_p",&fCompare_p,"comp_total_p/F");
58 
59  fTree->Branch("hHypHist_p",&fHypHist_p);
60  fTree->Branch("hSimHist_p",&fSimHist_p);
61  fTree->Branch("hCompareHist_p",&fCompareHist_p);
62 
63  fHypHist_l = h_h_l;
64  fSimHist_l = h_s_l;
65  fCompareHist_l = h_c_l;
66 
67  fHypHist_l->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
68  fSimHist_l->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
69  fCompareHist_l->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
70 
71  fHypHist_l->SetNameTitle("hHypHist_l","Hypothesis (Late);Opdet;PEs");
72  fSimHist_l->SetNameTitle("hSimHist_l","SimPhoton (Late);Opdet;PEs");
73  fCompareHist_l->SetNameTitle("hCompareHist_l","Comparison (Hyp - Sim) (Late);Opdet;PEs");
74 
75  fTree->Branch("hyp_PEs_l",&fHypPEs_l,"hyp_PEs_l/F");
76  fTree->Branch("hyp_PEsError_l",&fHypPEsError_l,"hyp_PEsError_l/F");
77  fTree->Branch("sim_PEs_l",&fSimPEs_l,"sim_PEs_l/F");
78 
79  fTree->Branch("hyp_Y_l",&fHypY_l,"hyp_Y_l/F");
80  fTree->Branch("sim_Y_l",&fSimY_l,"sim_Y_l/F");
81 
82  fTree->Branch("hyp_RMSY_l",&fHypRMSY_l,"hyp_RMSY_l/F");
83  fTree->Branch("sim_RMSY_l",&fSimRMSY_l,"sim_RMSY_l/F");
84 
85  fTree->Branch("hyp_Z_l",&fHypZ_l,"hyp_Z_l/F");
86  fTree->Branch("sim_Z_l",&fSimZ_l,"sim_Z_l/F");
87 
88  fTree->Branch("hyp_RMSZ_l",&fHypRMSZ_l,"hyp_RMSZ_l/F");
89  fTree->Branch("sim_RMSZ_l",&fSimRMSZ_l,"sim_RMSZ_l/F");
90 
91  fTree->Branch("comp_total_l",&fCompare_l,"comp_total_l/F");
92 
93  fTree->Branch("hHypHist_l",&fHypHist_l);
94  fTree->Branch("hSimHist_l",&fSimHist_l);
95  fTree->Branch("hCompareHist_l",&fCompareHist_l);
96 
97  fHypHist_t = h_h_t;
98  fSimHist_t = h_s_t;
99  fCompareHist_t = h_c_t;
100 
101  fHypHist_t->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
102  fSimHist_t->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
103  fCompareHist_t->SetBins(n_opdet,-0.5,(float)n_opdet - 0.5);
104 
105  fHypHist_t->SetNameTitle("hHypHist_t","Hypothesis (Total);Opdet;PEs");
106  fSimHist_t->SetNameTitle("hSimHist_t","SimPhoton (Total);Opdet;PEs");
107  fCompareHist_t->SetNameTitle("hCompareHist_t","Comparison (Hyp - Sim) (Total);Opdet;PEs");
108 
109  fTree->Branch("hyp_PEs_t",&fHypPEs_t,"hyp_PEs_t/F");
110  fTree->Branch("hyp_PEsError_t",&fHypPEsError_t,"hyp_PEsError_t/F");
111  fTree->Branch("sim_PEs_t",&fSimPEs_t,"sim_PEs_t/F");
112 
113  fTree->Branch("hyp_Y_t",&fHypY_t,"hyp_Y_t/F");
114  fTree->Branch("sim_Y_t",&fSimY_t,"sim_Y_t/F");
115 
116  fTree->Branch("hyp_RMSY_t",&fHypRMSY_t,"hyp_RMSY_t/F");
117  fTree->Branch("sim_RMSY_t",&fSimRMSY_t,"sim_RMSY_t/F");
118 
119  fTree->Branch("hyp_Z_t",&fHypZ_t,"hyp_Z_t/F");
120  fTree->Branch("sim_Z_t",&fSimZ_t,"sim_Z_t/F");
121 
122  fTree->Branch("hyp_RMSZ_t",&fHypRMSZ_t,"hyp_RMSZ_t/F");
123  fTree->Branch("sim_RMSZ_t",&fSimRMSZ_t,"sim_RMSZ_t/F");
124 
125  fTree->Branch("comp_total_t",&fCompare_t,"comp_total_t/F");
126 
127  fTree->Branch("hHypHist_t",&fHypHist_t);
128  fTree->Branch("hSimHist_t",&fSimHist_t);
129  fTree->Branch("hCompareHist_t",&fCompareHist_t);
130 }
131 
133  const unsigned int event,
134  const FlashHypothesisCollection& fhc,
135  const SimPhotonCounter& spc,
136  const std::vector<float>& posY,
137  const std::vector<float>& posZ)
138 {
139  if(fhc.GetVectorSize() != (unsigned int)fHypHist_p->GetNbinsX() ||
140  fhc.GetVectorSize() != spc.PromptPhotonVector().size() ||
141  fhc.GetVectorSize() != posY.size() ||
142  fhc.GetVectorSize() != posZ.size() ){
143  std::cout << (unsigned int)fHypHist_p->GetNbinsX() << " " << spc.PromptPhotonVector().size() << " " << posY.size() << " " << posZ.size() << std::endl;
144  throw std::runtime_error("ERROR in FlashHypothesisComparison: Mismatch in vector sizes.");
145  }
146  fRun = run;
147  fEvent = event;
148 
149  FillFlashHypothesisInfo(fhc,posY,posZ);
150  FillSimPhotonCounterInfo(spc,posY,posZ);
151  FillComparisonInfo(fhc,spc);
152 
153  if(fFillTree) fTree->Fill();
154 }
155 
157  const std::vector<float>& posY,
158  const std::vector<float>& posZ)
159 {
164 
165  for(size_t i=0; i<fhc.GetVectorSize(); i++)
166  fHypHist_p->SetBinContent(i+1,fhc.GetPromptHypothesis().GetHypothesis(i));
167 
172 
173  for(size_t i=0; i<fhc.GetVectorSize(); i++)
174  fHypHist_l->SetBinContent(i+1,fhc.GetLateHypothesis().GetHypothesis(i));
175 
180 
181  for(size_t i=0; i<fhc.GetVectorSize(); i++)
182  fHypHist_t->SetBinContent(i+1,fhc.GetLateHypothesis().GetHypothesis(i));
183 }
184 
186  const std::vector<float>& posY,
187  const std::vector<float>& posZ)
188 {
192 
193  for(size_t i=0; i<spc.PromptPhotonVector().size(); i++)
194  fSimHist_p->SetBinContent(i+1,spc.PromptPhotonVector()[i]);
195 
196  fSimPEs_l = spc.LatePhotonTotal();
199 
200  for(size_t i=0; i<spc.LatePhotonVector().size(); i++)
201  fSimHist_l->SetBinContent(i+1,spc.LatePhotonVector()[i]);
202 
206 
207  for(size_t i=0; i<spc.GetVectorSize(); i++)
208  fSimHist_t->SetBinContent(i+1,spc.TotalPhotonVector(i));
209 
210 }
211 
213  const SimPhotonCounter& spc)
214 {
215  std::vector<float> result_p,result_l,result_t;
219 
220  for(size_t i=0; i<result_p.size(); i++){
221  fCompareHist_p->SetBinContent(i+1,result_p[i]);
222  fCompareHist_l->SetBinContent(i+1,result_l[i]);
223  fCompareHist_t->SetBinContent(i+1,result_t[i]);
224  }
225 
226 }
void FillSimPhotonCounterInfo(const SimPhotonCounter &, const std::vector< float > &, const std::vector< float > &)
float GetTotalPEsError() const
const FlashHypothesis & GetLateHypothesis() const
const std::vector< float > & PromptPhotonVector() const
void FillComparisonInfo(const FlashHypothesisCollection &, const SimPhotonCounter &)
float GetTotalPEs() const
std::vector< float > TotalPhotonVector() const
void GetPosition(const std::vector< float > &, const std::vector< float > &, float &, float &)
void fill(const art::PtrVector< recob::Hit > &hits, int only_plane)
const FlashHypothesis & GetPromptHypothesis() const
void SetOutputObjects(TTree *, TH1F *, TH1F *, TH1F *, TH1F *, TH1F *, TH1F *, TH1F *, TH1F *, TH1F *, const unsigned int, bool fill=true)
float CompareByError(const FlashHypothesis &, const std::vector< float > &, std::vector< float > &)
const std::vector< float > & LatePhotonVector() const
float const & GetHypothesis(size_t i_opdet) const
float PromptPhotonTotal() const
void RunComparison(const unsigned int, const unsigned int, const FlashHypothesisCollection &, const SimPhotonCounter &, const std::vector< float > &, const std::vector< float > &)
const FlashHypothesis & GetTotalHypothesis() const
void FillFlashHypothesisInfo(const FlashHypothesisCollection &, const std::vector< float > &, const std::vector< float > &)
size_t GetVectorSize() const
std::vector< float > const & GetHypothesisVector() const
float LatePhotonTotal() const
Event finding and building.