FairRoot/PandaRoot
Functions
histos_ca.C File Reference

Go to the source code of this file.

Functions

int histos_ca ()
 

Function Documentation

int histos_ca ( )

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 histos_ca.C.

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

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