FairRoot/PandaRoot
Functions
outdated/run/QA_histos.C File Reference

Go to the source code of this file.

Functions

int QA_histos ()
 

Function Documentation

int QA_histos ( )

7 histos: Global efficiency (#true hits/#totaltruehits), for all the primary tracks, for all the detectors: Efficiency for MVD pixel Eff for MVD strips Eff for Stt parallel Eff for Stt skewed Eff for gem Eff for FTS (skewed? parallel?)

7 histos: Global Purity (#true hits/hits), for all the primary tracks, for all the detectors: Purity for MVD pixel Pur for MVD strips Pur for Stt parallel Pur for Stt skewed Pur for gem Pur for FTS (skewed? parallel?)

2 histos: Number of reco tracks correlated to MC track, for primaries (>1 -> clones) Number of MC tracks correlated to reco tracks for primaries (>1 -> broken tracks

7 histos: Resolution at the first point: delta_p delta_pz delta_pperp delta_theta delta_phi delta_position

7 histos: Resolution at the last point: delta_p delta_pz delta_pperp delta_theta delta_phi delta_position

delta_charge

efficiency vs theta

efficiency vs phi

efficiency vs mom

efficiency vs pt

efficiency vs pl

Write to output

Definition at line 10 of file outdated/run/QA_histos.C.

References cut, outfile, simtree, xmax, and xmin.

10  {
11 
12  TFile fileqa("recoqa_complete.root");
13  TTree *simtree = (TTree*) fileqa.Get("pndsim");
14  simtree->AddFriend("simtree", "sim_complete.root");
15 
16  TCut cut = "";
17 
28  double effmin = -0.1;
29  double effmax = 1.1;
30 
31  TH1F *eff_mvdpix = new TH1F("eff_mvdpix", "MVD pixel efficiency", 100, effmin, effmax);
32  TH1F *eff_mvdstr = new TH1F("eff_mvdstr", "MVD strip efficiency", 100, effmin, effmax);
33  TH1F *eff_stt = new TH1F("eff_stt", "STT efficiency", 100, effmin, effmax);
34  TH1F *eff_gem = new TH1F("eff_gem", "GEM efficiency", 100, effmin, effmax);
35  TH1F *eff_fts = new TH1F("eff_fts", "FTS efficiency", 100, effmin, effmax);
36 
37  simtree->Draw("RecoTrackInfo.GetMvdPixelEfficiency() >> eff_mvdpix", cut, "goff");
38  simtree->Draw("RecoTrackInfo.GetMvdStripEfficiency() >> eff_mvdstr", cut, "goff");
39  simtree->Draw("RecoTrackInfo.GetSttEfficiency() >> eff_stt", cut, "goff");
40  simtree->Draw("RecoTrackInfo.GetGemEfficiency() >> eff_gem", cut, "goff");
41  simtree->Draw("RecoTrackInfo.GetFtsEfficiency() >> eff_fts", cut, "goff");
42 
43  TH1F *eff_glo = new TH1F("eff_glo", "global efficiency", 100, effmin, effmax);
44  simtree->Draw("RecoTrackInfo.GetEfficiency() >> eff_glo", cut, "goff");
45 
55  double purmin = -0.1;
56  double purmax = 1.1;
57 
58  TH1F *pur_mvdpix = new TH1F("pur_mvdpix", "MVD pixel purity", 100, purmin, purmax);
59  TH1F *pur_mvdstr = new TH1F("pur_mvdstr", "MVD strip purity", 100, purmin, purmax);
60  TH1F *pur_stt = new TH1F("pur_stt", "STT purity", 100, purmin, purmax);
61  TH1F *pur_gem = new TH1F("pur_gem", "GEM purity", 100, purmin, purmax);
62  TH1F *pur_fts = new TH1F("pur_fts", "FTS purity", 100, purmin, purmax);
63 
64  simtree->Draw("RecoTrackInfo.GetMvdPixelPurity() >> pur_mvdpix", cut, "goff");
65  simtree->Draw("RecoTrackInfo.GetMvdStripPurity() >> pur_mvdstr", cut, "goff");
66  simtree->Draw("RecoTrackInfo.GetSttPurity() >> pur_stt", cut, "goff");
67  simtree->Draw("RecoTrackInfo.GetGemPurity() >> pur_gem", cut, "goff");
68  simtree->Draw("RecoTrackInfo.GetFtsPurity() >> pur_fts", cut, "goff");
69 
70 
71  TH1F *pur_glo = new TH1F("pur_glo", "global purity", 100, purmin, purmax);
72  simtree->Draw("RecoTrackInfo.GetPurity() >> pur_glo", cut, "goff");
73 
79  TH1F *hnofrecotracks = new TH1F("hnofrecotracks", "# of reco tracks associated to the same MC track", 20, 0, 20);
80  simtree->Draw("MCTrackInfo.GetNofRecoTracks() >> hnofrecotracks", cut, "goff");
81 
82  TH1F *hnofMCtracks = new TH1F("hnofMCtracks", "# of MC tracks associated to the same reco track", 20, 0, 20);
83  simtree->Draw("RecoTrackInfo.GetNofMCTracks() >> hnofMCtracks", cut, "goff");
84 
94  double pmin = -1;
95  double pmax = 1;
96  TH1F *hdelta_p_first = new TH1F("hdelta_p_first", "#Delta p @ first hit", 100, pmin, pmax);
97  TH1F *hdelta_pz_first = new TH1F("hdelta_pz_first", "#Delta pz @ first hit", 100, pmin, pmax);
98  TH1F *hdelta_pt_first = new TH1F("hdelta_pt_first", "#Delta pt @ first hit", 100, pmin, pmax);
99 
100  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Mag() - RecoTrackInfo.GetMCTrackInfo().GetMomentumFirst().Mag() >> hdelta_p_first", cut, "goff");
101  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Perp() - RecoTrackInfo.GetMCTrackInfo().GetMomentumFirst().Perp() >> hdelta_pt_first", cut, "goff");
102  simtree->Draw("RecoTrackInfo.GetMomentumFirst().Z() - RecoTrackInfo.GetMCTrackInfo().GetMomentumFirst().Z() >> hdelta_pz_first", cut, "goff");
103 
104  double thetamin = -5;
105  double thetamax = 5;
106  int ntheta = (thetamax - thetamin)*10;
107  TH1F *hdelta_theta_first = new TH1F("hdelta_theta_first", "#Delta #theta @ first hit", ntheta, thetamin, thetamax);
108  simtree->Draw("(RecoTrackInfo.GetMomentumFirst().Theta() - RecoTrackInfo.GetMCTrackInfo().GetMomentumFirst().Theta()) * TMath::RadToDeg() >> hdelta_theta_first", cut, "goff");
109 
110  double phimin = -5;
111  double phimax = 5;
112  int nphi = (phimax - phimin)*10;
113  TH1F *hdelta_phi_first = new TH1F("hdelta_phi_first", "#Delta #phi @ first hit", nphi, phimin, phimax);
114  simtree->Draw("(RecoTrackInfo.GetMomentumFirst().Phi() - RecoTrackInfo.GetMCTrackInfo().GetMomentumFirst().Phi()) * TMath::RadToDeg() >> hdelta_phi_first", cut, "goff");
115 
116 
117  double xmin = -0.2;
118  double xmax = 0.2;
119  TH1F *hdelta_x_first = new TH1F("hdelta_x_first", "#Delta x @ first hit", 100, xmin, xmax);
120  TH1F *hdelta_y_first = new TH1F("hdelta_y_first", "#Delta y @ first hit", 100, xmin, xmax);
121  simtree->Draw("RecoTrackInfo.GetPositionFirst().X() - RecoTrackInfo.GetMCTrackInfo().GetPositionFirst().X() >> hdelta_x_first", cut, "goff");
122  simtree->Draw("RecoTrackInfo.GetPositionFirst().Y() - RecoTrackInfo.GetMCTrackInfo().GetPositionFirst().Y() >> hdelta_y_first", cut, "goff");
123 
124  double zmin = -0.2;
125  double zmax = 0.2;
126  TH1F *hdelta_z_first = new TH1F("hdelta_z_first", "#Delta z @ first hit", 100, zmin, zmax);
127  simtree->Draw("RecoTrackInfo.GetPositionFirst().Z() - RecoTrackInfo.GetMCTrackInfo().GetPositionFirst().Z() >> hdelta_z_first", cut, "goff");
128 
129  double rmax = xmax;
130  TH1F *hdelta_r_first = new TH1F("hdelta_r_first", "#Delta r @ first hit", 100, -rmax, rmax);
131  simtree->Draw("RecoTrackInfo.GetPositionFirst().Perp() - RecoTrackInfo.GetMCTrackInfo().GetPositionFirst().Perp() >> hdelta_r_first", cut, "goff");
132 
142  TH1F *hdelta_p_last = new TH1F("hdelta_p_last", "#Delta p @ last hit", 100, pmin, pmax);
143  TH1F *hdelta_pz_last = new TH1F("hdelta_pz_last", "#Delta pz @ last hit", 100, pmin, pmax);
144  TH1F *hdelta_pt_last = new TH1F("hdelta_pt_last", "#Delta pt @ last hit", 100, pmin, pmax);
145 
146  simtree->Draw("RecoTrackInfo.GetMomentumLast().Mag() - RecoTrackInfo.GetMCTrackInfo().GetMomentumLast().Mag() >> hdelta_p_last", cut, "goff");
147  simtree->Draw("RecoTrackInfo.GetMomentumLast().Perp() - RecoTrackInfo.GetMCTrackInfo().GetMomentumLast().Perp() >> hdelta_pt_last", cut, "goff");
148  simtree->Draw("RecoTrackInfo.GetMomentumLast().Z() - RecoTrackInfo.GetMCTrackInfo().GetMomentumLast().Z() >> hdelta_pz_last", cut, "goff");
149 
150  TH1F *hdelta_theta_last = new TH1F("hdelta_theta_last", "#Delta #theta @ last hit", ntheta, thetamin, thetamax);
151  simtree->Draw("(RecoTrackInfo.GetMomentumLast().Theta() - RecoTrackInfo.GetMCTrackInfo().GetMomentumLast().Theta()) * TMath::RadToDeg() >> hdelta_theta_last", cut, "goff");
152 
153  TH1F *hdelta_phi_last = new TH1F("hdelta_phi_last", "#Delta #phi @ last hit", nphi, phimin, phimax);
154  simtree->Draw("(RecoTrackInfo.GetMomentumLast().Phi() - RecoTrackInfo.GetMCTrackInfo().GetMomentumLast().Phi()) * TMath::RadToDeg() >> hdelta_phi_last", cut, "goff");
155 
156  TH1F *hdelta_x_last = new TH1F("hdelta_x_last", "#Delta x @ last hit", 100, xmin, xmax);
157  TH1F *hdelta_y_last = new TH1F("hdelta_y_last", "#Delta y @ last hit", 100, xmin, xmax);
158  simtree->Draw("RecoTrackInfo.GetPositionLast().X() - RecoTrackInfo.GetMCTrackInfo().GetPositionLast().X() >> hdelta_x_last", cut, "goff");
159  simtree->Draw("RecoTrackInfo.GetPositionLast().Y() - RecoTrackInfo.GetMCTrackInfo().GetPositionLast().Y() >> hdelta_y_last", cut, "goff");
160 
161  double zminlast = -0.6;
162  double zmaxlast = 0.6;
163  TH1F *hdelta_z_last = new TH1F("hdelta_z_last", "#Delta z @ last hit", 100, zminlast, zmaxlast);
164  simtree->Draw("RecoTrackInfo.GetPositionLast().Z() - RecoTrackInfo.GetMCTrackInfo().GetPositionLast().Z() >> hdelta_z_last", cut, "goff");
165 
166  TH1F *hdelta_r_last = new TH1F("hdelta_r_last", "#Delta r @ last hit", 100, -rmax, rmax);
167  simtree->Draw("RecoTrackInfo.GetPositionLast().Perp() - RecoTrackInfo.GetMCTrackInfo().GetPositionLast().Perp() >> hdelta_r_last", cut, "goff");
168 
172  double chargemin = -3;
173  double chargemax = 3;
174  int ncharge = chargemax - chargemin;
175  TH1F *hdelta_charge = new TH1F("hdelta_charge", "#Delta #charge", ncharge, chargemin, chargemax);
176  simtree->Draw("RecoTrackInfo.GetCharge() - RecoTrackInfo.GetMCTrackInfo().GetCharge() >> hdelta_charge", cut, "goff");
177 
181  TH1F *hthetagen = new TH1F("hthetagen", "mc theta dist", 180, 0, 180);
182  TH1F *heffintheta = new TH1F("heffintheta", "efficiency vs #theta", 180, 0, 180);
183 
184  TCut cut_mc = "";
185  TCut cut_rec = cut_mc && "MCTrackInfo.GetRecoTrackID() != -1";
186  cut_rec = cut_rec && "RecoTrackInfo[MCTrackInfo.GetRecoTrackID()].GetEfficiency() > 0.8";
187  cut_rec = cut_rec && "RecoTrackInfo[MCTrackInfo.GetRecoTrackID()].IsClone() == 0";
188 
189  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> hthetagen", cut_mc, "goff");
190  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Theta() * TMath::RadToDeg() >> heffintheta", cut_rec, "goff");
191 
192  hthetagen->Sumw2();
193  heffintheta->Sumw2();
194  heffintheta->Divide(hthetagen);
195 
196 
200  TH1F *hphigen = new TH1F("hphigen", "mc phi dist", 180, -180, 180);
201  TH1F *heffinphi = new TH1F("heffinphi", "efficiency vs #phi", 180, -180, 180);
202  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> hphigen", cut_mc, "goff");
203  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Phi() * TMath::RadToDeg() >> heffinphi", cut_rec, "goff");
204 
205  hphigen->Sumw2();
206  heffinphi->Sumw2();
207  heffinphi->Divide(hphigen);
208 
212  TH1F *hmomgen = new TH1F("hmomgen", "mc mom dist", 100, 0, 3);
213  TH1F *heffinmom = new TH1F("heffinmom", "efficiency vs mom", 100, 0, 3);
214  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Mag() >> hmomgen", cut_mc, "goff");
215  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Mag() >> heffinmom", cut_rec, "goff");
216 
217  hmomgen->Sumw2();
218  heffinmom->Sumw2();
219  heffinmom->Divide(hmomgen);
220 
224  TH1F *hptgen = new TH1F("hptgen", "mc pt dist", 100, 0, 3);
225  TH1F *heffinpt = new TH1F("heffinpt", "efficiency vs pt", 100, 0, 3);
226  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> hptgen", cut_mc, "goff");
227  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Perp() >> heffinpt", cut_rec, "goff");
228 
229  hptgen->Sumw2();
230  heffinpt->Sumw2();
231  heffinpt->Divide(hptgen);
232 
236  TH1F *hplgen = new TH1F("hplgen", "mc pl dist", 100, 0, 3);
237  TH1F *heffinpl = new TH1F("heffinpl", "efficiency vs pl", 100, 0, 3);
238  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> hplgen", cut_mc, "goff");
239  simtree->Draw("MCTrack[MCTrackInfo.GetMCTrackID()].GetMomentum().Z() >> heffinpl", cut_rec, "goff");
240 
241  hplgen->Sumw2();
242  heffinpl->Sumw2();
243  heffinpl->Divide(hplgen);
244 
245 
246 
251  TFile outfile("QA_histograms.root", "RECREATE");
252  eff_mvdpix->Write();
253  eff_mvdstr->Write();
254  eff_stt->Write();
255  eff_gem->Write();
256  eff_fts->Write();
257  eff_glo->Write();
258 
259  pur_mvdpix->Write();
260  pur_mvdstr->Write();
261  pur_stt->Write();
262  pur_gem->Write();
263  pur_fts->Write();
264  pur_glo->Write();
265 
266  hnofrecotracks->Write();
267  hnofMCtracks->Write();
268 
269  hdelta_p_first->Write();
270  hdelta_pz_first->Write();
271  hdelta_pt_first->Write();
272  hdelta_theta_first->Write();
273  hdelta_phi_first->Write();
274  hdelta_x_first->Write();
275  hdelta_y_first->Write();
276  hdelta_z_first->Write();
277  hdelta_r_first->Write();
278 
279 
280  hdelta_p_last->Write();
281  hdelta_pz_last->Write();
282  hdelta_pt_last->Write();
283  hdelta_theta_last->Write();
284  hdelta_phi_last->Write();
285  hdelta_x_last->Write();
286  hdelta_y_last->Write();
287  hdelta_z_last->Write();
288  hdelta_r_last->Write();
289 
290  hdelta_charge->Write();
291 
292  heffintheta->Write();
293  heffinphi->Write();
294  heffinmom->Write();
295  heffinpt->Write();
296  heffinpl->Write();
297 
298  cout << " Test passed" << endl;
299  cout << " All ok " << endl;
300  return 0;
301 }
TTree * simtree
Definition: runPULL1.C:21
Double_t xmax
double cut[MAX]
Definition: autocutx.C:36
Double_t xmin
TString outfile