]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtest.C
Merging added and test macro
[u/mrichter/AliRoot.git] / TOF / AliTOFtest.C
1 Int_t AliTOFtest() 
2 {
3   //
4   // Test macro for the TOF code
5   // report bug to Fabrizio.Pierella@cern.ch
6   // Use case:
7   // start aliroot
8   // root [0] .L AliTOFtest.C
9   // root [1] AliTOFtest()
10
11   Int_t rc = 0;
12
13   // Initialize the test setup 
14
15   //gAlice->Init("$(ALICE_ROOT)/TOF/AliTOFconfig.C");
16   gAlice->Init("AliTOFconfig.C");
17
18   // Run one central Hijing event and create the hits (time required: 
19   // some minuts)
20
21   gAlice->SetDebug(2);
22   gAlice->Run(1);
23   
24   if (gAlice) delete gAlice;
25   TFile *file = (TFile *) gROOT->GetListOfFiles()->FindObject("TOF_test.root");
26   gAlice = (AliRun *) file->Get("gAlice");
27
28   // Analyze the TOF hits
29   if (rc = AliTOFanalyzeHits()) return rc;
30
31   return rc;
32
33 }
34
35 //_____________________________________________________________________________
36 Int_t AliTOFanalyzeHits()
37 {
38   //
39   // Analyzes the hits and fills QA-histograms 
40   //
41
42   Int_t rc = 0;
43
44   if (!gAlice) {
45     cout << "<AliTOFanalyzeHits> No AliRun object found" << endl;
46     rc = 1;
47     return rc;
48   }
49   gAlice->GetEvent(0);
50
51   // Get the pointer to the TOF detector 
52   AliTOF *tof = (AliTOF *) gAlice->GetDetector("TOF");
53   if (!tof) {
54     cout << "<AliTOFanalyzeHits> No TOF detector found" << endl;
55     rc = 2;
56     return rc;
57   }
58
59   // Define the histograms
60   // x,y,z, rho, tof, padx, padz, sector, plate, strip, (x vs y)
61
62   // hit-map in a plane
63   TH2F *h2hitMap = new TH2F("h2hitMap","Hit Map (projection on the plane)",2500,-12500.,12500.,800,-400.,400.);
64   // time of flight distribution for primaries and secondaries
65   TH1F *htofp    = new TH1F("htofp","Time of Flight (primaries)",800,0.,80.);
66   TH1F *htofs    = new TH1F("htofs","Time of Flight (secondaries)",800,0.,80.);
67   // momentum when striking the TOF for primaries and secondaries
68   TH1F *htofmomp = new TH1F("htofmomp","Momentum at TOF (primaries)",100,0.,10.);
69   TH1F *htofmoms = new TH1F("htofmoms","Momentum at TOF (secondaries)",100,0.,10.);
70   // TOF hit volumes
71   TH1F *hsector  = new TH1F("hsector","Sector",20,0.,20.);
72   TH1F *hplate   = new TH1F("hplate","Plate ", 6,0., 6.);
73   TH1F *hstrip   = new TH1F("hstrip","Strip ",25,0.,25.);
74   TH1F *hpadz    = new TH1F("hpadz","Pad along z ",3,0.,3.);
75   TH1F *hpadx    = new TH1F("hpadx","Pad along x",50,0.,50.);
76   // track length when striking the TOF (used by AliTOFT0)
77   TH1F *htrackLenp= new TH1F("htrackLenp","Track Length on TOF for Primaries",800,0.,800.);
78
79   // Get the pointer hit tree
80   TTree *hitTree = gAlice->TreeH();  
81   if (!hitTree) {
82     cout << "<AliTOFanalyzeHits> No hit tree found" << endl;
83     rc = 4;
84     return rc;
85   }
86
87   Int_t countHits = 0;
88   Int_t nBytes    = 0;
89
90   // Get the number of entries in the hit tree
91   // (Number of primary particles creating a hit somewhere)
92   Int_t nTrack = (Int_t) hitTree->GetEntries();
93   cout << "<AliTOFanalyzeHits> Found " << nTrack 
94        << " primary particles with hits" << endl;
95
96   Int_t nPrimaryOnTof = 0;
97   Int_t nSecondaryOnTof = 0;
98   Int_t nelectron  = 0;
99   Int_t npion      = 0;
100   Int_t nkaon      = 0;
101   Int_t nproton    = 0;    
102   Int_t nmuon      = 0;
103   
104   // Loop through all entries in the tree
105   for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
106
107     gAlice->ResetHits();
108     nBytes += hitTree->GetEvent(iTrack);
109
110
111     // Loop through the TOF hits  
112     Int_t iHit = 0;
113     AliTOFhitT0 *hit = (AliTOFhitT0 *) tof->FirstHit(-1);
114     while (hit) {
115
116       countHits++;
117       iHit++;
118
119       Float_t x     = hit->X();
120       Float_t y     = hit->Y();
121       Float_t z     = hit->Z();
122       Float_t circleLen=TMath::Sqrt(x*x+y*y)*TMath::ATan2(y,x);
123       h2hitMap->Fill(circleLen,z);
124
125       Float_t flightTime = hit->GetTof(); // [s]
126       flightTime*=1.e+09; // convert in [ns]
127       Float_t angle = hit->GetIncA();
128       Float_t tofmom= hit->GetMom(); // [GeV/c]
129       Float_t trackLen= hit->GetLen(); // [cm]
130
131       // TOF hit volumes
132       Int_t sector    = hit->GetSector(); // range [1-18]
133       Int_t plate     = hit->GetPlate();  // range [1- 5]
134       Int_t strip     = hit->GetStrip();  // range [1-20]
135       Int_t padz      = hit->GetPadz();   // range [1- 2]
136       Int_t padx      = hit->GetPadx();   // range [1-48]
137       // it is QA, then I perform QA!
138       Bool_t isHitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
139
140       if (isHitBad) {
141         cout << "<AliTOFanalyzeHits>  strange hit found" << endl;
142         rc = 3;
143         return rc;
144       }
145       // filling hit volume histos
146       hsector->Fill(sector);
147       hplate->Fill(plate);
148       hstrip->Fill(strip);
149       hpadx->Fill(padx);
150       hpadz->Fill(padz);
151
152
153       Int_t track     = hit->Track();
154       TParticle *part = gAlice->Particle(track);
155
156       // getting MC info for the current track
157       if (part->GetFirstMother()<0){
158         Int_t icod = TMath::Abs(part->GetPdgCode());
159         switch (icod) {
160         case 211:
161           npion++;
162           break ;
163         case 321:
164           nkaon++;
165           break ;
166         case 2212:
167           nproton++;
168           break ;
169         case 11:
170           nelectron++;
171           break ;
172         case 13:
173           nmuon++;
174           break ;
175         }
176         htofp->Fill(flightTime);
177         htofmomp->Fill(tofmom);
178         htrackLenp->Fill(trackLen);
179       } else {
180         htofs->Fill(flightTime);
181         htofmoms->Fill(tofmom);
182       }
183
184       // go to next hit
185       hit = (AliTOFhitT0 *) tof->NextHit();         
186
187     }
188
189   }
190
191   cout << "<AliTOFanalyzeHits> Found " << countHits << " hits in total" << endl;
192   cout << npion     << " primary pions reached the TOF detector"     << endl;
193   cout << nkaon     << " primary kaons reached the TOF detector"     << endl;
194   cout << nproton   << " primary protons reached the TOF detector"   << endl;
195   cout << nelectron << " primary electrons reached the TOF detector" << endl;
196   cout << nmuon     << " primary muons reached the TOF detector"     << endl;
197
198   
199   TCanvas *cHits = new TCanvas("cHits","AliTOFanalyzeHits hit volumes",50,50,900,900);
200   cHits->Divide(3,2);
201   cHits->cd(1);
202   hsector->Draw();
203   cHits->cd(2);
204   hplate->Draw();
205   cHits->cd(3);
206   hstrip->Draw();
207   cHits->cd(4);
208   hpadz->Draw();
209   cHits->cd(5);
210   hpadx->Draw();
211
212   TCanvas *chitmap = new TCanvas("chitmap","AliTOFanalyzeHits Hit Map",50,50,600,600);
213   chitmap->cd();
214   h2hitMap->Draw();
215
216   TCanvas *ctrackLen = new TCanvas("ctrackLen","AliTOFanalyzeHits Track Length for primaries on TOF",50,50,400,400);
217   ctrackLen->cd();
218   htrackLenp->Draw();
219
220   TCanvas *ctofmom = new TCanvas("ctofmom","AliTOFanalyzeHits flight times",50,50,700,700);
221   ctofmom->Divide(2,2);
222   ctofmom->cd(1);
223   gPad->SetLogy();
224   htofp->Draw();
225   ctofmom->cd(2);
226   gPad->SetLogy();
227   htofs->Draw();
228   ctofmom->cd(3);
229   gPad->SetLogy();
230   htofmomp->Draw();
231   ctofmom->cd(4);
232   gPad->SetLogy();
233   htofmoms->Draw();
234   
235
236   // save histos into file TOF_hitsQA.root
237   TFile *fout = new TFile("TOF_hitsQA.root","RECREATE");
238   h2hitMap->Write();
239   htofp->Write();
240   htofs->Write();
241   htofmomp->Write();
242   htofmoms->Write();
243   hsector->Write();
244   hplate->Write();
245   hstrip->Write();
246   hpadz->Write();
247   hpadx->Write();
248   htrackLenp->Write();
249   fout->Close(); 
250
251   return rc;
252
253 }