]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceMatch.cxx
Merge remote-tracking branch 'remotes/origin/master' into TPCdev
[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; continue; } 
263     
264      if(!fCutsRC->AcceptTrack(tpcTrack2)) { delete tpcTrack2; continue; }
265     // check matching
266     if (TMath::Abs(esdTrack->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; continue; }
267     if (TMath::Abs(esdTrack->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; continue; }
268     if (TMath::Abs(esdTrack->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; continue; }
269     
270     hasMatch=kTRUE;
271     break;
272   }
273   
274   FillHistograms(tpcTrack2,esdTrack,hasMatch);     
275   /*  if(tpcTrack2) { 
276     delete tpcTrack2;
277    
278     }*/
279 }
280
281 //_____________________________________________________________________________
282 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack)
283 {
284   //
285   // Match TPC and ITS min-bias tracks
286   // at radius between detectors
287   //
288   if(!esdTrack) return;
289   //   if(!esdFriendTrack) return;
290    
291   Bool_t isTPC = kFALSE;
292   Bool_t isMatch = kFALSE;
293
294
295   if(esdTrack->Charge()==0) return;
296   if(!esdTrack->GetTPCInnerParam()) return;
297   if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
298
299   if(!fCutsRC->AcceptTrack(esdTrack)) return;
300
301   isTPC = kTRUE;
302   
303   if( (esdTrack->GetStatus()&AliESDtrack::kITSrefit))
304     isMatch = kTRUE;
305   
306   if(isTPC){
307     Double_t vecTrackingEff[5] = { isMatch,esdTrack->Phi(), esdTrack->Pt(),esdTrack->Eta(),esdTrack->GetITSclusters(0) };
308     fTrackingEffHisto->Fill(vecTrackingEff);
309   }
310 }
311
312 //_____________________________________________________________________________
313 /*void AliPerformanceMatch::ProcessTPCTRD(AliStack* , AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
314 {
315   return;
316 }*/
317
318 void AliPerformanceMatch::ProcessTPCConstrain(AliStack* /*const stack*/, AliESDEvent *const esdEvent, AliESDtrack *const esdTrack){
319   //
320   // Contrain TPC inner track to the vertex
321   // then compare to the global tracks
322   //
323   if(!esdTrack) return;
324   if(esdTrack->Charge()==0) return;
325
326   if(!esdTrack->GetTPCInnerParam()) return;
327   if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) return;
328   if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
329   if(!fCutsRC->AcceptTrack(esdTrack)) return;
330
331   Double_t x[3]; esdTrack->GetXYZ(x);
332   Double_t b[3]; AliTracker::GetBxByBz(x,b);
333   Bool_t isOK = kFALSE;
334   const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
335   
336   AliExternalTrackParam * TPCinner = (AliExternalTrackParam *)esdTrack->GetTPCInnerParam();
337   if(!TPCinner) return;
338   
339   //  isOK = TPCinner->Rotate(esdTrack->GetAlpha());
340   //isOK = TPCinner->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
341
342   AliExternalTrackParam * TPCinnerC = new AliExternalTrackParam(*TPCinner);
343   if (TPCinnerC) {
344     isOK = TPCinnerC->ConstrainToVertex(vtxESD, b);
345
346     // transform to the track reference frame 
347     isOK = TPCinnerC->Rotate(esdTrack->GetAlpha());
348     isOK = TPCinnerC->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
349   }
350   if(!isOK) return;
351
352   Double_t sigmaPhi=0,deltaPhi=0,pullPhi=0;
353   deltaPhi = TPCinnerC->GetSnp() - esdTrack->GetSnp();
354   sigmaPhi = TMath::Sqrt(esdTrack->GetSigmaSnp2()+TPCinnerC->GetSigmaSnp2());
355   if(sigmaPhi!=0)
356     pullPhi = deltaPhi/sigmaPhi;
357
358   Double_t vTPCConstrain[4] = {pullPhi,esdTrack->Phi(),esdTrack->Pt(),esdTrack->Eta()};
359   fTPCConstrain->Fill(vTPCConstrain);  
360   
361   if(TPCinnerC)
362     delete TPCinnerC;
363
364   return;
365 }
366 //_____________________________________________________________________________
367 void AliPerformanceMatch::FillHistograms(AliESDtrack *const refParam, AliESDtrack *const param, Bool_t isRec) 
368 {
369   //
370   // fill performance histograms 
371   // (TPC always as reference)
372   //
373
374   if(!refParam) return;
375   if(!param) return;
376   if(!isRec) return;
377   
378   //
379   // Deltas (dy,dz,dphi,dtheta,dpt)
380   //
381   Float_t delta[5] = {0};
382   if(param && isRec) {
383     delta[0] = param->GetY()-refParam->GetY();
384     delta[1] = param->GetZ()-refParam->GetZ();
385     delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
386     delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
387     if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
388   }
389   // 
390   // Pulls (y,z,snp,tanl,1/pt)
391   //
392   Float_t sigma[5] = {0};
393   Float_t pull[5] = {0};
394   if(param && isRec) {
395     sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());   
396     sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
397     sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
398     sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
399     sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
400     if(sigma[0]) pull[0] = delta[0] / sigma[0]; 
401     if(sigma[1]) pull[1] = delta[1] / sigma[1]; 
402     if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2]; 
403     if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3]; 
404     if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4]; 
405   }
406
407   // Fill histograms
408   Double_t vResolHisto[9] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->Phi(),refParam->Eta(),refParam->Pt(),isRec};
409   if(fabs(pull[4])<5)
410     fResolHisto->Fill(vResolHisto);
411
412   Double_t vPullHisto[9] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->Phi(),refParam->Eta(),refParam->OneOverPt(),isRec};
413   if(fabs(pull[4])<5)
414     fPullHisto->Fill(vPullHisto);
415 }
416
417 //_____________________________________________________________________________
418 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend */*const esdFriend*/, const Bool_t bUseMC, const Bool_t /*bUseESDfriend*/)
419 {
420   // Process comparison information 
421   //
422   if(!esdEvent) 
423   {
424     Error("Exec","esdEvent not available");
425     return;
426   }
427   AliHeader* header = 0;
428   AliGenEventHeader* genHeader = 0;
429   AliStack* stack = 0;
430   TArrayF vtxMC(3);
431   
432   if(bUseMC)
433   {
434     if(!mcEvent) {
435       Error("Exec","mcEvent not available");
436       return;
437     }
438     // get MC event header
439     header = mcEvent->Header();
440     if (!header) {
441       Error("Exec","Header not available");
442       return;
443     }
444     // MC particle stack
445     stack = mcEvent->Stack();
446     if (!stack) {
447       Error("Exec","Stack not available");
448       return;
449     }
450     // get MC vertex
451     genHeader = header->GenEventHeader();
452     if (!genHeader) {
453       Error("Exec","Could not retrieve genHeader from Header");
454       return;
455     }
456     genHeader->PrimaryVertex(vtxMC);
457   } 
458   
459   // use ESD friends
460   /*  if(bUseESDfriend) {
461     if(!esdFriend) {
462       Error("Exec","esdFriend not available");
463       return;
464     }
465     }*/
466
467   // trigger
468   if(!bUseMC && GetTriggerClass()) {
469     Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
470     if(!isEventTriggered) return; 
471   }
472
473   // get TPC event vertex
474   const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
475   if(vtxESD && (vtxESD->GetStatus()<=0)) return;
476
477   //  Process events
478   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
479   { 
480     AliESDtrack *track = esdEvent->GetTrack(iTrack);
481     if(!track) continue;
482
483     /*AliESDfriendTrack *friendTrack=0;
484     if(bUseESDfriend) {
485       friendTrack=esdFriend->GetTrack(iTrack);
486       if(!friendTrack) continue;
487       }*/
488
489     if(GetAnalysisMode() == 0){
490       if(!IsUseTOFBunchCrossing()){
491         ProcessTPCITS(stack,track);
492       }
493       else
494         if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
495           ProcessTPCITS(stack,track);
496         }
497     }
498     /* else if(GetAnalysisMode() == 2) ProcessTPCTRD(stack,track,friendTrack);*/
499     else if(GetAnalysisMode() == 1) {ProcessITSTPC(iTrack,esdEvent,stack,track);}
500     else if(GetAnalysisMode() == 2){
501       if(!IsUseTOFBunchCrossing()){
502         ProcessTPCConstrain(stack,esdEvent,track);
503       }
504       else
505         if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
506           ProcessTPCConstrain(stack,esdEvent,track);
507         }
508     }
509     else {
510       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
511       return;
512     }
513   }
514
515 }
516
517 //_____________________________________________________________________________
518 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
519   // Create resolution histograms
520  
521   TH1F *hisr, *hism;
522   if (!gPad) new TCanvas;
523   hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
524   if (type) return hism;
525   else 
526     return hisr;
527 }
528
529 //_____________________________________________________________________________
530 void AliPerformanceMatch::Analyse() {
531   // Analyse comparison information and store output histograms
532   // in the folder "folderMatch"
533   //
534   TString selString;
535   /*
536   TH1::AddDirectory(kFALSE);
537   TH1F *h=0;
538   TH1F *h2=0;
539   TH2F *h2D=0;
540   */
541   TObjArray *aFolderObj = new TObjArray;
542
543   // write results in the folder 
544   // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
545   // c->cd();
546
547   // char name[256];
548   // char title[256];
549
550   if(GetAnalysisMode()==1) { 
551
552   fResolHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
553   fPullHisto->GetAxis(8)->SetRangeUser(1.0,2.0);  // only reconstructed
554   for(Int_t i=0; i<5; i++) 
555     {
556       for(Int_t j=5; j<8; j++) 
557         {
558           //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
559           if(j!=6) fResolHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
560           else fResolHisto->GetAxis(6)->SetRangeUser(-1.5,1.49);
561           fResolHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
562           
563           selString = "resol";
564           AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
565           
566           
567           if(j!=6) fPullHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
568           else  fPullHisto->GetAxis(6)->SetRangeUser(-1.5,1.49); // eta window
569           fPullHisto->GetAxis(7)->SetRangeUser(0.01,10.);  // pt threshold
570           selString = "pull";
571           AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
572           
573         }
574     }
575   }
576   // 
577   // TPC efficiency wrt ITS
578   //
579   if(GetAnalysisMode()==0) { 
580     selString = "trackingeff";
581     AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
582     
583     for(Int_t i=1; i<5; i++) 
584       {
585         // all ITS standalone tracks
586         fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
587         selString = "trackingeff_all";
588         AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
589         
590       // TPC tracks which has matching with TPC
591         fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
592         selString = "trackingeff_tpc";
593         AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
594       }
595   }
596
597   //
598   //TPC constrained to vertex
599   //
600   if(GetAnalysisMode()==2) { 
601     selString = "tpc";
602     //    for(Int_t i=1; i<4; i++)
603     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 2, &selString);
604     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 3, &selString);
605     AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 2, 3, &selString);
606   }
607   
608   printf("exportToFolder\n");
609   fAnalysisFolder = ExportToFolder(aFolderObj);
610   
611   // delete only TObjArray
612   if(fFolderObj) delete fFolderObj;
613   fFolderObj = aFolderObj;  
614   aFolderObj=0;  
615   
616   }
617
618
619 //_____________________________________________________________________________
620 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array) 
621 {
622   // recreate folder avery time and export objects to new one
623   //
624   AliPerformanceMatch * comp=this;
625   TFolder *folder = comp->GetAnalysisFolder();
626
627   TString name, title;
628   TFolder *newFolder = 0;
629   Int_t i = 0;
630   Int_t size = array->GetSize();
631
632   if(folder) { 
633      // get name and title from old folder
634      name = folder->GetName();  
635      title = folder->GetTitle();  
636
637          // delete old one
638      delete folder;
639
640          // create new one
641      newFolder = CreateFolder(name.Data(),title.Data());
642      newFolder->SetOwner();
643
644          // add objects to folder
645      while(i < size) {
646        newFolder->Add(array->At(i));
647        i++;
648      }
649   }
650   
651   return newFolder;
652 }
653  
654 //_____________________________________________________________________________
655 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) { 
656 // create folder for analysed histograms
657 //
658 TFolder *folder = 0;
659   folder = new TFolder(name.Data(),title.Data());
660
661   return folder;
662 }
663
664 //_____________________________________________________________________________
665 Long64_t AliPerformanceMatch::Merge(TCollection* const list) 
666 {
667   // Merge list of objects (needed by PROOF)
668
669   if (!list)
670   return 0;
671
672   if (list->IsEmpty())
673   return 1;
674   
675   Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
676
677   TIterator* iter = list->MakeIterator();
678   TObject* obj = 0;
679   TObjArray* objArrayList = 0;
680   objArrayList = new TObjArray();
681
682   // collection of generated histograms
683   Int_t count=0;
684   while((obj = iter->Next()) != 0) 
685   {
686     AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
687     if (entry == 0) continue; 
688     if (merge) {
689         if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
690         if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
691         if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }
692
693         if ((fTPCConstrain) && (entry->fTPCConstrain)) { fTPCConstrain->Add(entry->fTPCConstrain); }
694     }
695     // the analysisfolder is only merged if present
696     if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
697
698     count++;
699   }
700   if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
701   // to signal that track histos were not merged: reset
702   if (!merge) { fResolHisto->Reset(); fPullHisto->Reset(); fTrackingEffHisto->Reset(); fTPCConstrain->Reset(); }
703   // delete
704   if (objArrayList)  delete objArrayList;  objArrayList=0;
705 return count;
706 }