FairRoot/PandaRoot
Resolution_and_acceptance.cxx
Go to the documentation of this file.
1 /*
2  * test_main.cxx
3  *
4  * Created on: Feb 8, 2012
5  * Author: promme
6  */
7 
8 #include <iostream>
9 #include <sstream>
10 
11 #include<TApplication.h>
12 #include<TCanvas.h>
13 #include<TROOT.h>
14 #include<TString.h>
15 #include<TChain.h>
16 #include<TFile.h>
17 #include<TClonesArray.h>
18 #include<TSystem.h>
19 #include<TH1.h>
20 #include<TH2.h>
21 #include<TRotation.h>
22 #include<TVector3.h>
23 #include<TMath.h>
24 #include<TGaxis.h>
25 
26 #include<TStopwatch.h>
27 
28 #include<PndMCTrack.h>
29 #include<PndSdsMCPoint.h>
30 
31 // needed for geane backtracking
32 #include<FairRunAna.h>
33 #include<FairRootManager.h>
34 #include<FairGeane.h>
35 #include<FairRtdbRun.h>
36 #include<FairRuntimeDb.h>
37 #include<FairParRootFileIo.h>
39 #include<FairLogger.h>
40 #include<string>
41 
42 using namespace std;
43 
44 int main(int nargs, char** args) {
45  //TApplication myapp("myapp",0,0);
46  TString storePath = "./tmpOutput";
47  // set the storepath to somewhere else, if provided
48  if (nargs == 2) {
49  storePath = args[1];
50  } else {
51  if (nargs == 1)
52  cout << " Hint: you may change the standard path from "
53  << storePath.Data()
54  << " to any other folder by providing the path as an argument! "
55  << endl;
56  if (nargs > 2)
57  cout << " Error: too many arguments! " << endl;
58  }
59  TString startEvent = "0";
60  // ========================================================================
61  // Input file (MC events)
62  TString MCFile = storePath + "/Lumi_MC_";
63  MCFile += startEvent;
64  MCFile += ".root";
65  TString DigiFile = storePath + "/Lumi_digi_";
66  DigiFile += startEvent;
67  DigiFile += ".root";
68  // Digi file
69  TString RecoFile = storePath + "/Lumi_reco_";
70  RecoFile += startEvent;
71  RecoFile += ".root";
72  // TCand file
73  TString CandFile = storePath + "/Lumi_TCand_";
74  CandFile += startEvent;
75  CandFile += ".root";
76  // Parameter file << this one
77  TString parFile = storePath + "/Lumi_Params_";
78  parFile += startEvent;
79  parFile += ".root";
80  // Track file
81  TString TrkFile = storePath + "/Lumi_Track_";
82  TrkFile += startEvent;
83  TrkFile += ".root";
84 
85  // ---- Load libraries -------------------------------------------------
86  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
87  gSystem->Load("libSds");
88  gSystem->Load("libSdsReco");
89  gSystem->Load("libLmd");
90  gSystem->Load("libLmdReco");
91  gSystem->Load("libLmdTrk");
92  // ------------------------------------------------------------------------
93  // Output file
94  TString outFile = storePath + "/Lumi_Testing_";
95  outFile += startEvent;
96  outFile += ".root";
97 
98  // Output histogram file
99  TString outFilehist = storePath + "/Lumi_performance_histograms_";
100  outFilehist += startEvent;
101  outFilehist += ".root";
102 
103  std::cout << "MCFile : " << MCFile.Data() << std::endl;
104  std::cout << "DigiFile: " << DigiFile.Data() << std::endl;
105  std::cout << "RecoFile: " << RecoFile.Data() << std::endl;
106  std::cout << "TCandFile: " << CandFile.Data() << std::endl;
107  std::cout << "TrackFile: " << TrkFile.Data() << std::endl;
108  std::cout << "TestingFile: " << outFile.Data() << std::endl;
109  // --- Now choose concrete engines for the different tasks -------------
110  // ------------------------------------------------------------------------
111 
112 
113  // In general, the following parts need not be touched
114  // ========================================================================
115 
116 
117  // ----- Timer --------------------------------------------------------
118  TStopwatch timer;
119  timer.Start();
120  // ------------------------------------------------------------------------
121 
122 
123  // ----- Reconstruction run -------------------------------------------
124  FairRunAna *fRun = new FairRunAna();
125  fRun->SetInputFile(MCFile);
126  //fRun->AddFriend(DigiFile);
127  //fRun->AddFriend(RecoFile);
128  //fRun->AddFriend(CandFile);
129  //fRun->AddFriend(TrkFile);
130  fRun->SetOutputFile(outFile);
131  // ------------------------------------------------------------------------
132 
133  FairGeane *Geane = new FairGeane();
134  fRun->AddTask(Geane);
135 
136  // ----- Parameter database --------------------------------------------
137  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
138  FairParRootFileIo* parInput1 = new FairParRootFileIo(kTRUE);
139  parInput1->open(parFile.Data());
140  rtdb->setFirstInput(parInput1);
141 
142  // =========================================================================
143  // ====== Geane Back-Propagating ======
144  // =========================================================================
145 
146  // ----- SDS hit producer --------------------------------------------
147 
148  // PndEmcMapper *emcMap = PndEmcMapper::Instance(6);
149  Double_t fpBeam = 1.5;
150  TVector3 IP(0, 0, 0);
151  PndLmdPerformanceTask* lmdperformance = new PndLmdPerformanceTask(fpBeam,
152  IP);
153  lmdperformance->SetVerbose(0);
154  fRun->AddTask(lmdperformance);
155  rtdb->setOutput(parInput1);
156  lmdperformance->SetHistFilename(outFilehist);
157  rtdb->print();
158  // ===== End of Geane =====
159  // =========================================================================
160 
161  // ----- Intialise and run --------------------------------------------
162  fRun->Init();
163  // // PndEmcMapper *emcMap = PndEmcMapper::Instance(6);
164  // PndEmcMapper *emcMap = PndEmcMapper::Instance();
165  // //Geane->SetField(fRun->GetField());
166  FairRootManager* ioman = FairRootManager::Instance();
167  if (!ioman) {
168  std::cout << "Error: "
169  << "RootManager not instantiated!" << std::endl;
170  return 1;
171  }
172  fRun->Run(0, ioman->CheckMaxEventNo()); // hmm does not work yet, seems to return 0 but then is runs over all what was intended.
173  // ------------------------------------------------------------------------
174  rtdb->saveOutput();
175  rtdb->print();
176  // ----- Finish -------------------------------------------------------
177  timer.Stop();
178  Double_t rtime = timer.RealTime();
179  Double_t ctime = timer.CpuTime();
180  cout << endl << endl;
181  cout << "Macro finished succesfully." << endl;
182  cout << "Output file is " << outFile << endl;
183  cout << "Parameter file is " << parFile << endl;
184  cout << "Real time " << rtime << " s, CPU time " << ctime << " s" << endl;
185  cout << endl;
186  // ------------------------------------------------------------------------
187 
188  //myapp.Run();
189  return 0;
190 }
191 
192 /*
193  // some design constants of the LUMI detector
194  const int nplanes(4);
195  const int nsensors_per_plane(8);
196 
197  double last_percent(-1.);
198  // draw a progress bar only when the length changes significantly
199  void DrawProgressBar(int len, double percent);
200 
201  void Resolution_and_acceptance(){
202  cout << " analysis tool vers. 1.1 " << endl;
203 
204  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
205  gSystem->Load("libSds");
206  gSystem->Load("libSdsReco");
207  gSystem->Load("libLmd");
208  gSystem->Load("libLmdReco");
209  gSystem->Load("libLmdTrk");
210 
211  TStopwatch timer;
212  timer.Start();
213  // ------------------------------------------------------------------------
214 
215  TString storePath = "./tmpOutput";
216  TString startEvent = "0";
217  // ---- Input files --------------------------------------------------------
218  TString simMC = storePath + "/Lumi_MC_";
219  simMC += startEvent;
220  simMC += ".root";
221  TChain tMC("pndsim");
222  tMC.Add(simMC);
223 
224  TString DigiFile = storePath + "/Lumi_digi_";
225  DigiFile += startEvent;
226  DigiFile += ".root";
227  TChain tdigiHits("pndsim");
228  tdigiHits.Add(DigiFile);
229 
230  TString recHit = storePath + "/Lumi_reco_";
231  recHit += startEvent;
232  recHit += ".root";
233  TChain tHits("pndsim");
234  tHits.Add(recHit);
235 
236  TString trkCand = storePath + "/Lumi_TCand_";
237  trkCand += startEvent;
238  trkCand += ".root";
239  TChain tTrkCand("pndsim");
240  tTrkCand.Add(trkCand);
241 
242  TString recTrack = storePath + "/Lumi_Track_";
243  recTrack += startEvent;
244  recTrack += ".root";
245  TChain tTrkRec("pndsim");
246  tTrkRec.Add(recTrack);
247 
248  TString geaneFile = storePath + "/Lumi_Geane_";
249  geaneFile += startEvent;
250  geaneFile += ".root";
251  TChain tgeane("pndsim");
252  tgeane.Add(geaneFile);
253 
254  // Parameter file << needed for geane back tracking
255  TString parFile = storePath+"/Lumi_Params_";
256  parFile += startEvent;
257  parFile += ".root";
258 
259  // ---- Output file ----------------------------------------------------------------
260  TString fakefile = storePath + "/fake";
261  fakefile += startEvent;
262  fakefile += ".root";tmpOutput
263  // ---------------------------------------------------------------------------------
264 
265  // ---------------------------------------------------------------------------------
266 
267  // ---- needed to use the back propagation ----------------------------------------
268 
269  FairRunAna *fRun= new FairRunAna();
270  fRun->SetInputFile(simMC);
271  //fRun->AddFriend(DigiFile);
272  //fRun->AddFriend(recHit);
273  //fRun->AddFriend(trkCand);
274  //fRun->AddFriend(recTrack);
275  fRun->SetOutputFile(fakefile);
276  // ------------------------------------------------------------------------
277 
278  FairGeane *Geane = new FairGeane();
279  //fRun->AddTask(Geane);
280  FairGeanePro* fPro = Geane->Exec();
281  // new FairGeanePro();
282 
283  cout << " ok " << endl;
284 
285  // ----- Parameter database --------------------------------------------
286  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
287  FairParRootFileIo* parInput1 = new FairParRootFileIo(kTRUE);
288  parInput1->open(parFile.Data());
289  rtdb->setFirstInput(parInput1);
290  fRun->Init();
291 
292  // ---- Output file ----------------------------------------------------------------
293  TString out = storePath + "/Resolution_and_acceptance";
294  out += startEvent;
295  out += ".root";
296  TFile *output_histfile = new TFile(out, "RECREATE");
297  // ---------------------------------------------------------------------------------
298 
299  //--- MC info -----------------------------------------------------------------
300  TClonesArray* true_tracks = new TClonesArray("PndMCTrack");
301  tMC.SetBranchAddress("MCTrack", &true_tracks); //True Track to compare
302 
303  TClonesArray* true_points = new TClonesArray("PndSdsMCPoint");
304  tMC.SetBranchAddress("LMDPoint", &true_points); //True Points to compare
305  //----------------------------------------------------------------------------------
306 
307 
308  //--- Digitization info ------------------------------------------------------------
309  TClonesArray* fStripClusterArray = new TClonesArray("PndSdsClusterStrip");
310  tHits.SetBranchAddress("LMDStripClusterCand", &fStripClusterArray);
311 
312  TClonesArray* fStripDigiArray = new TClonesArray("PndSdsDigiStrip");
313  tdigiHits.SetBranchAddress("LMDStripDigis", &fStripDigiArray);
314  //----------------------------------------------------------------------------------
315 
316  //--- Real Hits --------------------------------------------------------------------
317  TClonesArray* rechit_array = new TClonesArray("PndSdsHit");
318  tHits.SetBranchAddress("LMDHitsStrip", &rechit_array); //Points for Tracks
319  //----------------------------------------------------------------------------------
320 
321 
322  //--- Track Candidate ---------------------------------------------------------------
323  TClonesArray* trkcand_array = new TClonesArray("PndTrackCand");
324  tTrkCand.SetBranchAddress("LMDTrackCand", &trkcand_array); //Points for Track Canidates
325  //-----------------------------------------------------------------------------------
326 
327  //--- Real tracks -------------------------------------------------------------------
328  TClonesArray* rec_trk = new TClonesArray("PndLinTrack");
329  tTrkRec.SetBranchAddress("LMDTrack", &rec_trk); //Tracks
330  //----------------------------------------------------------------------------------
331 
332  //--- Geane info ------------------------------------------------------------------
333  // TClonesArray* geaneArray =new TClonesArray("FairTrackParP");
334  //TClonesArray* geaneArray = new TClonesArray("FairTrackParH");
335  //tgeane.SetBranchAddress("GeaneTrackFinal", &geaneArray); //Tracks with parabolic parameterization
336  // angular acceptance in theta and phi
337 
338  TH2* hist_angular_distr_gen = new TH2F("hist_angular_distr_gen", "hist_angular_distr_gen", 100, 0., 20.e-3, 100, 0., 3.141);
339  TH2* hist_angular_distr_acc = new TH2F("hist_angular_distr_acc", "hist_angular_distr_acc", 100, 0., 20.e-3, 100, 0., 3.141);
340 
341  // spatial acceptance in x and y at the first lumi plane
342  TH2* hist_spatial_distr_gen = new TH2F("hist_spatial_distr_gen", "hist_spatial_distr_gen", 1000, -100., 100., 1000, -100., 100.);
343  TH2* hist_spatial_distr_acc = new TH2F("hist_spatial_distr_acc", "hist_spatial_distr_acc", 1000, -100., 100., 1000, -100., 100.);
344 
345  TTree* tree_results = new TTree("tree_results", "tree_results");
346  // initial momenta of the antiprotons
347  double px_init;
348  double py_init;
349  double pz_init;
350  double ptheta_init;
351  double pphi_init;
352  tree_results->Branch("px_init", &px_init);
353  tree_results->Branch("py_init", &py_init);
354  tree_results->Branch("pz_init", &pz_init);
355  tree_results->Branch("ptheta_init", &ptheta_init);
356  tree_results->Branch("pphi_init", &pphi_init);
357  // momenta at the entry of a detector plane
358  // in the reference system of the LUMI
359  double px_in;
360  double py_in;
361  double pz_in;
362  double ptheta_in;
363  double pphi_in;
364  tree_results->Branch("px_in", &px_in);
365  tree_results->Branch("py_in", &py_in);
366  tree_results->Branch("pz_in", &pz_in);
367  tree_results->Branch("ptheta_in", &ptheta_in);
368  tree_results->Branch("pphi_in", &pphi_in);
369  // position at the entry of a detector plane
370  // in the reference system of the LUMI
371  double x_in;
372  double y_in;
373  double z_in;
374  tree_results->Branch("x_in", &x_in);
375  tree_results->Branch("y_in", &y_in);
376  tree_results->Branch("z_in", &z_in);
377  // detector and sensor id of the mc hits
378  int det_id;
379  int sens_id;
380  tree_results->Branch("det_id", &det_id);
381  tree_results->Branch("sens_id", &sens_id);
382  // calculated plane and sensor in plane id
383  int plane;
384  int sensor;
385  tree_results->Branch("plane", &plane);
386  tree_results->Branch("sensor", &sensor);
387 
388  // Drawing directly from a tree is elegant but slow as every Draw call
389  // loops over the whole tree of events
390  // Better it is to create the desired histograms in advance
391  // and to fill those on event by event basis simultaneously
392 
393  // histograms per plane
394  TH2* hist_xy [nplanes];
395  TH1* hist_theta_init [nplanes];
396  TH1* hist_theta_in [nplanes];
397  TH1* hist_theta_diff [nplanes];
398  TH1* hist_theta_diff_rel [nplanes];
399 
400  // histograms per sensor
401  TH2* hists_xy [nplanes][nsensors_per_plane];
402  TH1* hists_theta_init [nplanes][nsensors_per_plane];
403  TH1* hists_theta_in [nplanes][nsensors_per_plane];
404  TH1* hists_theta_diff [nplanes][nsensors_per_plane];
405  TH1* hists_theta_diff_rel [nplanes][nsensors_per_plane];
406 
407  TCanvas temp_canvas("temp_canvas", "canvas for initialization", 100, 100);
408  temp_canvas.cd();
409  // loop over planes
410  cout << " constructing histograms for " << nplanes << " planes with " << nsensors_per_plane << " sensors per plane " << endl;
411 
412  for (int iplane = 0; iplane < nplanes; iplane++){
413  stringstream hist_name;
414  stringstream hist_title;
415 
416  hist_name << "hist_xy_plane_" << iplane;
417  hist_title << "xy hit distribution plane " << iplane;
418  hist_xy[iplane] = new TH2F(hist_name.str().c_str(), hist_title.str().c_str(),
419  100, -10, 10, 100, -10, 10);
420  hist_xy[iplane]->Draw();
421  hist_xy[iplane]->GetXaxis()->SetTitle("X [cm]");
422  hist_xy[iplane]->GetYaxis()->SetTitle("Y [cm]");
423 
424  hist_name .str("");
425  hist_title.str("");
426 
427  hist_name << "hist_theta_init_plane_" << iplane;
428  hist_title << "initial #Theta distribution plane " << iplane;
429  hist_theta_init[iplane] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
430  150, 0, 15e-3);
431  hist_theta_init[iplane]->Draw();
432  hist_theta_init[iplane]->GetXaxis()->SetTitle("#Theta [rad]");
433  hist_theta_init[iplane]->GetYaxis()->SetTitle("entries");
434 
435  hist_name .str("");
436  hist_title.str("");
437 
438  hist_name << "hist_theta_in_plane_" << iplane;
439  hist_title << "#Theta distribution at plane " << iplane;
440  hist_theta_in[iplane] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
441  200, 0, 20e-3);
442  hist_theta_in[iplane]->Draw();
443  hist_theta_in[iplane]->GetXaxis()->SetTitle("#Theta [rad]");
444  hist_theta_in[iplane]->GetYaxis()->SetTitle("entries");
445 
446  hist_name .str("");
447  hist_title.str("");
448 
449  hist_name << "hist_theta_diff_in_plane_" << iplane;
450  hist_title << "#Delta#Theta distribution at plane " << iplane;
451  hist_theta_diff[iplane] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
452  100, -20e-5, 20e-5);
453  hist_theta_diff[iplane]->Draw();
454  hist_theta_diff[iplane]->GetXaxis()->SetTitle("#Delta#Theta [rad]");
455  hist_theta_diff[iplane]->GetYaxis()->SetTitle("entries");
456 
457  hist_name .str("");
458  hist_title.str("");
459 
460  hist_name << "hist_theta_diff_rel_in_plane_" << iplane;
461  hist_title << "#Delta#Theta/#Theta distribution at plane " << iplane;
462  hist_theta_diff_rel[iplane] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
463  100, -1, 1);
464  hist_theta_diff_rel[iplane]->Draw();
465  hist_theta_diff_rel[iplane]->GetXaxis()->SetTitle("#Delta#Theta/#Theta");
466  hist_theta_diff_rel[iplane]->GetYaxis()->SetTitle("entries");
467 
468  // loop over sensors per plane which are enumerated linearly
469  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
470  hist_name .str("");
471  hist_title.str("");
472 
473  hist_name << "hist_xy_plane_" << iplane << "_sensor_" << isensor;
474  hist_title << "xy hit distribution plane " << iplane << " sensor " << isensor;
475  hists_xy[iplane][isensor] = new TH2F(hist_name.str().c_str(), hist_title.str().c_str(),
476  100, -10, 10, 100, -10, 10);
477  hists_xy[iplane][isensor]->Draw();
478  hists_xy[iplane][isensor]->GetXaxis()->SetTitle("X [cm]");
479  hists_xy[iplane][isensor]->GetYaxis()->SetTitle("Y [cm]");
480 
481  hist_name .str("");
482  hist_title.str("");
483 
484  hist_name << "hist_theta_init_plane_" << iplane << "_sensor_" << isensor;
485  hist_title << "initial #Theta distribution plane " << iplane << " sensor " << isensor;
486  hists_theta_init[iplane][isensor] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
487  150, 0, 15e-3);
488  hists_theta_init[iplane][isensor]->Draw();
489  hists_theta_init[iplane][isensor]->GetXaxis()->SetTitle("#Theta [rad]");
490  hists_theta_init[iplane][isensor]->GetYaxis()->SetTitle("entries");
491 
492  hist_name .str("");
493  hist_title.str("");
494 
495  hist_name << "hist_theta_in_plane_" << iplane << "_sensor_" << isensor;
496  hist_title << "#Theta distribution at plane " << iplane << " sensor " << isensor;
497  hists_theta_in[iplane][isensor] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
498  200, 0, 20e-3);
499  hists_theta_in[iplane][isensor]->Draw();
500  hists_theta_in[iplane][isensor]->GetXaxis()->SetTitle("#Theta [rad]");
501  hists_theta_in[iplane][isensor]->GetYaxis()->SetTitle("entries");
502 
503  hist_name .str("");
504  hist_title.str("");
505 
506  hist_name << "hist_theta_diff_in_plane_" << iplane << "_sensor_" << isensor;
507  hist_title << "#Delta#Theta distribution at plane " << iplane << " sensor " << isensor;
508  hists_theta_diff[iplane][isensor] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
509  100, -20e-5, 20e-5);
510  hists_theta_diff[iplane][isensor]->Draw();
511  hists_theta_diff[iplane][isensor]->GetXaxis()->SetTitle("#Delta#Theta [rad]");
512  hists_theta_diff[iplane][isensor]->GetYaxis()->SetTitle("entries");
513 
514  hist_name .str("");
515  hist_title.str("");
516 
517  hist_name << "hist_theta_diff_rel_in_plane_" << iplane << "_sensor_" << isensor;
518  hist_title << "#Delta#Theta/#Theta distribution at plane " << iplane << " sensor " << isensor;
519  hists_theta_diff_rel[iplane][isensor] = new TH1F(hist_name.str().c_str(), hist_title.str().c_str(),
520  100, -1, 1);
521  hists_theta_diff_rel[iplane][isensor]->Draw();
522  hists_theta_diff_rel[iplane][isensor]->GetXaxis()->SetTitle("#Delta#Theta/#Theta");
523  hists_theta_diff_rel[iplane][isensor]->GetYaxis()->SetTitle("entries");
524  }
525  }
526 
527 
528  // rotation into the lumi position according to the beamline file
529  // matrix tri6 - tr=1 rot=1 refl=0 scl=0
530  // 0.999192 0.000000 0.040200 Tx = 22.893063
531  // 0.000000 1.000000 -0.000000 Ty = 0.000000
532  // -0.040200 0.000000 0.999192 Tz = 1049.770390
533  // The matrix is a rotation around y of around 0.04025 radian
534  TRotation inv_lmdrotation;
535  //inv_lmdrotation.RotateY(-0.04025);
536  inv_lmdrotation.RotateY(-2.326/180.*TMath::Pi()); // as derived by mathias michel
537  TVector3 inv_lmdtranslation(
538  -22.893063, 0.000000, -1049.770390);
539 
540  int nEvents = tMC.GetEntries();
541  if (tTrkCand.GetEntries() != nEvents) {cout << "nMC <> nTrkCand! Exiting." << endl; return;}
542  if (tTrkRec.GetEntries() != nEvents) {cout << "nMC <> nTrkRec! Exiting." << endl; return;}
543  if (tHits.GetEntries() != nEvents) {cout << "nMC <> nHits! Exiting." << endl; return;}
544  if (tdigiHits.GetEntries() != nEvents) {cout << "nMC <> nTrkdigiHits! Exiting." << endl; return;}
545  //if (tgeane.GetEntries() != nEvents) {cout << "nMC <> TrkCand! Exiting." << endl; return;}
546  //TVector3 ip(0,0,0);
547  //PndLmdGeaneTask* lmdgeane = new PndLmdGeaneTask();
548  cout << " reading " << nEvents << " Events " << endl;
549  for (Int_t j = 0; j < nEvents; j++) {
550  DrawProgressBar(50, (j+1)/((double)nEvents));
551  tMC.GetEntry(j);
552  tTrkCand.GetEntry(j);
553  tTrkRec.GetEntry(j);
554  tHits.GetEntry(j);
555  tdigiHits.GetEntry(j);
556  const int nParticles = true_tracks->GetEntriesFast();
557 
558  int npbar = 0;
559  for (Int_t iParticle = 0; iParticle < nParticles; iParticle++) {
560  // if(nParticles!=4) continue;
561  PndMCTrack *mctrk = (PndMCTrack*) true_tracks->At(iParticle);
562  Int_t mcID = mctrk->GetPdgCode();
563  if (mctrk->IsGeneratorCreated()){//mcID == -2212) {
564  npbar++;
565  // if(nParticles!=4) cout<<"For event #"<<j<<" mcID="<<mcID<<endl;
566  TVector3 MomMC = mctrk->GetMomentum();
567  hist_angular_distr_gen->Fill(MomMC.Theta(), MomMC.Phi());
568  //hallTheta->Fill(MomMC.Theta());
569  //hallPhi->Fill(MomMC.Phi());
570  }
571  }
572  if (0&&npbar != 1)
573  cout << "found " << npbar << " anti-protons" << endl;
574 
575  int nHits = true_points->GetEntriesFast();
576  int nSdsHits = 0;
577  for (Int_t iHit = 0; iHit < nHits; iHit++) {
578  PndSdsMCPoint* mcpoint = (PndSdsMCPoint*) true_points->At(iHit);
579  PndMCTrack *mctrk = (PndMCTrack*) true_tracks->At(mcpoint->GetTrackID());
580  if (mctrk->IsGeneratorCreated()){ // take only anti-protons from the generator
581  if (mcpoint){
582  TVector3 _mcpoint((mcpoint->GetPosition())); // position at the sensor entrance
583  TVector3 _mctrack(mcpoint->GetPx(), mcpoint->GetPy(), mcpoint->GetPz()); // momentum of the track at the entrance
584  if (1) { // translate it into the reference system of the lumi monitor
585  _mcpoint = inv_lmdtranslation+_mcpoint;
586  _mcpoint = inv_lmdrotation*_mcpoint;
587  }
588  TVector3 momMC = mctrk->GetMomentum(); // momentum in the primary vertex
589  TVector3 posMC = mctrk->GetStartVertex(); // position of the primary vertex
590  if (1) { // use back propagation for the track
591  TVector3 initpos, initmom(0,0,momMC.Mag()),
592  StartPosErr, StartMomErr, \
593  ResPosErr, ResMomErr;
594  int charge = mctrk->GetPdgCode() > 0 ? 1 : mctrk->GetPdgCode() < 0 ? -1 : 0;
595  FairTrackParH *fStart =
596  new FairTrackParH(posMC, momMC, StartPosErr, StartMomErr, charge); // exchange fCharge by charge of PDGCode
597  FairTrackParH *fRes =
598  new FairTrackParH(initpos, initmom, ResPosErr, ResMomErr, charge);
599  Bool_t isProp = fPro->Propagate(fStart, fRes, mctrk->GetPdgCode());
600  if (!isProp) cout << "Propagate(): back propagation not possible! " << endl;
601  _mctrack = initmom;// inv_lmdrotation*_mctrack;
602 
603  }
604  det_id = mcpoint->GetDetectorID(); // store the detector id
605  sens_id = mcpoint->GetSensorID(); // store the sensor id
606  // calculate the plane and sensor on this plane
607  plane = sens_id/nsensors_per_plane;
608  sensor = sens_id-plane*nsensors_per_plane;
609 
610  // write those values out
611  px_init = momMC.X();
612  py_init = momMC.Y();
613  pz_init = momMC.Z();
614  ptheta_init = momMC.Theta();
615  pphi_init = momMC.Phi();
616  x_in = _mcpoint.X();
617  y_in = _mcpoint.Y();
618  z_in = _mcpoint.Z();
619  px_in = _mctrack.X();
620  py_in = _mctrack.Y();
621  pz_in = _mctrack.Z();
622  ptheta_in = _mctrack.Theta();
623  pphi_in = _mctrack.Phi();
624 
625  if (plane < 0 || plane >= nplanes || sensor < 0 || sensor >= nsensors_per_plane){
626  cout << " warning: calculated plane or sensor is wrong: plane " << plane << " sensor " << sensor << endl;
627  } else {
628  hist_xy [plane] ->Fill(x_in,y_in);
629  hists_xy[plane][sensor]->Fill(x_in,y_in);
630 
631  hist_theta_init [plane] ->Fill(ptheta_init);
632  hists_theta_init [plane][sensor]->Fill(ptheta_init);
633 
634  hist_theta_in [plane] ->Fill(ptheta_in);
635  hists_theta_in[plane][sensor]->Fill(ptheta_in);
636 
637  hist_theta_diff[plane] ->Fill(ptheta_in-ptheta_init);
638  hists_theta_diff[plane][sensor]->Fill(ptheta_in-ptheta_init);
639 
640  if (ptheta_init > 0){
641  hist_theta_diff_rel[plane] ->Fill((ptheta_in-ptheta_init)/ptheta_init);
642  hists_theta_diff_rel[plane][sensor]->Fill((ptheta_in-ptheta_init)/ptheta_init);
643  }
644 
645  tree_results->Fill();
646  }
647 
648  //cout << _mcpoint.X() << " " << _mcpoint.Y() << " "<< _mcpoint.Z() << endl;
649  hist_spatial_distr_acc->Fill(_mcpoint.X(), _mcpoint.Y());
650  hist_angular_distr_acc->Fill(momMC.Theta(), momMC.Phi());
651  nSdsHits++;
652  }
653  }
654  }
655  if (0 && nHits > 0)
656  cout << "found " << nHits << " hit(s) with " << nSdsHits << " sds hit(s)" << endl;
657  }
658 
659  hist_angular_distr_acc->Divide(hist_angular_distr_gen);
660 
661  TGaxis::SetMaxDigits(3);
662  TCanvas* canvas1 = new TCanvas("canvas1", "canvas1", 800, 800);
663  canvas1->Divide(2,2);
664  canvas1->cd(1);
665  hist_angular_distr_gen->Draw("COLZ");
666  hist_angular_distr_gen->SetXTitle("#Theta [rad]");
667  hist_angular_distr_gen->SetYTitle("#Phi [rad]");
668  gPad->Update();
669  canvas1->cd(2);
670  hist_angular_distr_acc->Draw("COLZ");
671  hist_angular_distr_acc->SetXTitle("#Theta [rad]");
672  hist_angular_distr_acc->SetYTitle("#Phi [rad]");
673  gPad->Update();
674  canvas1->cd(3);
675  hist_spatial_distr_gen->Draw("COLZ");
676  hist_spatial_distr_gen->SetXTitle("X [cm]");
677  hist_spatial_distr_gen->SetYTitle("Y [cm]");
678  hist_spatial_distr_gen->GetXaxis()->SetRangeUser(-10,10);
679  hist_spatial_distr_gen->GetYaxis()->SetRangeUser(-10,10);
680  gPad->Update();
681  canvas1->cd(4);
682  hist_spatial_distr_acc->Draw("COLZ");
683  hist_spatial_distr_acc->SetXTitle("X [cm]");
684  hist_spatial_distr_acc->SetYTitle("Y [cm]");
685  hist_spatial_distr_acc->GetXaxis()->SetRangeUser(-10,10);
686  hist_spatial_distr_acc->GetYaxis()->SetRangeUser(-10,10);
687  gPad->Update();
688 
689  // draw individual plane and sensor properties
690  //gStyle->SetPaperSize(25,25);
691  TCanvas canvas_properties_per_plane("canvas_properties_per_plane", "properties per plane", 800, 800);
692  canvas_properties_per_plane.cd();
693  // assuming to have 4 planes
694  if (nplanes != 4) cout << " warning: attempting to draw histograms for a number of planes not 4! " << endl;
695  canvas_properties_per_plane.Divide(2,2);
696 
697  TCanvas canvas_properties_per_sensor("canvas_properties_per_sensor", "properties per sensor", 800, 800);
698  canvas_properties_per_sensor.cd();
699  // assuming to have 8 sensors per plane
700  if (nsensors_per_plane != 8) cout << " warning: attempting to draw histograms for a number of sensors per plane not 8! " << endl;
701  canvas_properties_per_sensor.Divide(3,3);
702 
703  for (int iplane = 0; iplane < nplanes; iplane++){
704  canvas_properties_per_plane.cd(iplane+1);
705  hist_xy[iplane]->Draw("COLZ");
706  }
707  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps(");
708  for (int iplane = 0; iplane < nplanes; iplane++){
709  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
710  canvas_properties_per_sensor.cd(isensor+1);
711  hists_xy[iplane][isensor]->Draw("COLZ");
712  }
713  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
714  }
715 
716  for (int iplane = 0; iplane < nplanes; iplane++){
717  canvas_properties_per_plane.cd(iplane+1);
718  hist_theta_init[iplane]->Draw();
719  }
720  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps(");
721  for (int iplane = 0; iplane < nplanes; iplane++){
722  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
723  canvas_properties_per_sensor.cd(isensor+1);
724  hists_theta_init[iplane][isensor]->Draw();
725  canvas_properties_per_sensor.cd(9);
726  if (isensor > 0)
727  hists_theta_init[iplane][isensor]->Draw("same E");
728  else
729  hists_theta_init[iplane][isensor]->Draw("E");
730  hists_theta_init[iplane][isensor]->SetLineColor(kAzure-9+isensor);
731  }
732  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
733  }
734  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
735  canvas_properties_per_sensor.cd(isensor+1);
736  for (int iplane = 0; iplane < nplanes; iplane++){
737  if (iplane > 0)
738  hists_theta_init[iplane][isensor]->Draw("same");
739  else
740  hists_theta_init[iplane][isensor]->Draw();
741  hists_theta_init[iplane][isensor]->SetLineColor(kGray+iplane);
742  }
743  }
744  canvas_properties_per_sensor.cd(9);
745  gPad->Clear();
746  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
747 
748  for (int iplane = 0; iplane < nplanes; iplane++){
749  canvas_properties_per_plane.cd(iplane+1);
750  hist_theta_in[iplane]->Draw();
751  }
752  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps(");
753  for (int iplane = 0; iplane < nplanes; iplane++){
754  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
755  canvas_properties_per_sensor.cd(isensor+1);
756  hists_theta_in[iplane][isensor]->Draw();
757  canvas_properties_per_sensor.cd(9);
758  if (isensor > 0)
759  hists_theta_in[iplane][isensor]->Draw("same E");
760  else
761  hists_theta_in[iplane][isensor]->Draw("E");
762  hists_theta_in[iplane][isensor]->SetLineColor(kAzure-9+isensor);
763  }
764  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
765  }
766  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
767  canvas_properties_per_sensor.cd(isensor+1);
768  for (int iplane = 0; iplane < nplanes; iplane++){
769  if (iplane > 0)
770  hists_theta_in[iplane][isensor]->Draw("same");
771  else
772  hists_theta_in[iplane][isensor]->Draw();
773  hists_theta_in[iplane][isensor]->SetLineColor(kGray+iplane);
774  }
775  }
776  canvas_properties_per_sensor.cd(9);
777  gPad->Clear();
778  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
779 
780  for (int iplane = 0; iplane < nplanes; iplane++){
781  canvas_properties_per_plane.cd(iplane+1);
782  hist_theta_diff[iplane]->Draw();
783  }
784  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps(");
785  for (int iplane = 0; iplane < nplanes; iplane++){
786  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
787  canvas_properties_per_sensor.cd(isensor+1);
788  hists_theta_diff[iplane][isensor]->Draw();
789  canvas_properties_per_sensor.cd(9);
790  if (isensor > 0)
791  hists_theta_diff[iplane][isensor]->Draw("same E");
792  else
793  hists_theta_diff[iplane][isensor]->Draw("E");
794  hists_theta_diff[iplane][isensor]->SetLineColor(kAzure-9+isensor);
795  }
796  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
797  }
798  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
799  canvas_properties_per_sensor.cd(isensor+1);
800  for (int iplane = 0; iplane < nplanes; iplane++){
801  if (iplane > 0)
802  hists_theta_diff[iplane][isensor]->Draw("same");
803  else
804  hists_theta_diff[iplane][isensor]->Draw();
805  hists_theta_diff[iplane][isensor]->SetLineColor(kGray+iplane);
806  }
807  }
808  canvas_properties_per_sensor.cd(9);
809  gPad->Clear();
810  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
811 
812  for (int iplane = 0; iplane < nplanes; iplane++){
813  canvas_properties_per_plane.cd(iplane+1);
814  hist_theta_diff_rel[iplane]->Draw();
815  }
816  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps(");
817  for (int iplane = 0; iplane < nplanes; iplane++){
818  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
819  canvas_properties_per_sensor.cd(isensor+1);
820  hists_theta_diff_rel[iplane][isensor]->Draw();
821  canvas_properties_per_sensor.cd(9);
822  if (isensor > 0)
823  hists_theta_diff_rel[iplane][isensor]->Draw("same E");
824  else
825  hists_theta_diff_rel[iplane][isensor]->Draw("E");
826  hists_theta_diff_rel[iplane][isensor]->SetLineColor(kAzure-9+isensor);
827  }
828  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
829  }
830  for (int isensor = 0; isensor < nsensors_per_plane; isensor++){
831  canvas_properties_per_sensor.cd(isensor+1);
832  for (int iplane = 0; iplane < nplanes; iplane++){
833  if (iplane > 0)
834  hists_theta_diff_rel[iplane][isensor]->Draw("same");
835  else
836  hists_theta_diff_rel[iplane][isensor]->Draw();
837  hists_theta_diff_rel[iplane][isensor]->SetLineColor(kGray+iplane);
838  }
839  }
840  canvas_properties_per_sensor.cd(9);
841  gPad->Clear();
842  canvas_properties_per_sensor.Print("Resolution_acceptance_results.ps(");
843 
844  canvas_properties_per_plane.Clear();
845  canvas_properties_per_plane.Print("Resolution_acceptance_results.ps)");
846  // convert to pdf
847  cout << "converting to pdf " << system("ps2pdf Resolution_acceptance_results.ps")<< endl;
848  cout << "deleting ps file " << system("rm Resolution_acceptance_results.ps")<< endl;
849 
850  output_histfile->Write();
851 
852  hist_angular_distr_gen->SetDirectory(0);
853  hist_angular_distr_acc->SetDirectory(0);
854  hist_spatial_distr_gen->SetDirectory(0);
855  hist_spatial_distr_acc->SetDirectory(0);
856  output_histfile->Close();
857  }
858 
859  void DrawProgressBar(int len, double percent) {
860  if ((int)(last_percent*100) == (int)(percent*100)) return;
861  //cout << " drawing " << endl;
862  cout << "\x1B[2K"; // Erase the entire current line.
863  cout << "\x1B[0E"; // Move to the beginning of the current line.
864  string progress;
865  for (int i = 0; i < len; ++i) {
866  if (i < static_cast<int>(len * percent)) {
867  progress += "=";
868  } else {
869  progress += " ";
870  }
871  }
872  cout << "[" << progress << "] " << (static_cast<int>(100 * percent)) << "%";
873  flush(cout); // Required.
874  last_percent = percent;
875  }*/
876 
int main(int argc, char **argv)
TString RecoFile
TString outFile
Definition: hit_dirc.C:17
Int_t startEvent
TString storePath
FairGeane * Geane
FairRunAna * fRun
Definition: hit_dirc.C:58
TString DigiFile
Double_t
TString parFile
Definition: hit_dirc.C:14
TStopwatch timer
Definition: hit_dirc.C:51
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
void SetHistFilename(TString filename)
FairParRootFileIo * parInput1
Definition: hit_dirc.C:67
Double_t ctime
Definition: hit_dirc.C:114
TString MCFile
Double_t rtime
Definition: hit_dirc.C:113