]>
Commit | Line | Data |
---|---|---|
7f6ddf58 | 1 | /**************************************************************************** |
2 | * Very important, delicate and rather obscure macro. * | |
3 | * * | |
4 | * Creates list of "trackable" tracks, * | |
7f6ddf58 | 5 | * calculates efficiency, resolutions etc. * |
bf018fe4 | 6 | * There is a possibility to run this macro over several events. * |
7f6ddf58 | 7 | * * |
024a7fe9 | 8 | * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch * |
7f6ddf58 | 9 | * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de * |
10 | ****************************************************************************/ | |
11 | ||
c630aafd | 12 | #if !defined(__CINT__) || defined(__MAKECINT__) |
5102bab6 | 13 | #include <TMath.h> |
14 | #include <TError.h> | |
15 | #include <Riostream.h> | |
16 | #include <TH1F.h> | |
17 | #include <TH2F.h> | |
18 | #include <TTree.h> | |
19 | #include <TParticle.h> | |
5102bab6 | 20 | #include <TCanvas.h> |
21 | #include <TLine.h> | |
22 | #include <TText.h> | |
23 | #include <TBenchmark.h> | |
24 | #include <TStyle.h> | |
8b462fd8 | 25 | #include <TFile.h> |
bf018fe4 | 26 | #include <TROOT.h> |
5102bab6 | 27 | |
28 | #include "AliStack.h" | |
29 | #include "AliHeader.h" | |
30 | #include "AliTrackReference.h" | |
88cb7938 | 31 | #include "AliRunLoader.h" |
c630aafd | 32 | #include "AliRun.h" |
4514157e | 33 | #include "AliESDEvent.h" |
5102bab6 | 34 | |
35 | #include "AliSimDigits.h" | |
36 | #include "AliTPC.h" | |
37 | #include "AliTPCParamSR.h" | |
38 | #include "AliTPCClustersArray.h" | |
39 | #include "AliTPCClustersRow.h" | |
bf018fe4 | 40 | #include "AliTPCcluster.h" |
5102bab6 | 41 | #include "AliTPCLoader.h" |
be9c5115 | 42 | #endif |
43 | ||
4fa1e629 | 44 | Int_t GoodTracksTPC(const Char_t *dir="."); |
88cb7938 | 45 | |
c630aafd | 46 | extern AliRun *gAlice; |
5102bab6 | 47 | extern TBenchmark *gBenchmark; |
bf018fe4 | 48 | extern TROOT *gROOT; |
fbb58996 | 49 | |
4fa1e629 | 50 | static Int_t allgood=0; |
51 | static Int_t allselected=0; | |
52 | static Int_t allfound=0; | |
53 | ||
54 | Int_t AliTPCComparison | |
55 | (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") { | |
5102bab6 | 56 | gBenchmark->Start("AliTPCComparison"); |
57 | ||
58 | ::Info("AliTPCComparison.C","Doing comparison..."); | |
c630aafd | 59 | |
88cb7938 | 60 | |
bf018fe4 | 61 | |
62 | TH1F *hp=(TH1F*)gROOT->FindObject("hp"); | |
63 | if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.); | |
64 | hp->SetFillColor(4); | |
65 | ||
66 | TH1F *hl=(TH1F*)gROOT->FindObject("hl"); | |
67 | if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20); | |
68 | hl->SetFillColor(4); | |
69 | ||
70 | TH1F *hpt=(TH1F*)gROOT->FindObject("hpt"); | |
71 | if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); | |
72 | hpt->SetFillColor(2); | |
73 | ||
74 | TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt"); | |
75 | if (!hmpt) | |
76 | hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); | |
cc80f89e | 77 | hmpt->SetFillColor(6); |
8c555625 | 78 | |
bf018fe4 | 79 | |
80 | ||
81 | TH1F *hgood=(TH1F*)gROOT->FindObject("hgood"); | |
82 | if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1); | |
83 | ||
84 | TH1F *hfound=(TH1F*)gROOT->FindObject("hfound"); | |
85 | if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1); | |
86 | ||
87 | TH1F *hfake=(TH1F*)gROOT->FindObject("hfake"); | |
88 | if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1); | |
89 | ||
90 | TH1F *hg=(TH1F*)gROOT->FindObject("hg"); | |
91 | if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1); | |
8c555625 | 92 | hg->SetLineColor(4); hg->SetLineWidth(2); |
bf018fe4 | 93 | |
94 | TH1F *hf=(TH1F*)gROOT->FindObject("hf"); | |
95 | if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1); | |
8c555625 | 96 | hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2); |
97 | ||
bf018fe4 | 98 | TH1F *he=(TH1F*)gROOT->FindObject("he"); |
99 | if (!he) | |
100 | he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.); | |
101 | ||
102 | TH2F *hep=(TH2F*)gROOT->FindObject("hep"); | |
103 | if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.); | |
be9c5115 | 104 | hep->SetMarkerStyle(8); |
105 | hep->SetMarkerSize(0.4); | |
8c555625 | 106 | |
4fa1e629 | 107 | |
108 | Char_t fname[100]; | |
109 | sprintf(fname,"%s/GoodTracksTPC.root",dir); | |
110 | ||
111 | TFile *refFile=TFile::Open(fname,"old"); | |
112 | if (!refFile || !refFile->IsOpen()) { | |
113 | ::Info("AliTPCComparison.C","Marking good tracks (will take a while)..."); | |
114 | if (GoodTracksTPC(dir)) { | |
115 | ::Error("AliTPCComparison.C","Can't generate the reference file !"); | |
116 | return 1; | |
117 | } | |
118 | } | |
119 | refFile=TFile::Open(fname,"old"); | |
120 | if (!refFile || !refFile->IsOpen()) { | |
121 | ::Error("AliTPCComparison.C","Can't open the reference file !"); | |
122 | return 2; | |
123 | } | |
124 | ||
125 | TTree *tpcTree=(TTree*)refFile->Get("tpcTree"); | |
126 | if (!tpcTree) { | |
127 | ::Error("AliTPCComparison.C","Can't get the reference tree !"); | |
128 | return 3; | |
129 | } | |
130 | TBranch *branch=tpcTree->GetBranch("TPC"); | |
131 | if (!branch) { | |
132 | ::Error("AliTPCComparison.C","Can't get the TPC branch !"); | |
133 | return 4; | |
134 | } | |
135 | TClonesArray dummy("AliTrackReference",1000), *refs=&dummy; | |
136 | branch->SetAddress(&refs); | |
137 | ||
138 | ||
139 | sprintf(fname,"%s/AliESDs.root",dir); | |
140 | TFile *ef=TFile::Open(fname); | |
141 | if ((!ef)||(!ef->IsOpen())) { | |
142 | sprintf(fname,"%s/AliESDtpc.root",dir); | |
143 | ef=TFile::Open(fname); | |
144 | if ((!ef)||(!ef->IsOpen())) { | |
145 | ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !"); | |
146 | return 5; | |
147 | } | |
148 | } | |
4514157e | 149 | AliESDEvent* event = new AliESDEvent(); |
8b462fd8 | 150 | TTree* esdTree = (TTree*) ef->Get("esdTree"); |
151 | if (!esdTree) { | |
152 | ::Error("AliTPCComparison.C", "no ESD tree found"); | |
153 | return 6; | |
154 | } | |
4514157e | 155 | event->ReadFromTree(esdTree); |
4fa1e629 | 156 | |
157 | ||
158 | //******* Loop over events ********* | |
159 | Int_t e=0; | |
8b462fd8 | 160 | while (esdTree->GetEvent(e)) { |
4fa1e629 | 161 | cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n"; |
162 | ||
4fa1e629 | 163 | Int_t nentr=event->GetNumberOfTracks(); |
164 | allfound+=nentr; | |
165 | ||
166 | if (tpcTree->GetEvent(e++)==0) { | |
167 | cerr<<"No reconstructable tracks !\n"; | |
168 | continue; | |
169 | } | |
170 | ||
171 | Int_t ngood=refs->GetEntriesFast(); | |
172 | allgood+=ngood; | |
173 | ||
174 | const Int_t MAX=15000; | |
175 | //MI change | |
176 | Int_t track_notfound[MAX], itrack_notfound=0; | |
177 | Int_t track_fake[MAX], itrack_fake=0; | |
178 | Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0; | |
179 | ||
180 | Int_t i; | |
181 | for (Int_t k=0; k<ngood; k++) { | |
182 | AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k); | |
183 | Int_t lab=ref->Label(), tlab=-1; | |
184 | Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py()); | |
185 | ||
186 | if (ptg<1e-33) continue; // for those not crossing 0 pad row | |
187 | ||
188 | if (ptg<ptcutl) continue; | |
189 | if (ptg>ptcuth) continue; | |
190 | ||
191 | allselected++; | |
192 | ||
193 | hgood->Fill(ptg); | |
194 | ||
195 | AliESDtrack *track=0; | |
196 | for (i=0; i<nentr; i++) { | |
bf018fe4 | 197 | track=event->GetTrack(i); |
162085b7 | 198 | tlab=track->GetTPCLabel(); |
cc80f89e | 199 | if (lab==TMath::Abs(tlab)) break; |
4fa1e629 | 200 | } |
201 | if (i==nentr) { | |
202 | track_notfound[itrack_notfound++]=lab; | |
203 | continue; | |
204 | } | |
205 | ||
206 | //MI change - addition | |
207 | Int_t micount=0; | |
208 | Int_t mi; | |
209 | AliESDtrack * mitrack; | |
210 | for (mi=0; mi<nentr; mi++) { | |
211 | mitrack=event->GetTrack(mi); | |
212 | if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++; | |
213 | } | |
214 | if (micount>1) { | |
215 | track_multifound[itrack_multifound]=lab; | |
216 | track_multifound_n[itrack_multifound]=micount; | |
217 | itrack_multifound++; | |
218 | } | |
219 | if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue; | |
220 | if (lab==tlab) hfound->Fill(ptg); | |
221 | else { | |
222 | track_fake[itrack_fake++]=lab; | |
223 | hfake->Fill(ptg); | |
224 | } | |
cf98c13f | 225 | |
4fa1e629 | 226 | Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz); |
227 | Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]); | |
228 | if (phi<-TMath::Pi()) phi+=2*TMath::Pi(); | |
229 | if (phi>=TMath::Pi()) phi-=2*TMath::Pi(); | |
230 | Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]); | |
231 | Float_t lam=TMath::ATan2(pxpypz[2],pt); | |
232 | Float_t pt_1=1/pt; | |
233 | ||
234 | Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG ! | |
235 | ||
236 | if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons | |
237 | hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.); | |
238 | } else { | |
239 | Float_t phig=TMath::ATan2(ref->Py(),ref->Px()); | |
240 | hp->Fill((phi - phig)*1000.); | |
241 | ||
242 | Float_t lamg=TMath::ATan2(ref->Pz(),ptg); | |
243 | hl->Fill((lam - lamg)*1000.); | |
244 | ||
245 | hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.); | |
246 | } | |
247 | ||
248 | Float_t mom=pt/TMath::Cos(lam); | |
249 | Float_t dedx=track->GetTPCsignal(); | |
250 | hep->Fill(mom,dedx,1.); | |
251 | if (TMath::Abs(pdg)==211) //pions | |
252 | if (mom>0.4 && mom<0.5) { | |
253 | he->Fill(dedx,1.); | |
254 | } | |
cf98c13f | 255 | } |
4fa1e629 | 256 | |
257 | cout<<"\nList of Not found tracks :\n"; | |
258 | for ( i = 0; i< itrack_notfound; i++){ | |
259 | cout<<track_notfound[i]<<"\t"; | |
260 | if ((i%5)==4) cout<<"\n"; | |
cf98c13f | 261 | } |
4fa1e629 | 262 | cout<<"\nList of fake tracks :\n"; |
263 | for ( i = 0; i< itrack_fake; i++){ | |
264 | cout<<track_fake[i]<<"\t"; | |
265 | if ((i%5)==4) cout<<"\n"; | |
cf98c13f | 266 | } |
4fa1e629 | 267 | cout<<"\nList of multiple found tracks :\n"; |
268 | for ( i=0; i<itrack_multifound; i++) { | |
269 | cout<<"id. "<<track_multifound[i] | |
270 | <<" found - "<<track_multifound_n[i]<<"times\n"; | |
cc80f89e | 271 | } |
73042f01 | 272 | |
4fa1e629 | 273 | cout<<"Number of found tracks ="<<nentr<<endl; |
274 | cout<<"Number of \"good\" tracks ="<<ngood<<endl; | |
73042f01 | 275 | |
4fa1e629 | 276 | refs->Clear(); |
4fa1e629 | 277 | }// ***** End of the loop over events |
cc80f89e | 278 | |
8b462fd8 | 279 | delete event; |
4514157e | 280 | delete esdTree; |
4fa1e629 | 281 | ef->Close(); |
8c555625 | 282 | |
4fa1e629 | 283 | delete tpcTree; |
284 | refFile->Close(); | |
7f6ddf58 | 285 | |
4fa1e629 | 286 | Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries(); |
287 | if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n"; | |
288 | cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl; | |
289 | cout<<"Total number of found tracks ="<<allfound<<endl; | |
290 | cout<<"Total number of \"good\" tracks ="<<allgood<<endl; | |
291 | cout<<endl; | |
7f6ddf58 | 292 | |
8c555625 | 293 | gStyle->SetOptStat(111110); |
294 | gStyle->SetOptFit(1); | |
295 | ||
296 | TCanvas *c1=new TCanvas("c1","",0,0,700,850); | |
297 | ||
2d343a94 | 298 | Int_t minc=33; |
299 | ||
8c555625 | 300 | TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw(); |
301 | p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); | |
2d343a94 | 302 | hp->SetFillColor(4); hp->SetXTitle("(mrad)"); |
303 | if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd(); | |
8c555625 | 304 | |
305 | TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); | |
306 | p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10); | |
2d343a94 | 307 | hl->SetXTitle("(mrad)"); |
308 | if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd(); | |
8c555625 | 309 | |
310 | TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw(); | |
311 | p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); | |
2d343a94 | 312 | hpt->SetXTitle("(%)"); |
313 | if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd(); | |
8c555625 | 314 | |
315 | TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw(); | |
316 | p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10); | |
2d343a94 | 317 | hmpt->SetXTitle("(%)"); |
318 | if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd(); | |
8c555625 | 319 | |
320 | TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); | |
3c0f9266 | 321 | p5->SetFillColor(41); p5->SetFrameFillColor(10); |
8c555625 | 322 | hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2(); |
323 | hg->Divide(hfound,hgood,1,1.,"b"); | |
324 | hf->Divide(hfake,hgood,1,1.,"b"); | |
325 | hg->SetMaximum(1.4); | |
326 | hg->SetYTitle("Tracking efficiency"); | |
327 | hg->SetXTitle("Pt (GeV/c)"); | |
328 | hg->Draw(); | |
329 | ||
7f6ddf58 | 330 | TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4); |
8c555625 | 331 | line1->Draw("same"); |
7f6ddf58 | 332 | TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4); |
8c555625 | 333 | line2->Draw("same"); |
334 | ||
335 | hf->SetFillColor(1); | |
336 | hf->SetFillStyle(3013); | |
337 | hf->SetLineColor(2); | |
338 | hf->SetLineWidth(2); | |
339 | hf->Draw("histsame"); | |
340 | TText *text = new TText(0.461176,0.248448,"Fake tracks"); | |
341 | text->SetTextSize(0.05); | |
342 | text->Draw(); | |
343 | text = new TText(0.453919,1.11408,"Good tracks"); | |
344 | text->SetTextSize(0.05); | |
345 | text->Draw(); | |
3c0f9266 | 346 | |
347 | ||
348 | ||
349 | TCanvas *c2=new TCanvas("c2","",320,32,530,590); | |
350 | ||
351 | TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw(); | |
352 | p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); | |
353 | he->SetFillColor(2); he->SetFillStyle(3005); | |
354 | he->SetXTitle("Arbitrary Units"); | |
2d343a94 | 355 | if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd(); |
3c0f9266 | 356 | |
357 | TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); | |
358 | p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10); | |
359 | hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); | |
360 | hep->Draw(); c1->cd(); | |
361 | ||
2d343a94 | 362 | TFile fc("AliTPCComparison.root","RECREATE"); |
363 | c1->Write(); | |
364 | c2->Write(); | |
365 | fc.Close(); | |
366 | ||
4fa1e629 | 367 | gBenchmark->Stop("AliTPCComparison"); |
368 | gBenchmark->Show("AliTPCComparison"); | |
369 | ||
73042f01 | 370 | return 0; |
cc80f89e | 371 | } |
3c0f9266 | 372 | |
cc80f89e | 373 | |
4fa1e629 | 374 | Int_t GoodTracksTPC(const Char_t *dir) { |
375 | if (gAlice) { | |
376 | delete gAlice->GetRunLoader(); | |
377 | delete gAlice;//if everything was OK here it is already NULL | |
378 | gAlice = 0x0; | |
379 | } | |
380 | ||
381 | Char_t fname[100]; | |
382 | sprintf(fname,"%s/galice.root",dir); | |
cc80f89e | 383 | |
4fa1e629 | 384 | AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON"); |
385 | if (!rl) { | |
386 | ::Error("GoodTracksTPC","Can't start session !"); | |
387 | return 1; | |
c630aafd | 388 | } |
4fa1e629 | 389 | |
390 | rl->LoadgAlice(); | |
391 | rl->LoadHeader(); | |
392 | rl->LoadKinematics(); | |
393 | rl->LoadTrackRefs(); | |
394 | ||
c630aafd | 395 | AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader"); |
396 | if (tpcl == 0x0) { | |
4fa1e629 | 397 | ::Error("GoodTracksTPC","Can not find TPCLoader !"); |
398 | delete rl; | |
399 | return 2; | |
c630aafd | 400 | } |
4fa1e629 | 401 | |
88cb7938 | 402 | AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC"); |
73042f01 | 403 | Int_t ver = TPC->IsVersion(); |
4fa1e629 | 404 | cout<<"TPC version "<<ver<<" has been found !\n"; |
405 | if (ver==1) tpcl->LoadRecPoints(); | |
406 | else if (ver==2) tpcl->LoadDigits(); | |
407 | else { | |
408 | ::Error("GoodTracksTPC","Wrong TPC version !"); | |
409 | delete rl; | |
410 | return 3; | |
411 | } | |
73042f01 | 412 | |
88cb7938 | 413 | rl->CdGAFile(); |
414 | AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60"); | |
5102bab6 | 415 | if (!digp) { |
4fa1e629 | 416 | ::Error("AliTPCHits2Digits.C","TPC parameters have not been found !"); |
417 | delete rl; | |
418 | return 4; | |
5102bab6 | 419 | } |
73042f01 | 420 | TPC->SetParam(digp); |
421 | ||
422 | Int_t nrow_up=digp->GetNRowUp(); | |
423 | Int_t nrows=digp->GetNRowLow()+nrow_up; | |
424 | Int_t zero=digp->GetZeroSup(); | |
7f6ddf58 | 425 | Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap); |
cc80f89e | 426 | Int_t good_number=Int_t(0.4*nrows); |
427 | ||
4fa1e629 | 428 | Int_t nev=rl->GetNumberOfEvents(); |
429 | ::Info("GoodTracksTPC","Number of events : %d\n",nev); | |
430 | ||
431 | ||
432 | sprintf(fname,"%s/GoodTracksTPC.root",dir); | |
433 | TFile *file=TFile::Open(fname,"recreate"); | |
434 | ||
435 | TClonesArray dummy("AliTrackReference",1000), *refs=&dummy; | |
436 | TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks"); | |
437 | tpcTree.Branch("TPC",&refs); | |
438 | ||
439 | //******** Loop over generated events | |
440 | for (Int_t e=0; e<nev; e++) { | |
441 | Int_t nt=0; | |
442 | refs->Clear(); | |
443 | ||
444 | Int_t i; | |
445 | ||
446 | rl->GetEvent(e); file->cd(); | |
447 | ||
448 | Int_t np = rl->GetHeader()->GetNtrack(); | |
449 | cout<<"Event "<<e<<" Number of particles: "<<np<<endl; | |
450 | ||
451 | //******** Fill the "good" masks | |
452 | Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0; | |
453 | switch (ver) { | |
454 | case 1: | |
455 | { | |
456 | AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca; | |
457 | ca.Setup(digp); | |
458 | ca.SetClusterType("AliTPCcluster"); | |
459 | ca.ConnectTree(tpcl->TreeR()); | |
460 | Int_t nrows=Int_t(ca.GetTree()->GetEntries()); | |
461 | for (Int_t n=0; n<nrows; n++) { | |
7f6ddf58 | 462 | AliSegmentID *s=ca.LoadEntry(n); |
463 | Int_t sec,row; | |
464 | digp->AdjustSectorRow(s->GetID(),sec,row); | |
465 | AliTPCClustersRow &clrow = *ca.GetRow(sec,row); | |
466 | Int_t ncl=clrow.GetArray()->GetEntriesFast(); | |
467 | while (ncl--) { | |
468 | AliTPCcluster *c=(AliTPCcluster*)clrow[ncl]; | |
469 | Int_t lab=c->GetLabel(0); | |
470 | if (lab<0) continue; //noise cluster | |
471 | lab=TMath::Abs(lab); | |
472 | ||
473 | if (sec>=digp->GetNInnerSector()) | |
474 | if (row==nrow_up-1) good[lab]|=0x4000; | |
475 | if (sec>=digp->GetNInnerSector()) | |
476 | if (row==nrow_up-1-gap) good[lab]|=0x1000; | |
477 | ||
478 | if (sec>=digp->GetNInnerSector()) | |
479 | if (row==nrow_up-1-shift) good[lab]|=0x2000; | |
480 | if (sec>=digp->GetNInnerSector()) | |
481 | if (row==nrow_up-1-gap-shift) good[lab]|=0x800; | |
482 | ||
483 | good[lab]++; | |
484 | } | |
485 | ca.ClearRow(sec,row); | |
4fa1e629 | 486 | } |
487 | delete pca; | |
488 | } | |
489 | break; | |
490 | case 2: | |
491 | { | |
492 | TTree *TD=tpcl->TreeD(); | |
88cb7938 | 493 | |
4fa1e629 | 494 | AliSimDigits da, *digits=&da; |
495 | TD->GetBranch("Segment")->SetAddress(&digits); | |
7f6ddf58 | 496 | |
4fa1e629 | 497 | Int_t *count = new Int_t[np]; |
498 | Int_t i; | |
499 | for (i=0; i<np; i++) count[i]=0; | |
7f6ddf58 | 500 | |
4fa1e629 | 501 | Int_t sectors_by_rows=(Int_t)TD->GetEntries(); |
502 | for (i=0; i<sectors_by_rows; i++) { | |
7f6ddf58 | 503 | if (!TD->GetEvent(i)) continue; |
504 | Int_t sec,row; | |
505 | digp->AdjustSectorRow(digits->GetID(),sec,row); | |
4fa1e629 | 506 | //cerr<<sec<<' '<<row<<'\r'; |
7f6ddf58 | 507 | digits->First(); |
508 | do { //Many thanks to J.Chudoba who noticed this | |
509 | Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn(); | |
510 | Short_t dig = digits->GetDigit(it,ip); | |
511 | Int_t idx0=digits->GetTrackID(it,ip,0); | |
512 | Int_t idx1=digits->GetTrackID(it,ip,1); | |
513 | Int_t idx2=digits->GetTrackID(it,ip,2); | |
8ded1b7a | 514 | if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1; |
515 | if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1; | |
516 | if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1; | |
afc42102 | 517 | } while (digits->Next()); |
7f6ddf58 | 518 | for (Int_t j=0; j<np; j++) { |
519 | if (count[j]>1) { | |
520 | if (sec>=digp->GetNInnerSector()) | |
521 | if (row==nrow_up-1 ) good[j]|=0x4000; | |
522 | if (sec>=digp->GetNInnerSector()) | |
523 | if (row==nrow_up-1-gap) good[j]|=0x1000; | |
524 | ||
525 | if (sec>=digp->GetNInnerSector()) | |
526 | if (row==nrow_up-1-shift) good[j]|=0x2000; | |
527 | if (sec>=digp->GetNInnerSector()) | |
528 | if (row==nrow_up-1-gap-shift) good[j]|=0x800; | |
529 | good[j]++; | |
530 | } | |
531 | count[j]=0; | |
532 | } | |
4fa1e629 | 533 | } |
534 | delete[] count; | |
535 | } | |
536 | break; | |
afc42102 | 537 | } |
7f6ddf58 | 538 | |
4fa1e629 | 539 | |
540 | //****** select tracks which are "good" enough | |
541 | AliStack* stack = rl->Stack(); | |
542 | for (i=0; i<np; i++) { | |
543 | if ((good[i]&0x5000) != 0x5000) | |
544 | if ((good[i]&0x2800) != 0x2800) continue; | |
545 | if ((good[i]&0x7FF ) < good_number) continue; | |
546 | ||
547 | TParticle *p = (TParticle*)stack->Particle(i); | |
548 | if (p == 0x0) { | |
88cb7938 | 549 | cerr<<"Can not get particle "<<i<<endl; |
550 | continue; | |
024a7fe9 | 551 | } |
4fa1e629 | 552 | if (p->Pt()<0.100) continue; |
553 | if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue; | |
bf018fe4 | 554 | |
4fa1e629 | 555 | Double_t vx=p->Vx(),vy=p->Vy(); |
556 | if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue; | |
7f6ddf58 | 557 | |
4fa1e629 | 558 | AliTrackReference *ref=new((*refs)[nt]) AliTrackReference(); |
559 | ||
560 | ref->SetLabel(i); | |
561 | Int_t pdg=p->GetPdgCode(); | |
562 | ref->SetLength(pdg); //This will the particle's PDG ! | |
563 | ref->SetMomentum(0.,0.,0.); | |
564 | ref->SetPosition(0.,0.,0.); | |
565 | ||
566 | nt++; | |
567 | } | |
568 | ||
569 | //**** check if there is also information at the entrance of the TPC | |
570 | TTree *TR=rl->TreeTR(); | |
571 | TBranch *branch=TR->GetBranch("TPC"); | |
572 | if (branch==0) { | |
573 | ::Error("GoodTracksTPC","No track references !"); | |
574 | delete rl; | |
575 | return 5; | |
576 | } | |
577 | TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy; | |
578 | branch->SetAddress(&tpcRefs); | |
579 | Int_t nr=(Int_t)TR->GetEntries(); | |
580 | for (Int_t r=0; r<nr; r++) { | |
581 | //cerr<<r<<' '<<nr<<'\r'; | |
582 | TR->GetEvent(r); | |
583 | if (tpcRefs->GetEntriesFast()==0) continue; | |
584 | AliTrackReference *tpcRef=(AliTrackReference*)tpcRefs->UncheckedAt(0); | |
585 | Int_t j; | |
586 | AliTrackReference *ref=0; | |
587 | for (j=0; j<nt; j++) { | |
588 | ref=(AliTrackReference *)refs->UncheckedAt(j); | |
589 | if (ref->Label()==tpcRef->Label()) break; | |
590 | ref=0; | |
591 | } | |
592 | if (ref==0) continue; | |
593 | ||
594 | ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz()); | |
595 | ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z()); | |
596 | ||
597 | tpcRefs->Clear(); | |
598 | } | |
599 | ||
600 | tpcTree.Fill(); | |
601 | ||
602 | delete[] good; | |
603 | } //****** end of the loop over generated events | |
7f6ddf58 | 604 | |
88cb7938 | 605 | |
4fa1e629 | 606 | tpcTree.Write(); |
607 | file->Close(); | |
7f6ddf58 | 608 | |
4fa1e629 | 609 | delete rl; |
610 | return 0; | |
8c555625 | 611 | } |