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