FairRoot/PandaRoot
tracking/trackingQA/QA_histos.C
Go to the documentation of this file.
1 
11  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(1), "Spurious found");
12  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(2), "Partially found");
13  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(3), "Fully found");
14  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(-1), "Not found secondary");
15  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(-2), "Not found primary");
16  hist->GetXaxis()->SetBinLabel(hist->FindFixBin(-3), "Found");
17 }
18 
19 int QA_histos(TString prefix, TString trackBranch, Bool_t forward = kFALSE) {
20 
21  PndFileNameCreator creator(prefix.Data());
22  TString extension;
23  if (forward == kFALSE)
24  extension = "trackingQA_";
25  else extension = "trackingQAfwd_";
26  extension.Append(trackBranch);
27  TString simFile = creator.GetSimFileName();
28  TString trackQAFile = creator.GetCustomFileName(extension.Data());
29  extension.Append("_histos");
30  TString outputFile = creator.GetCustomFileName(extension.Data());
31 
32  std::cout << "InFileName: " << simFile.Data() << std::endl;
33  std::cout << "TrackQAFile: " << trackQAFile.Data() << std::endl;
34  std::cout << "OutputFileName: " << outputFile.Data() << std::endl;
35 
36  TFile fileqa(trackQAFile.Data());
37  TTree *simtree = (TTree*) fileqa.Get("pndsim");
38  std::cout << "SimTree: " << simtree << std::endl;
39  std::cout << "SimTree.GetEntries() " << simtree->GetEntries() << std::endl;
40  simtree->AddFriend("pndsim", simFile.Data());
41 
42  TCut cut = "RecoTrackInfo.GetQuality() > 0";
43 
54  double effmin = -0.1;
55  double effmax = 1.1;
56 
57  TH1F *eff_mvdpix = new TH1F("eff_mvdpix", "MVD pixel efficiency", 100, effmin, effmax);
58  TH1F *eff_mvdstr = new TH1F("eff_mvdstr", "MVD strip efficiency", 100, effmin, effmax);
59  TH1F *eff_stt = new TH1F("eff_stt", "STT efficiency", 100, effmin, effmax);
60  TH1F *eff_gem = new TH1F("eff_gem", "GEM efficiency", 100, effmin, effmax);
61  TH1F *eff_fts = new TH1F("eff_fts", "FTS efficiency", 100, effmin, effmax);
62 
63 
64  std::cout << "Before first Draw!" << std::endl;
65 
66  simtree->Draw("RecoTrackInfo.GetMvdPixelEfficiency() >> eff_mvdpix", cut, "goff");
67  simtree->Draw("RecoTrackInfo.GetMvdStripEfficiency() >> eff_mvdstr", cut, "goff");
68  simtree->Draw("RecoTrackInfo.GetSttEfficiency() >> eff_stt", cut, "goff");
69  simtree->Draw("RecoTrackInfo.GetGemEfficiency() >> eff_gem", cut, "goff");
70  simtree->Draw("RecoTrackInfo.GetFtsEfficiency() >> eff_fts", cut, "goff");
71 //
72  TH1F *eff_glo = new TH1F("eff_glo", "global efficiency", 100, effmin, effmax);
73  simtree->Draw("RecoTrackInfo.GetEfficiency() >> eff_glo", cut, "goff");
74 
75  std::cout << "Efficiencies filled!" << std::endl;
85  double purmin = -0.1;
86  double purmax = 1.1;
87 
88  TH1F *pur_mvdpix = new TH1F("pur_mvdpix", "MVD pixel purity", 100, purmin, purmax);
89  TH1F *pur_mvdstr = new TH1F("pur_mvdstr", "MVD strip purity", 100, purmin, purmax);
90  TH1F *pur_stt = new TH1F("pur_stt", "STT purity", 100, purmin, purmax);
91  TH1F *pur_gem = new TH1F("pur_gem", "GEM purity", 100, purmin, purmax);
92  TH1F *pur_fts = new TH1F("pur_fts", "FTS purity", 100, purmin, purmax);
93 
94  simtree->Draw("RecoTrackInfo.GetMvdPixelPurity() >> pur_mvdpix", cut, "goff");
95  simtree->Draw("RecoTrackInfo.GetMvdStripPurity() >> pur_mvdstr", cut, "goff");
96  simtree->Draw("RecoTrackInfo.GetSttPurity() >> pur_stt", cut, "goff");
97  simtree->Draw("RecoTrackInfo.GetGemPurity() >> pur_gem", cut, "goff");
98  simtree->Draw("RecoTrackInfo.GetFtsPurity() >> pur_fts", cut, "goff");
99 
100 
101  TH1F *pur_glo = new TH1F("pur_glo", "global purity", 100, purmin, purmax);
102  simtree->Draw("RecoTrackInfo.GetPurity() >> pur_glo", cut, "goff");
103 
104  std::cout << "Purities filled!" << std::endl;
110  TH1F *hnofrecotracks = new TH1F("hnofrecotracks", "# of reco tracks associated to the same MC track", 20, 0, 20);
111  simtree->Draw("MCTrackInfo.GetNofRecoTracks() >> hnofrecotracks", cut, "goff");
112 
113  TH1F *hnofMCtracks = new TH1F("hnofMCtracks", "# of MC tracks associated to the same reco track", 20, 0, 20);
114  simtree->Draw("RecoTrackInfo.GetNofMCTracks() >> hnofMCtracks", cut, "goff");
115 
116  std::cout << "RecoTracks filled!" << std::endl;
117 
127  double pmin = -1;
128  double pmax = 1;
129  TH1F *hdelta_p_first = new TH1F("hdelta_p_first", "#Delta p @ first hit", 100, pmin, pmax);
130  TH1F *hdelta_pz_first = new TH1F("hdelta_pz_first", "#Delta pz @ first hit", 100, pmin, pmax);
131  TH1F *hdelta_pt_first = new TH1F("hdelta_pt_first", "#Delta pt @ first hit", 100, pmin, pmax);
132 
133  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Mag() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumFirst().Mag() >> hdelta_p_first", cut, "goff");
134  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Perp() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumFirst().Perp() >> hdelta_pt_first", cut, "goff");
135  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Z() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumFirst().Z() >> hdelta_pz_first", cut, "goff");
136 
137  double thetamin = -5;
138  double thetamax = 5;
139  int ntheta = (thetamax - thetamin)*10;
140  TH1F *hdelta_theta_first = new TH1F("hdelta_theta_first", "#Delta #theta @ first hit", ntheta, thetamin, thetamax);
141  simtree->Draw("(RecoTrackInfo.GetMomentumFirst().Theta() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumFirst().Theta()) * TMath::RadToDeg() >> hdelta_theta_first", cut, "goff");
142 
143  double phimin = -5;
144  double phimax = 5;
145  int nphi = (phimax - phimin)*10;
146  TH1F *hdelta_phi_first = new TH1F("hdelta_phi_first", "#Delta #phi @ first hit", nphi, phimin, phimax);
147  simtree->Draw("(RecoTrackInfo.GetMomentumFirst().Phi() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumFirst().Phi()) * TMath::RadToDeg() >> hdelta_phi_first", cut, "goff");
148 
149 
150  double xmin = -0.2;
151  double xmax = 0.2;
152  TH1F *hdelta_x_first = new TH1F("hdelta_x_first", "#Delta x @ first hit", 100, xmin, xmax);
153  TH1F *hdelta_y_first = new TH1F("hdelta_y_first", "#Delta y @ first hit", 100, xmin, xmax);
154  simtree->Draw("RecoTrackInfo.GetPositionFirst().X() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionFirst().X() >> hdelta_x_first", cut, "goff");
155  simtree->Draw("RecoTrackInfo.GetPositionFirst().Y() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionFirst().Y() >> hdelta_y_first", cut, "goff");
156 
157  double zmin = -0.2;
158  double zmax = 0.2;
159  TH1F *hdelta_z_first = new TH1F("hdelta_z_first", "#Delta z @ first hit", 100, zmin, zmax);
160  simtree->Draw("RecoTrackInfo.GetPositionFirst().Z() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionFirst().Z() >> hdelta_z_first", cut, "goff");
161 
162  double rmax = xmax;
163  TH1F *hdelta_r_first = new TH1F("hdelta_r_first", "#Delta r @ first hit", 100, -rmax, rmax);
164  simtree->Draw("RecoTrackInfo.GetPositionFirst().Perp() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionFirst().Perp() >> hdelta_r_first", cut, "goff");
165 
166  std::cout << "DeltaPFirst filled!" << std::endl;
167 
177  TH1F *hdelta_p_last = new TH1F("hdelta_p_last", "#Delta p @ last hit", 100, pmin, pmax);
178  TH1F *hdelta_pz_last = new TH1F("hdelta_pz_last", "#Delta pz @ last hit", 100, pmin, pmax);
179  TH1F *hdelta_pt_last = new TH1F("hdelta_pt_last", "#Delta pt @ last hit", 100, pmin, pmax);
180 
181  simtree->Draw("RecoTrackInfo.GetMomentumLast().Mag() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumLast().Mag() >> hdelta_p_last", cut, "goff");
182  simtree->Draw("RecoTrackInfo.GetMomentumLast().Perp() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumLast().Perp() >> hdelta_pt_last", cut, "goff");
183  simtree->Draw("RecoTrackInfo.GetMomentumLast().Z() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumLast().Z() >> hdelta_pz_last", cut, "goff");
184 
185  TH1F *hdelta_theta_last = new TH1F("hdelta_theta_last", "#Delta #theta @ last hit", ntheta, thetamin, thetamax);
186  simtree->Draw("(RecoTrackInfo.GetMomentumLast().Theta() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumLast().Theta()) * TMath::RadToDeg() >> hdelta_theta_last", cut, "goff");
187 
188  TH1F *hdelta_phi_last = new TH1F("hdelta_phi_last", "#Delta #phi @ last hit", nphi, phimin, phimax);
189  simtree->Draw("(RecoTrackInfo.GetMomentumLast().Phi() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetMomentumLast().Phi()) * TMath::RadToDeg() >> hdelta_phi_last", cut, "goff");
190 
191  TH1F *hdelta_x_last = new TH1F("hdelta_x_last", "#Delta x @ last hit", 100, xmin, xmax);
192  TH1F *hdelta_y_last = new TH1F("hdelta_y_last", "#Delta y @ last hit", 100, xmin, xmax);
193  simtree->Draw("RecoTrackInfo.GetPositionLast().X() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionLast().X() >> hdelta_x_last", cut, "goff");
194  simtree->Draw("RecoTrackInfo.GetPositionLast().Y() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionLast().Y() >> hdelta_y_last", cut, "goff");
195 
196  double zminlast = -0.6;
197  double zmaxlast = 0.6;
198  TH1F *hdelta_z_last = new TH1F("hdelta_z_last", "#Delta z @ last hit", 100, zminlast, zmaxlast);
199  simtree->Draw("RecoTrackInfo.GetPositionLast().Z() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionLast().Z() >> hdelta_z_last", cut, "goff");
200 
201  TH1F *hdelta_r_last = new TH1F("hdelta_r_last", "#Delta r @ last hit", 100, -rmax, rmax);
202  simtree->Draw("RecoTrackInfo.GetPositionLast().Perp() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetPositionLast().Perp() >> hdelta_r_last", cut, "goff");
203 
204  std::cout << "DeltaPLast filled!" << std::endl;
208  double chargemin = -3;
209  double chargemax = 3;
210  int ncharge = chargemax - chargemin;
211  TH1F *hdelta_charge = new TH1F("hdelta_charge", "#Delta #charge", ncharge, chargemin, chargemax);
212  simtree->Draw("RecoTrackInfo.GetCharge() - MCTrackInfo[RecoTrackInfo.GetIdealTrackId()].GetCharge() >> hdelta_charge", cut, "goff");
213 
217  TH1F *hthetagen = new TH1F("hthetagen", "mc theta dist", 180, 0, 180);
218  TH1F *heffintheta = new TH1F("heffintheta", "hit efficiency vs #theta", 180, 0, 180);
219 
220  TCut cut_mc = "";
221  TCut cut_rec = cut_mc && "MCTrackInfo.GetRecoTrackID() != -1";
222  cut_rec = cut_rec && "RecoTrackInfo[MCTrackInfo.GetRecoTrackID()].GetEfficiency() > 0.8";
223  cut_rec = cut_rec && "RecoTrackInfo[MCTrackInfo.GetRecoTrackID()].IsClone() == 0";
224 
225  TCut cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
226  TCut cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
227  TCut cut_rec_prim = cut_rec && cut_mc_prim;
228  TCut cut_rec_sec = cut_rec && cut_mc_sec;
229 
230  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> hthetagen", cut_mc, "goff");
231  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> heffintheta", cut_rec, "goff");
232 
233  hthetagen->Sumw2();
234  heffintheta->Sumw2();
235  heffintheta->Divide(hthetagen);
236 
237 
241  TH1F *hphigen = new TH1F("hphigen", "mc phi dist", 180, -180, 180);
242  TH1F *heffinphi = new TH1F("heffinphi", "hit efficiency vs #phi", 180, -180, 180);
243  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> hphigen", cut_mc, "goff");
244  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> heffinphi", cut_rec, "goff");
245 
246  hphigen->Sumw2();
247  heffinphi->Sumw2();
248  heffinphi->Divide(hphigen);
249 
254  TH1F *hMcPDist = new TH1F("hMcPDist","Momentum Distribution", 30,0,15);
255  TH1F *hMcPtDist = new TH1F("hMcPtDist","Transv. Momentum Distribution", 30,0,15);
256  TH1F *hMcPlDist = new TH1F("hMcPlDist","Long. Momentum Distribution", 30,0,15);
257 
258  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Mag()>>hMcPDist","","goff");
259  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp()>>hMcPtDist","","goff");
260  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z()>>hMcPlDist","","goff");
261 
262  Float_t maxP = hMcPDist->GetBinLowEdge(hMcPDist->FindLastBinAbove(0)) + hMcPDist->GetBinWidth(hMcPDist->FindLastBinAbove(0));
263  Float_t maxPt = hMcPtDist->GetBinLowEdge(hMcPtDist->FindLastBinAbove(0)) + hMcPtDist->GetBinWidth(hMcPtDist->FindLastBinAbove(0));
264  Float_t maxPl = hMcPlDist->GetBinLowEdge(hMcPlDist->FindLastBinAbove(0)) + hMcPlDist->GetBinWidth(hMcPlDist->FindLastBinAbove(0));
265  Float_t minPl = hMcPlDist->GetBinLowEdge(hMcPlDist->FindFirstBinAbove(0)) + hMcPlDist->GetBinWidth(hMcPlDist->FindFirstBinAbove(0));
266 
267  if (minPl > 0) minPl = 0;
268  Int_t Pres = 20;
269 
270  std::cout << "MaxP = " << maxP << " : MaxPt = " << maxPt << " : minPl = " << minPl << " : maxPl = " << maxPl << " : P res = " << Pres << std::endl;
271 
272 
276  TH1F *hmomgen = new TH1F("hmomgen", "mc mom dist", maxP*Pres, 0, maxP);
277  TH1F *heffinmom = new TH1F("heffinmom", "hit efficiency vs mom", maxP*Pres, 0, maxP);
278  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Mag() >> hmomgen", cut_mc, "goff");
279  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Mag() >> heffinmom", cut_rec, "goff");
280 
281  hmomgen->Sumw2();
282  heffinmom->Sumw2();
283  heffinmom->Divide(hmomgen);
284 
288  TH1F *hptgen = new TH1F("hptgen", "mc pt dist", maxPt * Pres, 0, maxPt);
289  TH1F *heffinpt = new TH1F("heffinpt", "hit efficiency vs pt", maxPt * Pres, 0, maxPt);
290 
291  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hptgen", cut_mc, "goff");
292  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> heffinpt", cut_rec, "goff");
293 
294  hptgen->Sumw2();
295  heffinpt->Sumw2();
296  heffinpt->Divide(hptgen);
297 
298  TH1F *hptgenprim = new TH1F("hptgenprim", "mc pt dist prim", maxPt * Pres, 0, maxPt);
299  TH1F *heffinptprim = new TH1F("heffinptprim", "hit efficiency vs pt prim", maxPt * Pres, 0, maxPt);
300 
301  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hptgenprim", cut_mc_prim, "goff");
302  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> heffinptprim", cut_rec_prim, "goff");
303 
304  hptgenprim->Sumw2();
305  heffinptprim->Sumw2();
306  heffinptprim->Divide(hptgenprim);
307 
308  TH1F *hptgensec = new TH1F("hptgensec", "mc pt dist sec", maxPt * Pres, 0, maxPt);
309  TH1F *heffinptsec = new TH1F("heffinptsec", "hit efficiency vs pt sec", maxPt * Pres, 0, maxPt);
310 
311  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hptgensec", cut_mc_sec, "goff");
312  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> heffinptsec", cut_rec_sec, "goff");
313 
314  hptgensec->Sumw2();
315  heffinptsec->Sumw2();
316  heffinptsec->Divide(hptgensec);
317 
321  TH1F *hplgen = new TH1F("hplgen", "mc pl dist", (maxPl - minPl)*Pres, minPl, maxPl);
322  TH1F *heffinpl = new TH1F("heffinpl", "hit efficiency vs pl", (maxPl - minPl)*Pres, minPl, maxPl);
323  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> hplgen", cut_mc, "goff");
324  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> heffinpl", cut_rec, "goff");
325 
326  hplgen->Sumw2();
327  heffinpl->Sumw2();
328  heffinpl->Divide(hplgen);
329 
334  TH1F* hQualityPrim = new TH1F("hQualityPrim", "Quality for possible primary tracks", 8, -3.5, 4.5);
335  TH1F* hQualitySec = new TH1F("hQualitySec", "Quality for possible secondary tracks", 8, -3.5, 4.5);
336  hQualityPrim->SetStats(kFALSE);
337  hQualitySec->SetStats(kFALSE);
338 
339 
340  LabelQualyHistogram(hQualityPrim);
341  LabelQualyHistogram(hQualitySec);
342 
343  simtree->Draw("MCTrackInfo.GetQuality()>>hQualityPrim", cut_mc_prim,"goff");
344  simtree->Draw("MCTrackInfo.GetQuality()>>hQualitySec", cut_mc_sec,"goff");
345 
346  // std::cout << "Primary Tracks: " << nPrim << " Secondary Tracks: " << nSec << std::endl;
347  hQualityPrim->Scale(1.0/hQualityPrim->Integral());
348  hQualitySec->Scale(1.0/hQualitySec->Integral());
349 // NormalizeHistogram(hQualityPrim);
350 // NormalizeHistogram(hQualitySec);
351  hQualityPrim->SetBinContent(1, 1.0 - hQualityPrim->GetBinContent(2));
352  hQualitySec->SetBinContent(1, 1.0 - hQualitySec->GetBinContent(3));
353 
354 
355 
356 
357  // constant->SetParameter(0, nPrim);
358  // hQualityPrim->Divide(constant, nPrim);
359  // constant->SetParameter(0, nSec);
360  // hQualitySec->Divide(constant, nSec);
361 
365  cut_rec = "MCTrackInfo.GetQuality() > 0";
366 
367  cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
368  cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
369  cut_rec_prim = cut_rec && cut_mc_prim;
370  cut_rec_sec = cut_rec && cut_mc_sec;
371 
372  TH1F* htfeffintheta = new TH1F("htfeffintheta", "Track finding efficiency vs theta", 180, 0, 180);
373  TH1F* hmctfeffintheta = new TH1F("hmctfeffintheta", "Track finding efficiency vs theta", 180, 0, 180);
374  htfeffintheta->Sumw2();
375  hmctfeffintheta->Sumw2();
376  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> htfeffintheta", cut_rec, "goff");
377  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> hmctfeffintheta", cut_mc, "goff");
378  htfeffintheta->Divide(hmctfeffintheta);
379 
380  TH1F *htfeffinthetaprim = new TH1F("htfeffinthetaprim", "Track finding efficiency vs theta prim", 180, 0, 180);
381  TH1F *hmctfeffinthetaprim = new TH1F("hmctfeffinthetaprim", "Track finding efficiency vs theta prim", 180, 0, 180);
382  htfeffinthetaprim->Sumw2();
383  hmctfeffinthetaprim->Sumw2();
384  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> htfeffinthetaprim", cut_rec_prim, "goff");
385  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> hmctfeffinthetaprim", cut_mc_prim, "goff");
386  htfeffinthetaprim->Divide(hmctfeffinthetaprim);
387 
388  TH1F *htfeffinthetasec = new TH1F("htfeffinthetasec", "Track finding efficiency vs theta sec", 180, 0, 180);
389  TH1F *hmctfeffinthetasec = new TH1F("hmctfeffinthetasec", "Track finding efficiency vs theta sec", 180, 0, 180);
390  htfeffinthetasec->Sumw2();
391  hmctfeffinthetasec->Sumw2();
392  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> htfeffinthetasec", cut_rec_sec, "goff");
393  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> hmctfeffinthetasec", cut_mc_sec, "goff");
394  htfeffinthetasec->Divide(hmctfeffinthetasec);
395 
399  cut_rec = "MCTrackInfo.GetQuality() > 0";
400 
401  cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
402  cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
403  cut_rec_prim = cut_rec && cut_mc_prim;
404  cut_rec_sec = cut_rec && cut_mc_sec;
405 
406  TH1F* htfeffinphi = new TH1F("htfeffinphi", "Track finding efficiency vs phi", 180, -180, 180);
407  TH1F* hmctfeffinphi = new TH1F("hmctfeffinphi", "Track finding efficiency vs phi", 180, -180, 180);
408  htfeffinphi->Sumw2();
409  hmctfeffinphi->Sumw2();
410  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> htfeffinphi", cut_rec, "goff");
411  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> hmctfeffinphi", cut_mc, "goff");
412  htfeffinphi->Divide(hmctfeffinphi);
413 
414  TH1F *htfeffinphiprim = new TH1F("htfeffinphiprim", "Track finding efficiency vs phi prim", 180, -180, 180);
415  TH1F *hmctfeffinphiprim = new TH1F("hmctfeffinphiprim", "Track finding efficiency vs phi prim", 180, -180, 180);
416  htfeffinphiprim->Sumw2();
417  hmctfeffinphiprim->Sumw2();
418  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> htfeffinphiprim", cut_rec_prim, "goff");
419  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> hmctfeffinphiprim", cut_mc_prim, "goff");
420  htfeffinphiprim->Divide(hmctfeffinphiprim);
421 
422  TH1F *htfeffinphisec = new TH1F("htfeffinphisec", "Track finding efficiency vs phi sec", 180, -180, 180);
423  TH1F *hmctfeffinphisec = new TH1F("hmctfeffinphisec", "Track finding efficiency vs phi sec", 180, -180, 180);
424  htfeffinphisec->Sumw2();
425  hmctfeffinphisec->Sumw2();
426  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> htfeffinphisec", cut_rec_sec, "goff");
427  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> hmctfeffinphisec", cut_mc_sec, "goff");
428  htfeffinphisec->Divide(hmctfeffinphisec);
429 
433  cut_rec = "MCTrackInfo.GetQuality() > 0";
434 
435  cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
436  cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
437  cut_rec_prim = cut_rec && cut_mc_prim;
438  cut_rec_sec = cut_rec && cut_mc_sec;
439 
440  TH1F* htfeffinpt = new TH1F("htfeffinpt", "Track finding efficiency vs pt", maxPt * Pres, 0, maxPt);
441  TH1F* hmctfeffinpt = new TH1F("hmctfeffinpt", "Track finding efficiency vs pt", maxPt * Pres, 0, maxPt);
442  htfeffinpt->Sumw2();
443  hmctfeffinpt->Sumw2();
444  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> htfeffinpt", cut_rec, "goff");
445  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hmctfeffinpt", cut_mc, "goff");
446  htfeffinpt->Divide(hmctfeffinpt);
447 
448  TH1F *htfeffinptprim = new TH1F("htfeffinptprim", "Track finding efficiency vs pt prim", maxPt * Pres, 0, maxPt);
449  TH1F *hmctfeffinptprim = new TH1F("hmctfeffinptprim", "Track finding efficiency vs pt prim", maxPt * Pres, 0, maxPt);
450  htfeffinptprim->Sumw2();
451  hmctfeffinptprim->Sumw2();
452  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> htfeffinptprim", cut_rec_prim, "goff");
453  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hmctfeffinptprim", cut_mc_prim, "goff");
454  htfeffinptprim->Divide(hmctfeffinptprim);
455 
456  TH1F *htfeffinptsec = new TH1F("htfeffinptsec", "Track finding efficiency vs pt sec", maxPt * Pres, 0, maxPt);
457  TH1F *hmctfeffinptsec = new TH1F("hmctfeffinptsec", "Track finding efficiency vs pt sec", maxPt * Pres, 0, maxPt);
458  htfeffinptsec->Sumw2();
459  hmctfeffinptsec->Sumw2();
460  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> htfeffinptsec", cut_rec_sec, "goff");
461  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hmctfeffinptsec", cut_mc_sec, "goff");
462  htfeffinptsec->Divide(hmctfeffinptsec);
463 
467  cut_rec = "MCTrackInfo.GetQuality() > 0";
468 
469  cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
470  cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
471  cut_rec_prim = cut_rec && cut_mc_prim;
472  cut_rec_sec = cut_rec && cut_mc_sec;
473 
474  TH1F* htfeffinpl = new TH1F("htfeffinpl", "Track finding efficiency vs pl", (maxPl - minPl)*Pres, minPl, maxPl);
475  TH1F* hmctfeffinpl = new TH1F("hmctfeffinpl", "Track finding efficiency vs pl", (maxPl - minPl)*Pres, minPl, maxPl);
476  htfeffinpl->Sumw2();
477  hmctfeffinpl->Sumw2();
478  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> htfeffinpl", cut_rec, "goff");
479  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> hmctfeffinpl", cut_mc, "goff");
480  htfeffinpl->Divide(hmctfeffinpl);
481 
482  TH1F *htfeffinplprim = new TH1F("htfeffinplprim", "Track finding efficiency vs pl prim", (maxPl - minPl)*Pres, minPl, maxPl);
483  TH1F *hmctfeffinplprim = new TH1F("hmctfeffinplprim", "Track finding efficiency vs pl prim", (maxPl - minPl)*Pres, minPl, maxPl);
484  htfeffinplprim->Sumw2();
485  hmctfeffinplprim->Sumw2();
486  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> htfeffinplprim", cut_rec_prim, "goff");
487  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> hmctfeffinplprim", cut_mc_prim, "goff");
488  htfeffinplprim->Divide(hmctfeffinplprim);
489 
490  TH1F *htfeffinplsec = new TH1F("htfeffinplsec", "Track finding efficiency vs pl sec", (maxPl - minPl)*Pres, minPl, maxPl);
491  TH1F *hmctfeffinplsec = new TH1F("hmctfeffinplsec", "Track finding efficiency vs pl sec", (maxPl - minPl)*Pres, minPl, maxPl);
492  htfeffinplsec->Sumw2();
493  hmctfeffinplsec->Sumw2();
494  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> htfeffinplsec", cut_rec_sec, "goff");
495  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> hmctfeffinplsec", cut_mc_sec, "goff");
496  htfeffinplsec->Divide(hmctfeffinplsec);
497 
501  cut_rec = "MCTrackInfo.GetQuality() > 0";
502 
503  cut_mc_prim = cut_mc && "MCTrackInfo.GetMCQuality() == -2";
504  cut_mc_sec = cut_mc && "MCTrackInfo.GetMCQuality() == -1";
505  cut_rec_prim = cut_rec && cut_mc_prim;
506  cut_rec_sec = cut_rec && cut_mc_sec;
507 
508  TH1F* hstartZRange = new TH1F("hstartZRange","Range of StartZ", 2000, -1000, 1000);
509  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> hstartZRange", "", "goff");
510  Float_t maxStartZ = hstartZRange->GetBinLowEdge(hstartZRange->FindLastBinAbove(0)) + hstartZRange->GetBinWidth(hstartZRange->FindLastBinAbove(0));
511  Float_t minStartZ = hstartZRange->GetBinLowEdge(hstartZRange->FindFirstBinAbove(0)) + hstartZRange->GetBinWidth(hstartZRange->FindFirstBinAbove(0));
512  Float_t startZRes = 2;
513 
514  std::cout << "MinStartZ = " << minStartZ << " : MaxStartZ " << maxStartZ << " : StartZRes = " << startZRes << std::endl;
515 
516  TH1F* htfeffinstartZ = new TH1F("htfeffinstartZ", "Track finding efficiency vs startZ", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
517  TH1F* hmctfeffinstartZ = new TH1F("hmctfeffinstartZ", "Track finding efficiency vs startZ", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
518  htfeffinstartZ->Sumw2();
519  hmctfeffinstartZ->Sumw2();
520  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> htfeffinstartZ", cut_rec, "goff");
521  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> hmctfeffinstartZ", cut_mc, "goff");
522  htfeffinstartZ->Divide(hmctfeffinstartZ);
523 
524  TH1F *htfeffinstartZprim = new TH1F("htfeffinstartZprim", "Track finding efficiency vs startZ prim", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
525  TH1F *hmctfeffinstartZprim = new TH1F("hmctfeffinstartZprim", "Track finding efficiency vs startZ prim", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
526  htfeffinstartZprim->Sumw2();
527  hmctfeffinstartZprim->Sumw2();
528  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> htfeffinstartZprim", cut_rec_prim, "goff");
529  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> hmctfeffinstartZprim", cut_mc_prim, "goff");
530  htfeffinstartZprim->Divide(hmctfeffinstartZprim);
531 
532  TH1F *htfeffinstartZsec = new TH1F("htfeffinstartZsec", "Track finding efficiency vs startZ sec", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
533  TH1F *hmctfeffinstartZsec = new TH1F("hmctfeffinstartZsec", "Track finding efficiency vs startZ sec", (maxStartZ - minStartZ)*startZRes, minStartZ, maxStartZ);
534  htfeffinstartZsec->Sumw2();
535  hmctfeffinstartZsec->Sumw2();
536  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> htfeffinstartZsec", cut_rec_sec, "goff");
537  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].fStartZ >> hmctfeffinstartZsec", cut_mc_sec, "goff");
538  htfeffinstartZsec->Divide(hmctfeffinstartZsec);
539 
544  TFile outfile(outputFile.Data(), "RECREATE");
545  eff_mvdpix->Write();
546  eff_mvdstr->Write();
547  eff_stt->Write();
548  eff_gem->Write();
549  eff_fts->Write();
550  eff_glo->Write();
551 
552  pur_mvdpix->Write();
553  pur_mvdstr->Write();
554  pur_stt->Write();
555  pur_gem->Write();
556  pur_fts->Write();
557  pur_glo->Write();
558 
559  hnofrecotracks->Write();
560  hnofMCtracks->Write();
561 
562  hdelta_p_first->Write();
563  hdelta_pz_first->Write();
564  hdelta_pt_first->Write();
565  hdelta_theta_first->Write();
566  hdelta_phi_first->Write();
567  hdelta_x_first->Write();
568  hdelta_y_first->Write();
569  hdelta_z_first->Write();
570  hdelta_r_first->Write();
571 
572 
573  hdelta_p_last->Write();
574  hdelta_pz_last->Write();
575  hdelta_pt_last->Write();
576  hdelta_theta_last->Write();
577  hdelta_phi_last->Write();
578  hdelta_x_last->Write();
579  hdelta_y_last->Write();
580  hdelta_z_last->Write();
581  hdelta_r_last->Write();
582 
583  hdelta_charge->Write();
584 
585  heffintheta->Write();
586  heffinphi->Write();
587  heffinmom->Write();
588  heffinpt->Write();
589  heffinptprim->Write();
590  heffinptsec->Write();
591  heffinpl->Write();
592 
593  htfeffintheta->Write();
594  htfeffinthetaprim->Write();
595  htfeffinthetasec->Write();
596 
597  htfeffinphi->Write();
598  htfeffinphiprim->Write();
599  htfeffinphisec->Write();
600 
601  htfeffinpt->Write();
602  htfeffinptprim->Write();
603  htfeffinptsec->Write();
604 
605  htfeffinpl->Write();
606  htfeffinplprim->Write();
607  htfeffinplsec->Write();
608 
609  htfeffinstartZ->Write();
610  htfeffinstartZprim->Write();
611  htfeffinstartZsec->Write();
612 
613  hQualityPrim->Write();
614  hQualitySec->Write();
615 
616 
617  cout << " Test passed" << endl;
618  cout << " All ok " << endl;
619  return 0;
620 }
TTree * simtree
Definition: runPULL1.C:21
Double_t xmax
double cut[MAX]
Definition: autocutx.C:36
TString simFile
Definition: bump_emc.C:11
A simple class which adds the corresponding file extensions to a given base class.
PndMvdCreateDefaultApvMap * creator
void LabelQualyHistogram(TH1 *hist)
Double_t xmin
TH1F * hist
int QA_histos()
TString outfile