]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceMatch.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGPP / TPC / AliPerformanceMatch.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliPerformanceMatch class. It checks the track matching 
3 // quality (residuals, pulls) bewteen tracking detectors TPC-ITS and TPC-TRD 
4 // with TPC as the reference detector. In addition, the ITS and TRD track 
5 // reconstruction efficiency with respect to TPC is calculated.
6 //
7 // The matchinig quality parameters are stored in the THnSparse histograms. 
8 // The analysis of these histograms to extract reference information is done in 
9 // the Analyse() function (post-processing).  
10 //  
11 // The histograms with reference information can be exported to the ROOT folder.
12 //
13 // Author: J.Otwinowski 17/10/2009 
14 // Changes by M.Knichel 22/10/2010 
15 //------------------------------------------------------------------------------
16
17 /*
18  
19   // after running comparison task, read the file, and get component
20   gROOT->LoadMacro("$ALICE_ROOT/PWGPP/Macros/LoadMyLibs.C");
21   LoadMyLibs();
22
23   TFile f("Output.root");
24   AliPerformanceMatch * compObj = (AliPerformanceMatch*)coutput->FindObject("AliPerformanceMatchTPCTRD");
25  
26   // analyse comparison data
27   compObj->Analyse();
28
29   // the output histograms/graphs will be stored in the folder "folderMatch" 
30   compObj->GetAnalysisFolder()->ls("*");
31
32   // user can save whole comparison object (or only folder with anlysed histograms) 
33   // in the seperate output file (e.g.)
34   TFile fout("Analysed_Match.root","recreate");
35   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
36   fout.Close();
37
38 */
39
40 #include "TCanvas.h"
41 #include "TH1.h"
42 #include "TH2.h"
43 #include "TAxis.h"
44 #include "TF1.h"
45
46 #include "AliPerformanceMatch.h" 
47 #include "AliESDEvent.h" 
48 #include "AliESDVertex.h"
49 #include "AliESDtrack.h"
50 #include "AliESDfriendTrack.h"
51 #include "AliESDfriend.h"
52 #include "AliLog.h" 
53 #include "AliMCEvent.h" 
54 #include "AliMCParticle.h" 
55 #include "AliHeader.h" 
56 #include "AliGenEventHeader.h" 
57 #include "AliStack.h" 
58 #include "AliMCInfoCuts.h" 
59 #include "AliRecInfoCuts.h" 
60 #include "AliTracker.h" 
61 #include "AliTRDtrackV1.h" 
62 #include "AliTreeDraw.h" 
63
64 using namespace std;
65
66 ClassImp(AliPerformanceMatch)
67
68 Bool_t AliPerformanceMatch::fgMergeTHnSparse = kFALSE;
69 Bool_t AliPerformanceMatch::fgUseMergeTHnSparse = kFALSE;
70
71 //_____________________________________________________________________________
72 AliPerformanceMatch::AliPerformanceMatch(const Char_t* name, const Char_t* title, Int_t analysisMode, Bool_t hptGenerator):
73   AliPerformanceObject(name,title),
74   fResolHisto(0),
75   fPullHisto(0),
76   fTrackingEffHisto(0),
77   fTPCConstrain(0),
78   fFolderObj(0),
79   // Cuts 
80   fCutsRC(0),  
81   fCutsMC(0),  
82
83   // histogram folder 
84   fAnalysisFolder(0),
85   fUseHLT(0)
86 {
87   // named constructor  
88   // 
89   SetAnalysisMode(analysisMode);
90   SetHptGenerator(hptGenerator);
91
92   Init();
93 }
94
95 //_____________________________________________________________________________
96 AliPerformanceMatch::~AliPerformanceMatch()
97 {
98   // destructor
99    
100   if(fResolHisto) delete fResolHisto; fResolHisto=0;     
101   if(fPullHisto)  delete fPullHisto;  fPullHisto=0;
102   if(fTrackingEffHisto) delete fTrackingEffHisto; fTrackingEffHisto = 0x0;
103   if(fTPCConstrain) delete fTPCConstrain; fTPCConstrain = 0x0;
104   
105   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
106   if(fFolderObj) delete fFolderObj; fFolderObj=0;
107 }
108
109 //_____________________________________________________________________________
110 void AliPerformanceMatch::Init(){
111
112   //
113   // Make performance histogrms
114   //
115
116   // set pt bins
117   Int_t nPtBins = 50;
118   Double_t ptMin = 1.e-2, ptMax = 20.;
119  
120   Double_t *binsPt = 0;
121
122   if (IsHptGenerator())  { 
123         ptMax = 100.;
124   } 
125    binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
126
127
128   // Double_t yMin = -0.02, yMax = 0.02;
129   // Double_t zMin = -12.0, zMax = 12.0;
130   // if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) { 
131   //   yMin = -100.; yMax = 100.; 
132   //   zMin = -100.; zMax = 100.; 
133   // }
134
135   //
136   //init ITS TPC Mactching
137   //
138   if(GetAnalysisMode()==1){
139     // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
140     Int_t binsResolHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
141     Double_t minResolHisto[9]={-1.,-1.,-0.03,-0.03,-0.2, 0., -1.5, ptMin,0};
142     Double_t maxResolHisto[9]={ 1., 1., 0.03, 0.03, 0.2, 2.*TMath::Pi(), 1.5, ptMax,2};
143     
144     fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:phi:eta:pt:isRec",9,binsResolHisto,minResolHisto,maxResolHisto);
145     fResolHisto->SetBinEdges(7,binsPt);
146     
147     fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
148     fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
149     fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
150     fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
151     fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
152     fResolHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
153     fResolHisto->GetAxis(6)->SetTitle("#eta_{ref}");
154     fResolHisto->GetAxis(7)->SetTitle("p_{Tref} (GeV/c)");
155     fResolHisto->GetAxis(8)->SetTitle("isReconstructed");
156     fResolHisto->Sumw2();
157     
158     //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
159     Int_t binsPullHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
160     Double_t minPullHisto[9]={-5.,-5.,-5.,-5.,-5., 0.,-1.5, ptMin,0};
161     Double_t maxPullHisto[9]={ 5., 5., 5., 5., 5., 2.*TMath::Pi(), 1.5, ptMax,2};
162     fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:phi:eta:1pt:isRec",9,binsPullHisto,minPullHisto,maxPullHisto);
163     
164     fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
165     fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
166     fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
167     fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
168     fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
169     fPullHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
170     fPullHisto->GetAxis(6)->SetTitle("#eta_{ref}");
171     fPullHisto->GetAxis(7)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
172     fPullHisto->GetAxis(8)->SetTitle("isReconstructed");
173     fPullHisto->Sumw2();
174   }
175   
176   //
177   //TPC  ITS  Mactching
178   // 
179   if(GetAnalysisMode()==0){
180     // -> has match:y:z:snp:tgl:phi:pt:ITSclusters
181     Int_t binsTrackingEffHisto[5]    = { 2,   90,       nPtBins, 30,    7   };
182     Double_t minTrackingEffHisto[5]  = {-0.5, 0.,          ptMin, -1.5, -0.5 };
183     Double_t maxTrackingEffHisto[5]  = { 1.5, 2*TMath::Pi(), ptMax, 1.5,  6.5 };
184     
185     fTrackingEffHisto = new THnSparseF("fTrackingEffHisto","has match:phi:pt:eta:ITSclusters",5,binsTrackingEffHisto,minTrackingEffHisto,maxTrackingEffHisto);
186     fTrackingEffHisto->SetBinEdges(2,binsPt);
187     fTrackingEffHisto->GetAxis(0)->SetTitle("IsMatching");
188     fTrackingEffHisto->GetAxis(1)->SetTitle("phi (rad)");
189     fTrackingEffHisto->GetAxis(2)->SetTitle("p_{T}");
190     fTrackingEffHisto->GetAxis(3)->SetTitle("eta");
191     fTrackingEffHisto->GetAxis(4)->SetTitle("number of ITS clusters");
192     fTrackingEffHisto->Sumw2();
193   }
194
195   //
196   //TPC constraining to vertex
197   //
198   if(GetAnalysisMode()==2){
199     //initilization of fTPCConstrain
200     Int_t  binsTPCConstrain[4] = {100, 90,            nPtBins, 30};
201     Double_t minTPCConstrain[4] = {-5, 0,             ptMin,   -1.5};
202     Double_t maxTPCConstrain[4] = {5,  2*TMath::Pi(), ptMax,  1.5};
203     
204     fTPCConstrain = new THnSparseF("fTPCConstrain","pull_phi:phi:pt:eta",4, binsTPCConstrain,minTPCConstrain,maxTPCConstrain);
205     fTPCConstrain->SetBinEdges(2,binsPt);
206     fTPCConstrain->GetAxis(0)->SetTitle("(#phi-#phi_{ref})/#sigma");
207     fTPCConstrain->GetAxis(1)->SetTitle("phi (rad)");
208     fTPCConstrain->GetAxis(2)->SetTitle("p_{T}");
209     fTPCConstrain->GetAxis(3)->SetTitle("eta");
210     fTPCConstrain->Sumw2();
211   }
212
213   // Init cuts 
214   if(!fCutsMC) 
215     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
216   if(!fCutsRC) 
217     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
218
219   // init folder
220   fAnalysisFolder = CreateFolder("folderMatch","Analysis Matching Folder");
221   
222    // save merge status in object
223   fMergeTHnSparseObj = fgMergeTHnSparse;
224
225   delete [] binsPt;
226
227 }
228
229 //_____________________________________________________________________________
230 void AliPerformanceMatch::ProcessITSTPC(Int_t iTrack, AliESDEvent *const esdEvent, AliStack* /*const stack*/, AliESDtrack *const esdTrack)
231 {
232   //
233   // addition to standard analysis - check if ITS stand-alone tracks have a match in the TPC
234   // Origin: A. Kalwait
235   // Modified: J. Otwinowski
236   if(!esdEvent) return;
237   if(!esdTrack) return;
238   //  if(!esdFriendTrack) return;
239
240   // ITS stand alone tracks with SPD layers 
241   if(!(esdTrack->GetStatus() & AliESDtrack::kITSpureSA)) return;
242   if(!(esdTrack->GetStatus() & AliESDtrack::kITSrefit)) return;
243   if(esdTrack->GetNcls(0)<4) return;
244    if(!esdTrack->HasPointOnITSLayer(0) && !esdTrack->HasPointOnITSLayer(1)) return;
245   
246   const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
247   //   const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
248  AliESDtrack* tpcTrack2 = NULL;
249   Bool_t hasMatch = kFALSE;
250   for (Int_t jTrack = 0; jTrack < esdEvent->GetNumberOfTracks(); jTrack++) {
251     // loop for all those tracks and check if the corresponding TPC track is found
252     if (jTrack==iTrack) continue;
253     AliESDtrack *trackTPC = esdEvent->GetTrack(jTrack);
254     if (!trackTPC) continue;
255     if (!trackTPC->GetTPCInnerParam()) continue;
256     if(!(trackTPC->GetStatus() & AliESDtrack::kTPCrefit)) continue;
257     
258     // TPC nClust/track after first tracking pass
259     // if(trackTPC->GetTPCNclsIter1()<fCutsRC->GetMinNClustersTPC()) continue;
260     tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);
261     if(!tpcTrack2) continue;
262     if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; tpcTrack2=0; continue; } 
263     
264      if(!fCutsRC->AcceptTrack(tpcTrack2)) { delete tpcTrack2; tpcTrack2=0; continue; }
265     // check matching
266     if (TMath::Abs(esdTrack->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; tpcTrack2=0; continue; }
267     if (TMath::Abs(esdTrack->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
268     if (TMath::Abs(esdTrack->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
269     
270     hasMatch=kTRUE;
271     break;
272   }
273   
274   FillHistograms(tpcTrack2,esdTrack,hasMatch);     
275   delete tpcTrack2;
276 }
277
278 //_____________________________________________________________________________
279 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack)
280 {
281   //
282   // Match TPC and ITS min-bias tracks
283   // at radius between detectors
284   //
285   if(!esdTrack) return;
286   //   if(!esdFriendTrack) return;
287    
288   Bool_t isTPC = kFALSE;
289   Bool_t isMatch = kFALSE;
290
291
292   if(esdTrack->Charge()==0) return;
293   if(!esdTrack->GetTPCInnerParam()) return;
294   if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
295
296   if(!fCutsRC->AcceptTrack(esdTrack)) return;
297
298   isTPC = kTRUE;
299   
300   if( (esdTrack->GetStatus()&AliESDtrack::kITSrefit))
301     isMatch = kTRUE;
302   
303   if(isTPC){
304     Double_t vecTrackingEff[5] = { static_cast<Double_t>(isMatch),esdTrack->Phi(), esdTrack->Pt(),esdTrack->Eta(),static_cast<Double_t>(esdTrack->GetITSclusters(0)) };
305     fTrackingEffHisto->Fill(vecTrackingEff);
306   }
307 }
308
309 //_____________________________________________________________________________
310 /*void AliPerformanceMatch::ProcessTPCTRD(AliStack* , AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
311 {
312   return;
313 }*/
314
315 void AliPerformanceMatch::ProcessTPCConstrain(AliStack* /*const stack*/, AliESDEvent *const esdEvent, AliESDtrack *const esdTrack){
316   //
317   // Contrain TPC inner track to the vertex
318   // then compare to the global tracks
319   //
320   if(!esdTrack) return;
321   if(esdTrack->Charge()==0) return;
322
323   if(!esdTrack->GetTPCInnerParam()) return;
324   if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) return;
325   if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
326   if(!fCutsRC->AcceptTrack(esdTrack)) return;
327
328   Double_t x[3]; esdTrack->GetXYZ(x);
329   Double_t b[3]; AliTracker::GetBxByBz(x,b);
330   Bool_t isOK = kFALSE;
331   const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
332   
333   AliExternalTrackParam * TPCinner = (AliExternalTrackParam *)esdTrack->GetTPCInnerParam();
334   if(!TPCinner) return;
335   
336   //  isOK = TPCinner->Rotate(esdTrack->GetAlpha());
337   //isOK = TPCinner->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
338
339   AliExternalTrackParam * TPCinnerC = new AliExternalTrackParam(*TPCinner);
340   if (TPCinnerC) {
341     isOK = TPCinnerC->ConstrainToVertex(vtxESD, b);
342
343     // transform to the track reference frame 
344     isOK = TPCinnerC->Rotate(esdTrack->GetAlpha());
345     isOK = TPCinnerC->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
346   }
347   if(!isOK) return;
348
349   Double_t sigmaPhi=0,deltaPhi=0,pullPhi=0;
350   deltaPhi = TPCinnerC->GetSnp() - esdTrack->GetSnp();
351   sigmaPhi = TMath::Sqrt(esdTrack->GetSigmaSnp2()+TPCinnerC->GetSigmaSnp2());
352   if(sigmaPhi!=0)
353     pullPhi = deltaPhi/sigmaPhi;
354
355   Double_t vTPCConstrain[4] = {pullPhi,esdTrack->Phi(),esdTrack->Pt(),esdTrack->Eta()};
356   fTPCConstrain->Fill(vTPCConstrain);  
357   
358   if(TPCinnerC)
359     delete TPCinnerC;
360
361   return;
362 }
363 //_____________________________________________________________________________
364 void AliPerformanceMatch::FillHistograms(AliESDtrack *const refParam, AliESDtrack *const param, Bool_t isRec) 
365 {
366   //
367   // fill performance histograms 
368   // (TPC always as reference)
369   //
370
371   if(!refParam) return;
372   if(!param) return;
373   if(!isRec) return;
374   
375   //
376   // Deltas (dy,dz,dphi,dtheta,dpt)
377   //
378   Float_t delta[5] = {0};
379   if(param && isRec) {
380     delta[0] = param->GetY()-refParam->GetY();
381     delta[1] = param->GetZ()-refParam->GetZ();
382     delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
383     delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
384     if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
385   }
386   // 
387   // Pulls (y,z,snp,tanl,1/pt)
388   //
389   Float_t sigma[5] = {0};
390   Float_t pull[5] = {0};
391   if(param && isRec) {
392     sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());   
393     sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
394     sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
395     sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
396     sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
397     if(sigma[0]) pull[0] = delta[0] / sigma[0]; 
398     if(sigma[1]) pull[1] = delta[1] / sigma[1]; 
399     if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2]; 
400     if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3]; 
401     if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4]; 
402   }
403
404   // Fill histograms
405   Double_t vResolHisto[9] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->Phi(),refParam->Eta(),refParam->Pt(),static_cast<Double_t>(isRec)};
406   if(fabs(pull[4])<5)
407     fResolHisto->Fill(vResolHisto);
408
409   Double_t vPullHisto[9] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->Phi(),refParam->Eta(),refParam->OneOverPt(),static_cast<Double_t>(isRec)};
410   if(fabs(pull[4])<5)
411     fPullHisto->Fill(vPullHisto);
412 }
413
414 //_____________________________________________________________________________
415 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend */*const esdFriend*/, const Bool_t bUseMC, const Bool_t /*bUseESDfriend*/)
416 {
417   // Process comparison information 
418   //
419   if(!esdEvent) 
420   {
421     Error("Exec","esdEvent not available");
422     return;
423   }
424   AliHeader* header = 0;
425   AliGenEventHeader* genHeader = 0;
426   AliStack* stack = 0;
427   TArrayF vtxMC(3);
428   
429   if(bUseMC)
430   {
431     if(!mcEvent) {
432       Error("Exec","mcEvent not available");
433       return;
434     }
435     // get MC event header
436     header = mcEvent->Header();
437     if (!header) {
438       Error("Exec","Header not available");
439       return;
440     }
441     // MC particle stack
442     stack = mcEvent->Stack();
443     if (!stack) {
444       Error("Exec","Stack not available");
445       return;
446     }
447     // get MC vertex
448     genHeader = header->GenEventHeader();
449     if (!genHeader) {
450       Error("Exec","Could not retrieve genHeader from Header");
451       return;
452     }
453     genHeader->PrimaryVertex(vtxMC);
454   } 
455   
456   // use ESD friends
457   /*  if(bUseESDfriend) {
458     if(!esdFriend) {
459       Error("Exec","esdFriend not available");
460       return;
461     }
462     }*/
463
464   // trigger
465   if(!bUseMC && GetTriggerClass()) {
466     Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
467     if(!isEventTriggered) return; 
468   }
469
470   // get TPC event vertex
471   const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
472   if(vtxESD && (vtxESD->GetStatus()<=0)) return;
473
474   //  Process events
475   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
476   { 
477     AliESDtrack *track = esdEvent->GetTrack(iTrack);
478     if(!track) continue;
479
480     /*AliESDfriendTrack *friendTrack=0;
481     if(bUseESDfriend) {
482       friendTrack=esdFriend->GetTrack(iTrack);
483       if(!friendTrack) continue;
484       }*/
485
486     if(GetAnalysisMode() == 0){
487       if(!IsUseTOFBunchCrossing()){
488         ProcessTPCITS(stack,track);
489       }
490       else
491         if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
492           ProcessTPCITS(stack,track);
493         }
494     }
495     /* else if(GetAnalysisMode() == 2) ProcessTPCTRD(stack,track,friendTrack);*/
496     else if(GetAnalysisMode() == 1) {ProcessITSTPC(iTrack,esdEvent,stack,track);}
497     else if(GetAnalysisMode() == 2){
498       if(!IsUseTOFBunchCrossing()){
499         ProcessTPCConstrain(stack,esdEvent,track);
500       }
501       else
502         if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
503           ProcessTPCConstrain(stack,esdEvent,track);
504         }
505     }
506     else {
507       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
508       return;
509     }
510   }
511
512 }
513
514 //_____________________________________________________________________________
515 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
516   // Create resolution histograms
517  
518   TH1F *hisr, *hism;
519   if (!gPad) new TCanvas;
520   hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
521   if (type) return hism;
522   else 
523     return hisr;
524 }
525
526 //_____________________________________________________________________________
527 void AliPerformanceMatch::Analyse() {
528   // Analyse comparison information and store output histograms
529   // in the folder "folderMatch"
530   //
531   TString selString;
532   /*
533   TH1::AddDirectory(kFALSE);
534   TH1F *h=0;
535   TH1F *h2=0;
536   TH2F *h2D=0;
537   */
538   TObjArray *aFolderObj = new TObjArray;
539
540   // write results in the folder 
541   // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
542   // c->cd();
543
544   // char name[256];
545   // char title[256];
546
547   if(GetAnalysisMode()==1) { 
548
549   fResolHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
550   fPullHisto->GetAxis(8)->SetRangeUser(1.0,2.0);  // only reconstructed
551   for(Int_t i=0; i<5; i++) 
552     {
553       for(Int_t j=5; j<8; j++) 
554         {
555           //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
556           if(j!=6) fResolHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
557           else fResolHisto->GetAxis(6)->SetRangeUser(-1.5,1.49);
558           fResolHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
559           
560           selString = "resol";
561           AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
562           
563           
564           if(j!=6) fPullHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
565           else  fPullHisto->GetAxis(6)->SetRangeUser(-1.5,1.49); // eta window
566           fPullHisto->GetAxis(7)->SetRangeUser(0.01,10.);  // pt threshold
567           selString = "pull";
568           AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
569           
570         }
571     }
572   }
573   // 
574   // TPC efficiency wrt ITS
575   //
576   if(GetAnalysisMode()==0) { 
577     selString = "trackingeff";
578     AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
579     
580     for(Int_t i=1; i<5; i++) 
581       {
582         // all ITS standalone tracks
583         fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
584         selString = "trackingeff_all";
585         AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
586         
587       // TPC tracks which has matching with TPC
588         fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
589         selString = "trackingeff_tpc";
590         AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
591       }
592   }
593
594   //
595   //TPC constrained to vertex
596   //
597   if(GetAnalysisMode()==2) { 
598     selString = "tpc";
599     //    for(Int_t i=1; i<4; i++)
600     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 2, &selString);
601     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 3, &selString);
602     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 2, 3, &selString);
603   }
604   
605   printf("exportToFolder\n");
606   fAnalysisFolder = ExportToFolder(aFolderObj);
607   
608   // delete only TObjArray
609   if(fFolderObj) delete fFolderObj;
610   fFolderObj = aFolderObj;  
611   aFolderObj=0;  
612   
613   }
614
615
616 //_____________________________________________________________________________
617 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array) 
618 {
619   // recreate folder avery time and export objects to new one
620   //
621   AliPerformanceMatch * comp=this;
622   TFolder *folder = comp->GetAnalysisFolder();
623
624   TString name, title;
625   TFolder *newFolder = 0;
626   Int_t i = 0;
627   Int_t size = array->GetSize();
628
629   if(folder) { 
630      // get name and title from old folder
631      name = folder->GetName();  
632      title = folder->GetTitle();  
633
634          // delete old one
635      delete folder;
636
637          // create new one
638      newFolder = CreateFolder(name.Data(),title.Data());
639      newFolder->SetOwner();
640
641          // add objects to folder
642      while(i < size) {
643        newFolder->Add(array->At(i));
644        i++;
645      }
646   }
647   
648   return newFolder;
649 }
650  
651 //_____________________________________________________________________________
652 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) { 
653 // create folder for analysed histograms
654 //
655 TFolder *folder = 0;
656   folder = new TFolder(name.Data(),title.Data());
657
658   return folder;
659 }
660
661 //_____________________________________________________________________________
662 Long64_t AliPerformanceMatch::Merge(TCollection* const list) 
663 {
664   // Merge list of objects (needed by PROOF)
665
666   if (!list)
667   return 0;
668
669   if (list->IsEmpty())
670   return 1;
671   
672   Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
673
674   TIterator* iter = list->MakeIterator();
675   TObject* obj = 0;
676   TObjArray* objArrayList = 0;
677   objArrayList = new TObjArray();
678
679   // collection of generated histograms
680   Int_t count=0;
681   while((obj = iter->Next()) != 0) 
682   {
683     AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
684     if (entry == 0) continue; 
685     if (merge) {
686         if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
687         if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
688         if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }
689
690         if ((fTPCConstrain) && (entry->fTPCConstrain)) { fTPCConstrain->Add(entry->fTPCConstrain); }
691     }
692     // the analysisfolder is only merged if present
693     if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
694
695     count++;
696   }
697   if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
698   // to signal that track histos were not merged: reset
699   if (!merge) { fResolHisto->Reset(); fPullHisto->Reset(); fTrackingEffHisto->Reset(); fTPCConstrain->Reset(); }
700   // delete
701   if (objArrayList)  delete objArrayList;  objArrayList=0;
702 return count;
703 }