]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliPerformanceTPC.cxx
extentions by Jochen Thaeder
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliPerformanceTPC.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliPerformanceTPC class. It keeps information from 
3 // comparison of reconstructed and MC particle tracks. In addtion, 
4 // it keeps selection cuts used during comparison. The comparison 
5 // information is stored in the ROOT histograms. Analysis of these 
6 // histograms can be done by using Analyse() class function. The result of 
7 // the analysis (histograms/graphs) are stored in the folder which is
8 // a data member of AliPerformanceTPC.
9 //
10 // Author: J.Otwinowski 04/02/2008 
11 // Changes by M.Knichel 15/10/2010
12 //------------------------------------------------------------------------------
13
14 /*
15  
16   // after running comparison task, read the file, and get component
17   gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
18   LoadMyLibs();
19
20   TFile f("Output.root");
21   AliPerformanceTPC * compObj = (AliPerformanceTPC*)coutput->FindObject("AliPerformanceTPC");
22  
23   // analyse comparison data
24   compObj->Analyse();
25
26   // the output histograms/graphs will be stored in the folder "folderTPC" 
27   compObj->GetAnalysisFolder()->ls("*");
28
29   // user can save whole comparison object (or only folder with anlysed histograms) 
30   // in the seperate output file (e.g.)
31   TFile fout("Analysed_TPC.root","recreate");
32   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
33   fout.Close();
34
35 */
36
37 #include "TCanvas.h"
38 #include "TH1.h"
39 #include "TH2.h"
40 #include "TH3.h"
41 #include "TAxis.h"
42 #include "TPostScript.h"
43 #include "TString.h"
44 #include "TUUID.h"
45 #include "TTree.h"
46 #include "TChain.h"
47 #include "AliTPCPerformanceSummary.h"
48 #include "TSystem.h"
49
50 #include "AliPerformanceTPC.h" 
51 #include "AliESDEvent.h" 
52 #include "AliESDVertex.h"
53 #include "AliESDtrack.h"
54 #include "AliLog.h" 
55 #include "AliMCEvent.h" 
56 #include "AliHeader.h" 
57 #include "AliGenEventHeader.h" 
58 #include "AliStack.h" 
59 #include "AliMCInfoCuts.h" 
60 #include "AliRecInfoCuts.h" 
61 #include "AliTracker.h" 
62 #include "AliTreeDraw.h" 
63 #include "AliTPCTransform.h" 
64 #include "AliTPCseed.h" 
65 #include "AliTPCcalibDB.h" 
66 #include "AliESDfriend.h" 
67 #include "AliESDfriendTrack.h" 
68 #include "AliTPCclusterMI.h" 
69
70 using namespace std;
71
72 ClassImp(AliPerformanceTPC)
73
74 Bool_t AliPerformanceTPC::fgMergeTHnSparse = kFALSE;
75
76
77 //_____________________________________________________________________________
78 /*
79 AliPerformanceTPC::AliPerformanceTPC():
80   AliPerformanceObject("AliPerformanceTPC"),
81   fTPCClustHisto(0),
82   fTPCEventHisto(0),
83   fTPCTrackHisto(0),
84   fFolderObj(0),
85
86   // Cuts 
87   fCutsRC(0),  
88   fCutsMC(0),  
89
90   // histogram folder 
91   fAnalysisFolder(0),
92   
93   fUseHLT(kFALSE)
94
95 {
96   Init();
97 }
98 */
99
100 //_____________________________________________________________________________
101 AliPerformanceTPC::AliPerformanceTPC(Char_t* name, Char_t* title,Int_t analysisMode,Bool_t hptGenerator, Int_t run):
102   AliPerformanceObject(name,title,run),
103   fTPCClustHisto(0),
104   fTPCEventHisto(0),
105   fTPCTrackHisto(0),
106   fFolderObj(0),
107
108   // Cuts 
109   fCutsRC(0),  
110   fCutsMC(0),  
111
112   // histogram folder 
113   fAnalysisFolder(0),
114   
115   fUseHLT(kFALSE)
116
117 {
118   // named constructor  
119   // 
120   SetAnalysisMode(analysisMode);
121   SetHptGenerator(hptGenerator);
122
123   Init();
124 }
125
126
127 //_____________________________________________________________________________
128 AliPerformanceTPC::~AliPerformanceTPC()
129 {
130   // destructor
131    
132   if(fTPCClustHisto) delete fTPCClustHisto; fTPCClustHisto=0;     
133   if(fTPCEventHisto) delete fTPCEventHisto; fTPCEventHisto=0;     
134   if(fTPCTrackHisto) delete fTPCTrackHisto; fTPCTrackHisto=0;   
135   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
136   if(fFolderObj) delete fFolderObj; fFolderObj=0;
137 }
138
139
140 //_____________________________________________________________________________
141 void AliPerformanceTPC::Init()
142 {
143   //
144   // histogram bining
145   //
146
147   // set pt bins
148   Int_t nPtBins = 50;
149   Double_t ptMin = 1.e-2, ptMax = 10.;
150
151   Double_t *binsPt = 0;
152   if (IsHptGenerator())  { 
153     nPtBins = 100; ptMax = 100.;
154     binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
155   } else {
156     binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
157   }
158
159   /*
160   Int_t nPtBins = 31;
161   Double_t binsPt[32] = {0.,0.05,0.1,0.15,0.2,0.25,0.3,0.35,0.4,0.45,0.5,0.55,0.6,0.7,0.8,0.9,1.0,1.2,1.4,1.6,1.8,2.0,2.25,2.5,2.75,3.,3.5,4.,5.,6.,8.,10.};
162   Double_t ptMin = 0., ptMax = 10.; 
163
164   if(IsHptGenerator() == kTRUE) {
165     nPtBins = 100;
166     ptMin = 0.; ptMax = 100.; 
167   }
168   */
169   //
170
171   //
172   //padRow:phi:TPCSide:pad:detector
173   Int_t binsTPCClustHisto[5] =   {160,  180,  2, 256, 512};
174   Double_t minTPCClustHisto[5] = {0.,   0.,   0., -128, 0};
175   Double_t maxTPCClustHisto[5] = {160., 2.*TMath::Pi(), 2., 128, 512};
176
177   fTPCClustHisto = new THnSparseF("fTPCClustHisto","padRow:phi:TPCSide:pad:detector",5,binsTPCClustHisto,minTPCClustHisto,maxTPCClustHisto);
178   fTPCClustHisto->GetAxis(0)->SetTitle("padRow");
179   fTPCClustHisto->GetAxis(1)->SetTitle("phi (rad)");
180   fTPCClustHisto->GetAxis(2)->SetTitle("TPCSide");
181   fTPCClustHisto->GetAxis(3)->SetTitle("pad");
182   fTPCClustHisto->GetAxis(4)->SetTitle("detector");
183   //fTPCClustHisto->Sumw2();
184   
185   Int_t maxMult;
186   if (fHighMultiplicity) { maxMult = 4001; } else { maxMult = 151; }
187   // Xv:Yv:Zv:mult:multP:multN:vertStatus
188   Int_t binsTPCEventHisto[7]=  {100,  100,   100,  maxMult,  maxMult,  maxMult, 2   };
189   Double_t minTPCEventHisto[7]={-10., -10., -30., -0.5,  -0.5,  -0.5, -0.5  };
190   Double_t maxTPCEventHisto[7]={ 10.,  10.,  30.,  maxMult-0.5,  maxMult-0.5, maxMult-0.5, 1.5 };
191
192   fTPCEventHisto = new THnSparseF("fTPCEventHisto","Xv:Yv:Zv:mult:multP:multN:vertStatus",7,binsTPCEventHisto,minTPCEventHisto,maxTPCEventHisto);
193   fTPCEventHisto->GetAxis(0)->SetTitle("Xv (cm)");
194   fTPCEventHisto->GetAxis(1)->SetTitle("Yv (cm)");
195   fTPCEventHisto->GetAxis(2)->SetTitle("Zv (cm)");
196   fTPCEventHisto->GetAxis(3)->SetTitle("mult");
197   fTPCEventHisto->GetAxis(4)->SetTitle("multP");
198   fTPCEventHisto->GetAxis(5)->SetTitle("multN");
199   fTPCEventHisto->GetAxis(6)->SetTitle("vertStatus");
200   //fTPCEventHisto->Sumw2();
201
202
203   // nTPCClust:chi2PerTPCClust:nTPCClustFindRatio:DCAr:DCAz:eta:phi:pt:charge:vertStatus
204    Int_t binsTPCTrackHisto[10]=  { 160,  20,  60,  30, 30,  30,   144,             nPtBins,   3, 2 };
205    Double_t minTPCTrackHisto[10]={ 0.,   0.,  0., -3,  -3., -1.5, 0.,             ptMin,   -1.5, -0.5 };
206    Double_t maxTPCTrackHisto[10]={ 160., 5., 1.2, 3,   3.,  1.5, 2.*TMath::Pi(), ptMax,    1.5,  1.5 };
207   
208   // nTPCClust:chi2PerTPCClust:nTPCClustFindRatio:DCAr:DCAz:eta:phi:pt:charge:vertStatus
209 //   Int_t binsTPCTrackHisto[10]=  { 160,  50,  60,  100, 100,  30,   144,            nPtBins,    3, 2 };
210 //   Double_t minTPCTrackHisto[10]={ 0.,   0.,  0., -10,  -10., -1.5, 0.,             ptMin,   -1.5, 0 };
211 //   Double_t maxTPCTrackHisto[10]={ 160., 10., 1.2, 10,   10.,  1.5, 2.*TMath::Pi(), ptMax,    1.5, 2 };
212
213
214
215   fTPCTrackHisto = new THnSparseF("fTPCTrackHisto","nClust:chi2PerClust:nClust/nFindableClust:DCAr:DCAz:eta:phi:pt:charge:vertStatus",10,binsTPCTrackHisto,minTPCTrackHisto,maxTPCTrackHisto);
216   fTPCTrackHisto->SetBinEdges(7,binsPt);
217
218   fTPCTrackHisto->GetAxis(0)->SetTitle("nClust");
219   fTPCTrackHisto->GetAxis(1)->SetTitle("chi2PerClust");
220   fTPCTrackHisto->GetAxis(2)->SetTitle("nClust/nFindableClust");
221   fTPCTrackHisto->GetAxis(3)->SetTitle("DCAr (cm)");
222   fTPCTrackHisto->GetAxis(4)->SetTitle("DCAz (cm)");
223   fTPCTrackHisto->GetAxis(5)->SetTitle("#eta");
224   fTPCTrackHisto->GetAxis(6)->SetTitle("#phi (rad)");
225   fTPCTrackHisto->GetAxis(7)->SetTitle("p_{T} (GeV/c)");
226   fTPCTrackHisto->GetAxis(8)->SetTitle("charge");
227   fTPCTrackHisto->GetAxis(9)->SetTitle("vertStatus");
228   //fTPCTrackHisto->Sumw2();
229
230   // Init cuts 
231   if(!fCutsMC) {
232     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
233   }
234   if(!fCutsRC) {
235     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object"); 
236   }
237
238   // init folder
239   fAnalysisFolder = CreateFolder("folderTPC","Analysis Resolution Folder");
240   
241 }
242
243
244 //_____________________________________________________________________________
245 void AliPerformanceTPC::ProcessTPC(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent *const esdEvent, Bool_t vertStatus)
246 {
247 //
248 // fill TPC QA info
249 //
250   if(!esdEvent) return;
251   if(!esdTrack) return;
252
253   if( IsUseTrackVertex() ) 
254   { 
255     // Relate TPC inner params to prim. vertex
256     const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
257     Double_t x[3]; esdTrack->GetXYZ(x);
258     Double_t b[3]; AliTracker::GetBxByBz(x,b);
259     Bool_t isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
260     if(!isOK) return;
261
262     /*
263       // JMT -- recaluclate DCA for HLT if not present
264       if ( dca[0] == 0. && dca[1] == 0. ) {
265         track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
266       }
267     */
268   }
269
270   // Fill TPC only resolution comparison information 
271   const AliExternalTrackParam *track = esdTrack->GetTPCInnerParam();
272   if(!track) return;
273
274   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
275   esdTrack->GetImpactParametersTPC(dca,cov);
276
277   Float_t q = esdTrack->Charge();
278   Float_t pt = track->Pt();
279   Float_t eta = track->Eta();
280   Float_t phi = track->Phi();
281   Int_t nClust = esdTrack->GetTPCclusters(0);
282   Int_t nFindableClust = esdTrack->GetTPCNclsF();
283
284   Float_t chi2PerCluster = 0.;
285   if(nClust>0.) chi2PerCluster = esdTrack->GetTPCchi2()/Float_t(nClust);
286
287   Float_t clustPerFindClust = 0.;
288   if(nFindableClust>0.) clustPerFindClust = Float_t(nClust)/nFindableClust;
289   
290   //
291   // select primaries
292   //
293   Double_t dcaToVertex = -1;
294   if( fCutsRC->GetDCAToVertex2D() ) 
295   {
296       dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
297   }
298   if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
299   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
300   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
301
302   Double_t vTPCTrackHisto[10] = {nClust,chi2PerCluster,clustPerFindClust,dca[0],dca[1],eta,phi,pt,q,vertStatus};
303   fTPCTrackHisto->Fill(vTPCTrackHisto); 
304  
305   //
306   // Fill rec vs MC information
307   //
308   if(!stack) return;
309
310 }
311
312
313 //_____________________________________________________________________________
314 void AliPerformanceTPC::ProcessTPCITS(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent* const esdEvent, Bool_t vertStatus)
315 {
316   // Fill comparison information (TPC+ITS) 
317   if(!esdTrack) return;
318   if(!esdEvent) return;
319
320   if( IsUseTrackVertex() ) 
321   { 
322     // Relate TPC inner params to prim. vertex
323     const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
324     Double_t x[3]; esdTrack->GetXYZ(x);
325     Double_t b[3]; AliTracker::GetBxByBz(x,b);
326     Bool_t isOK = esdTrack->RelateToVertexBxByBz(vtxESD, b, kVeryBig);
327     if(!isOK) return;
328
329     /*
330       // JMT -- recaluclate DCA for HLT if not present
331       if ( dca[0] == 0. && dca[1] == 0. ) {
332         track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
333       }
334     */
335   }
336
337   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
338   esdTrack->GetImpactParameters(dca,cov);
339
340   if ((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0) return; // ITS refit
341   if ((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return; // TPC refit
342   if (esdTrack->GetITSclusters(0)<fCutsRC->GetMinNClustersITS()) return;  // min. nb. ITS clusters
343
344   Float_t q = esdTrack->Charge();
345   Float_t pt = esdTrack->Pt();
346   Float_t eta = esdTrack->Eta();
347   Float_t phi = esdTrack->Phi();
348   Int_t nClust = esdTrack->GetTPCclusters(0);
349   Int_t nFindableClust = esdTrack->GetTPCNclsF();
350
351   Float_t chi2PerCluster = 0.;
352   if(nClust>0.) chi2PerCluster = esdTrack->GetTPCchi2()/Float_t(nClust);
353
354   Float_t clustPerFindClust = 0.;
355   if(nFindableClust>0.) clustPerFindClust = Float_t(nClust)/nFindableClust;
356   
357   //
358   // select primaries
359   //
360   Double_t dcaToVertex = -1;
361   if( fCutsRC->GetDCAToVertex2D() ) 
362   {
363       dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
364   }
365   if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
366   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
367   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
368
369   Double_t vTPCTrackHisto[10] = {nClust,chi2PerCluster,clustPerFindClust,dca[0],dca[1],eta,phi,pt,q,vertStatus};
370   fTPCTrackHisto->Fill(vTPCTrackHisto); 
371  
372   //
373   // Fill rec vs MC information
374   //
375   if(!stack) return;
376 }
377
378
379 //_____________________________________________________________________________
380 void AliPerformanceTPC::ProcessConstrained(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/, AliESDEvent* const /*esdEvent*/)
381 {
382   // Fill comparison information (constarained parameters) 
383   AliDebug(AliLog::kWarning, "Warning: Not implemented");
384 }
385
386
387 //_____________________________________________________________________________
388 void AliPerformanceTPC::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
389 {
390   // Process comparison information 
391   //
392   if(!esdEvent) 
393   {
394     Error("Exec","esdEvent not available");
395     return;
396   }
397   AliHeader* header = 0;
398   AliGenEventHeader* genHeader = 0;
399   AliStack* stack = 0;
400   TArrayF vtxMC(3);
401   
402   if(bUseMC)
403   {
404     if(!mcEvent) {
405       Error("Exec","mcEvent not available");
406       return;
407     }
408     // get MC event header
409     header = mcEvent->Header();
410     if (!header) {
411       Error("Exec","Header not available");
412       return;
413     }
414     // MC particle stack
415     stack = mcEvent->Stack();
416     if (!stack) {
417       Error("Exec","Stack not available");
418       return;
419     }
420     // get MC vertex
421     genHeader = header->GenEventHeader();
422     if (!genHeader) {
423       Error("Exec","Could not retrieve genHeader from Header");
424       return;
425     }
426     genHeader->PrimaryVertex(vtxMC);
427   } 
428   
429   // trigger
430   if(!bUseMC &&GetTriggerClass()) {
431     Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
432     if(!isEventTriggered) return; 
433   }
434
435   // get TPC event vertex
436   const AliESDVertex *vtxESD = NULL; 
437   if(fUseTrackVertex) {
438     vtxESD = esdEvent->GetPrimaryVertexTracks();
439   } else {
440     vtxESD = esdEvent->GetPrimaryVertexTPC();
441   }
442   if(!vtxESD) return;
443
444   //  events with rec. vertex
445   Int_t mult=0; Int_t multP=0; Int_t multN=0;
446   
447   // changed to take all events but store vertex status
448 //  if(vtxESD->GetStatus() >0)
449 //  {
450   // store vertex status
451   Bool_t vertStatus = vtxESD->GetStatus();
452   //  Process ESD events
453   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
454   { 
455     AliESDtrack *track = esdEvent->GetTrack(iTrack);
456     if(!track) continue;
457
458     // if not fUseKinkDaughters don't use tracks with kink index > 0
459     if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;
460     
461     if(bUseESDfriend && esdFriend && esdFriend->TestSkipBit()==kFALSE) 
462     {
463       AliESDfriendTrack *friendTrack=esdFriend->GetTrack(iTrack);
464       if(friendTrack) 
465       {
466         //
467         TObject *calibObject=0;
468         AliTPCseed *seed=0;
469         for (Int_t j=0;(calibObject=friendTrack->GetCalibObject(j));++j) {
470             if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) {
471             break;
472           }
473         }
474
475         // 
476         for (Int_t irow=0;irow<159;irow++) {
477         if(!seed) continue;
478           
479           AliTPCclusterMI *cluster=seed->GetClusterPointer(irow);
480           if (!cluster) continue;
481
482              Float_t gclf[3];
483              cluster->GetGlobalXYZ(gclf);
484
485              //Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
486              //Int_t i[1]={cluster->GetDetector()};
487              //transform->Transform(x,i,0,1);
488              //printf("gx %f gy  %f  gz %f \n", cluster->GetX(), cluster->GetY(),cluster->GetZ());
489              //printf("gclf[0] %f gclf[1]  %f  gclf[2] %f \n", gclf[0], gclf[1],  gclf[2]);
490      
491              Int_t TPCside; 
492              if(gclf[2]>0.) TPCside=0; // A side 
493              else TPCside=1;
494
495              //
496              //Double_t vTPCClust1[3] = { gclf[0], gclf[1],  TPCside };
497              //fTPCClustHisto1->Fill(vTPCClust1);
498
499              //  
500              Float_t pad = cluster->GetPad();
501              Int_t detector = cluster->GetDetector();
502              Double_t phi = TMath::ATan2(gclf[1],gclf[0]);
503              if(phi < 0) phi += 2.*TMath::Pi();
504             
505              Double_t vTPCClust[5] = { irow, phi, TPCside, pad, detector };
506              fTPCClustHisto->Fill(vTPCClust);
507         }
508       }
509     }
510
511     if(GetAnalysisMode() == 0) ProcessTPC(stack,track,esdEvent,vertStatus);
512     else if(GetAnalysisMode() == 1) ProcessTPCITS(stack,track,esdEvent,vertStatus);
513     else if(GetAnalysisMode() == 2) ProcessConstrained(stack,track,esdEvent);
514     else {
515       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
516       return;
517     }
518
519     // TPC only
520     if(!fUseHLT){
521       AliESDtrack *tpcTrack = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
522       if(!tpcTrack) continue;
523       
524       // track selection
525       if( fCutsRC->AcceptTrack(tpcTrack) ) { 
526         mult++;
527         if(tpcTrack->Charge()>0.) multP++;
528         if(tpcTrack->Charge()<0.) multN++;
529       }
530       
531       if(tpcTrack) delete tpcTrack;
532     }
533     else
534       if( fCutsRC->AcceptTrack(track) ) { 
535         //Printf("Still here for HLT");
536         mult++;
537         if(track->Charge()>0.) multP++;
538         if(track->Charge()<0.) multN++;
539       }
540     
541   }
542 //  }
543   //
544   
545   Double_t vTPCEvent[7] = {vtxESD->GetXv(),vtxESD->GetYv(),vtxESD->GetZv(),mult,multP,multN,vtxESD->GetStatus()};
546   fTPCEventHisto->Fill(vTPCEvent);
547 }
548
549
550 //_____________________________________________________________________________
551 void AliPerformanceTPC::Analyse()
552 {
553     //
554     // Analyse comparison information and store output histograms
555     // in the folder "folderTPC"
556     //
557     TH1::AddDirectory(kFALSE);
558     TH1::SetDefaultSumw2(kFALSE);
559     TObjArray *aFolderObj = new TObjArray;
560     //aFolderObj->SetOwner(); // objects are owned by fanalysisFolder
561     TString selString;
562
563     //
564     // Cluster histograms
565     //
566     AddProjection(aFolderObj, "clust", fTPCClustHisto, 0, 1, 2);
567     
568     selString = "all";
569     for(Int_t i=0; i <= 2; i++) {
570         AddProjection(aFolderObj, "clust", fTPCClustHisto, i, &selString);
571     }
572     
573     //fTPCClustHisto->GetAxis(2)->SetRange(1,1); // A-side
574     //selString = "A_side";
575     //AddProjection(aFolderObj, fTPCClustHisto, 0, 1, &selString);
576     
577     //fTPCClustHisto->GetAxis(2)->SetRange(2,2); // C-side
578     //selString = "C_side";
579     //AddProjection(aFolderObj, fTPCClustHisto, 0, 1, &selString);
580     
581     //reset range
582     fTPCClustHisto->GetAxis(2)->SetRange(1,2); 
583     
584     //
585     // event histograms
586     //
587     for(Int_t i=0; i<=6; i++) {
588       AddProjection(aFolderObj, "event", fTPCEventHisto, i);
589     }    
590     AddProjection(aFolderObj, "event", fTPCEventHisto, 4, 5);
591     AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 1);
592     AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 3);
593     AddProjection(aFolderObj, "event", fTPCEventHisto, 1, 3);
594     AddProjection(aFolderObj, "event", fTPCEventHisto, 2, 3);
595
596     // reconstructed vertex status > 0
597     fTPCEventHisto->GetAxis(6)->SetRange(2,2);
598     selString = "recVertex";
599     for(Int_t i=0; i<=5; i++) {
600       AddProjection(aFolderObj, "event", fTPCEventHisto, i, &selString);
601     }
602     AddProjection(aFolderObj, "event", fTPCEventHisto, 4, 5, &selString);
603     AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 1, &selString);
604     AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 3, &selString);
605     AddProjection(aFolderObj, "event", fTPCEventHisto, 1, 3, &selString);
606     AddProjection(aFolderObj, "event", fTPCEventHisto, 2, 3, &selString);
607
608     // reset cuts
609     fTPCEventHisto->GetAxis(6)->SetRange(1,2);
610
611     //
612     // Track histograms 
613     // 
614     // all with vertex
615     fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,1.5);
616     fTPCTrackHisto->GetAxis(9)->SetRangeUser(0.5,1.5);
617     selString = "all_recVertex";
618     for(Int_t i=0; i <= 9; i++) {
619         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);        
620     }
621     for(Int_t i=0; i <= 4; i++) {
622         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, 5, 7, &selString);        
623     }    
624
625
626
627     // Track histograms (pos with vertex)
628     fTPCTrackHisto->GetAxis(8)->SetRangeUser(0,1.5);
629     selString = "pos_recVertex";
630     for(Int_t i=0; i <= 9; i++) {
631         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);
632     }
633     for(Int_t i=0; i <= 4; i++) { for(Int_t j=5; j <= 5; j++) { for(Int_t k=j+1; k <= 7; k++) {
634         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, j, k, &selString);
635     }  }  }
636     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 2, &selString);
637     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 5, &selString);
638     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 2, 5, &selString);
639     AddProjection(aFolderObj, "track", fTPCTrackHisto, 1, 2, 5, &selString);
640     AddProjection(aFolderObj, "track", fTPCTrackHisto, 3, 4, 5, &selString);
641     AddProjection(aFolderObj, "track", fTPCTrackHisto, 5, 6, 7, &selString);
642   
643     // Track histograms (neg with vertex)
644     fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,0);
645     selString = "neg_recVertex";
646     for(Int_t i=0; i <= 9; i++) {
647         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);
648     }
649     for(Int_t i=0; i <= 4; i++) { for(Int_t j=5; j <= 5; j++) { for(Int_t k=j+1; k <= 7; k++) {
650         AddProjection(aFolderObj, "track", fTPCTrackHisto, i, j, k, &selString);
651     }  }  }
652     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 2, &selString);
653     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 5, &selString);
654     AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 2, 5, &selString);
655     AddProjection(aFolderObj, "track", fTPCTrackHisto, 1, 2, 5, &selString);
656     AddProjection(aFolderObj, "track", fTPCTrackHisto, 3, 4, 5, &selString);
657
658     //restore cuts
659     fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,1.5);
660     fTPCTrackHisto->GetAxis(9)->SetRangeUser(-0.5,1.5);
661   
662   
663     printf("exportToFolder\n");
664     // export objects to analysis folder
665     fAnalysisFolder = ExportToFolder(aFolderObj);
666     if (fFolderObj) delete fFolderObj;
667     fFolderObj = aFolderObj;
668     aFolderObj=0;
669 }
670
671
672 //_____________________________________________________________________________
673 TFolder* AliPerformanceTPC::ExportToFolder(TObjArray * array) 
674 {
675   // recreate folder avery time and export objects to new one
676   //
677   AliPerformanceTPC * comp=this;
678   TFolder *folder = comp->GetAnalysisFolder();
679
680   TString name, title;
681   TFolder *newFolder = 0;
682   Int_t i = 0;
683   Int_t size = array->GetSize();
684
685   if(folder) { 
686      // get name and title from old folder
687      name = folder->GetName();  
688      title = folder->GetTitle();  
689
690          // delete old one
691      delete folder;
692
693          // create new one
694      newFolder = CreateFolder(name.Data(),title.Data());
695      newFolder->SetOwner();
696
697          // add objects to folder
698      while(i < size) {
699            newFolder->Add(array->At(i));
700            i++;
701          }
702   }
703
704 return newFolder;
705 }
706
707 //_____________________________________________________________________________
708 Long64_t AliPerformanceTPC::Merge(TCollection* const list) 
709 {
710   // Merge list of objects (needed by PROOF)
711
712   if (!list)
713   return 0;
714
715   if (list->IsEmpty())
716   return 1;
717
718   TIterator* iter = list->MakeIterator();
719   TObject* obj = 0;
720   TObjArray* objArrayList = 0;
721   objArrayList = new TObjArray();
722
723   // collection of generated histograms
724   Int_t count=0;
725   while((obj = iter->Next()) != 0) 
726   {
727     AliPerformanceTPC* entry = dynamic_cast<AliPerformanceTPC*>(obj);
728     if (entry == 0) continue; 
729     if (fgMergeTHnSparse) {
730         if ((fTPCClustHisto) && (entry->fTPCClustHisto)) { fTPCClustHisto->Add(entry->fTPCClustHisto); }
731         if ((fTPCEventHisto) && (entry->fTPCEventHisto)) { fTPCEventHisto->Add(entry->fTPCEventHisto); }
732         if ((fTPCTrackHisto) && (entry->fTPCTrackHisto)) { fTPCTrackHisto->Add(entry->fTPCTrackHisto); }
733     }
734     // the analysisfolder is only merged if present
735     if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
736
737     count++;
738   }
739   if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
740   // to signal that track histos were not merged: reset
741   if (!fgMergeTHnSparse) { fTPCTrackHisto->Reset(); fTPCClustHisto->Reset(); fTPCEventHisto->Reset(); }
742   // delete
743   if (objArrayList)  delete objArrayList;  objArrayList=0;
744 return count;
745 }
746
747
748 //_____________________________________________________________________________
749 TFolder* AliPerformanceTPC::CreateFolder(TString name, TString title) 
750
751 // create folder for analysed histograms
752 //
753 TFolder *folder = 0;
754   folder = new TFolder(name.Data(),title.Data());
755
756   return folder;
757 }
758
759 //_____________________________________________________________________________
760 TTree* AliPerformanceTPC::CreateSummary()
761 {
762     // implementaion removed, switched back to use AliPerformanceSummary (now called in AliPerformanceTask)
763     return 0;
764 }
765