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