]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliPerformanceMatch.cxx
move all TPC related classes
[u/mrichter/AliRoot.git] / PWG1 / 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 //------------------------------------------------------------------------------
15
16 /*
17  
18   // after running comparison task, read the file, and get component
19   gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
20   LoadMyLibs();
21
22   TFile f("Output.root");
23   AliPerformanceMatch * compObj = (AliPerformanceMatch*)coutput->FindObject("AliPerformanceMatchTPCTRD");
24  
25   // analyse comparison data
26   compObj->Analyse();
27
28   // the output histograms/graphs will be stored in the folder "folderMatch" 
29   compObj->GetAnalysisFolder()->ls("*");
30
31   // user can save whole comparison object (or only folder with anlysed histograms) 
32   // in the seperate output file (e.g.)
33   TFile fout("Analysed_Match.root","recreate");
34   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
35   fout.Close();
36
37 */
38
39 #include "TCanvas.h"
40 #include "TH1.h"
41 #include "TH2.h"
42 #include "TAxis.h"
43 #include "TF1.h"
44
45 #include "AliPerformanceMatch.h" 
46 #include "AliESDEvent.h" 
47 #include "AliESDVertex.h"
48 #include "AliESDtrack.h"
49 #include "AliESDfriendTrack.h"
50 #include "AliESDfriend.h"
51 #include "AliLog.h" 
52 #include "AliMCEvent.h" 
53 #include "AliMCParticle.h" 
54 #include "AliHeader.h" 
55 #include "AliGenEventHeader.h" 
56 #include "AliStack.h" 
57 #include "AliMCInfoCuts.h" 
58 #include "AliRecInfoCuts.h" 
59 #include "AliTracker.h" 
60 #include "AliTRDtrackV1.h" 
61 #include "AliTreeDraw.h" 
62
63 using namespace std;
64
65 ClassImp(AliPerformanceMatch)
66
67 //_____________________________________________________________________________
68 AliPerformanceMatch::AliPerformanceMatch():
69   AliPerformanceObject("AliPerformanceMatch"),
70   fResolHisto(0),
71   fPullHisto(0),
72
73   // Cuts 
74   fCutsRC(0),  
75   fCutsMC(0),  
76
77   // histogram folder 
78   fAnalysisFolder(0)
79 {
80   Init();
81 }
82
83 //_____________________________________________________________________________
84 AliPerformanceMatch::AliPerformanceMatch(Char_t* name="AliPerformanceMatch", Char_t* title="AliPerformanceMatch",Int_t analysisMode=0,Bool_t hptGenerator=kFALSE):
85   AliPerformanceObject(name,title),
86   fResolHisto(0),
87   fPullHisto(0),
88
89   // Cuts 
90   fCutsRC(0),  
91   fCutsMC(0),  
92
93   // histogram folder 
94   fAnalysisFolder(0)
95 {
96   // named constructor  
97   // 
98   SetAnalysisMode(analysisMode);
99   SetHptGenerator(hptGenerator);
100
101   Init();
102 }
103
104 //_____________________________________________________________________________
105 AliPerformanceMatch::~AliPerformanceMatch()
106 {
107   // destructor
108    
109   if(fResolHisto) delete fResolHisto; fResolHisto=0;     
110   if(fPullHisto)  delete fPullHisto;  fPullHisto=0;     
111   
112   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
113 }
114
115 //_____________________________________________________________________________
116 void AliPerformanceMatch::Init(){
117
118   //
119   // Make performance histogrms
120   //
121
122   // set pt bins
123   Int_t nPtBins = 50;
124   Double_t ptMin = 1.e-2, ptMax = 10.;
125
126   Double_t *binsPt = 0;
127   if (IsHptGenerator())  { 
128     nPtBins = 100; ptMax = 100.;
129     binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
130   } else {
131     binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
132   }
133
134   Double_t yMin = -0.02, yMax = 0.02;
135   Double_t zMin = -12.0, zMax = 12.0;
136   if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) { 
137     yMin = -100.; yMax = 100.; 
138     zMin = -100.; zMax = 100.; 
139   }
140
141   // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
142   Int_t binsResolHisto[11]={100,100,100,100,100,25,50,90,30,nPtBins,2};
143   Double_t minResolHisto[11]={-1.,-1.,-0.03,-0.03,-0.2, yMin, zMin, 0., -1.5, ptMin,0};
144   Double_t maxResolHisto[11]={ 1., 1., 0.03, 0.03, 0.2, yMax, zMax, 2.*TMath::Pi(), 1.5, ptMax,2};
145
146   fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:y:z:phi:eta:pt:isRec",11,binsResolHisto,minResolHisto,maxResolHisto);
147   fResolHisto->SetBinEdges(9,binsPt);
148
149   fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
150   fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
151   fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
152   fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
153   fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
154   fResolHisto->GetAxis(5)->SetTitle("y_{ref} (cm)");
155   fResolHisto->GetAxis(6)->SetTitle("z_{ref} (cm)");
156   fResolHisto->GetAxis(7)->SetTitle("#phi_{ref} (rad)");
157   fResolHisto->GetAxis(8)->SetTitle("#eta_{ref}");
158   fResolHisto->GetAxis(9)->SetTitle("p_{Tref} (GeV/c)");
159   fResolHisto->GetAxis(10)->SetTitle("isReconstructed");
160   fResolHisto->Sumw2();
161
162   //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
163   Int_t binsPullHisto[11]={100,100,100,100,100,50,50,50,50,nPtBins,2};
164   Double_t minPullHisto[11]={-5.,-5.,-5.,-5.,-5., yMin, zMin,-1.,-2.0, ptMin,0};
165   Double_t maxPullHisto[11]={ 5., 5., 5., 5., 5., yMax, zMax, 1., 2.0, ptMax,2};
166   fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec",11,binsPullHisto,minPullHisto,maxPullHisto);
167
168   fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
169   fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
170   fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
171   fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
172   fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
173   fPullHisto->GetAxis(5)->SetTitle("y_{ref} (cm)");
174   fPullHisto->GetAxis(6)->SetTitle("z_{ref} (cm)");
175   fPullHisto->GetAxis(7)->SetTitle("sin#phi_{ref}");
176   fPullHisto->GetAxis(8)->SetTitle("tan#lambda_{ref}");
177   fPullHisto->GetAxis(9)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
178   fPullHisto->GetAxis(10)->SetTitle("isReconstructed");
179   fPullHisto->Sumw2();
180
181   // Init cuts 
182   if(!fCutsMC) 
183     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
184   if(!fCutsRC) 
185     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
186
187   // init folder
188   fAnalysisFolder = CreateFolder("folderRes","Analysis Resolution Folder");
189 }
190
191 //_____________________________________________________________________________
192 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
193 {
194   //
195   // Match TPC and ITS min-bias tracks
196   // at radius between detectors
197   //
198   if(!esdTrack) return;
199   if(!esdFriendTrack) return;
200    
201   //
202   // Propagate tracks to the radius between TPC-ITS
203   // using B-field and material budget
204   //
205   Double_t radius = fCutsRC->GetTPCITSMatchingRadius();
206   Double_t mass = esdTrack->GetMass();
207   Double_t step=1.0; // cm
208
209
210   //
211   // Propagate TPCinner (reference detector)
212   //
213   Bool_t isTPCOK=kFALSE;
214   AliExternalTrackParam *innerTPC=NULL;
215
216   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
217   esdTrack->GetImpactParametersTPC(dca,cov);
218
219   if( (esdTrack->GetNcls(1)>fCutsRC->GetMinNClustersTPC()) && 
220       (TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ()) && 
221       (esdTrack->GetTPCInnerParam()) &&
222       (innerTPC=new AliExternalTrackParam(*(esdTrack->GetTPCInnerParam())))) 
223   {
224      isTPCOK = AliTracker::PropagateTrackToBxByBz(innerTPC,radius,mass,step,kTRUE);
225   }
226   if(!isTPCOK) return;
227
228   //
229   // Propagate ITSouter
230   //
231   Bool_t isITSOK=kFALSE;
232   AliExternalTrackParam *outerITS=NULL;
233
234   if( (esdTrack->GetNcls(0)>fCutsRC->GetMinNClustersITS()) &&
235       (esdFriendTrack->GetITSOut()) &&
236       (outerITS=new AliExternalTrackParam(*(esdFriendTrack->GetITSOut())))) 
237   {
238     isITSOK = AliTracker::PropagateTrackToBxByBz(outerITS,radius,mass,step,kTRUE);
239   }
240
241   //
242   // Fill histograms (TPC reference detector)
243   //
244   if(isTPCOK)
245     FillHistograms(innerTPC,outerITS,isITSOK);
246
247   if(outerITS) delete outerITS;  
248   if(innerTPC) delete innerTPC;  
249 }
250
251 //_____________________________________________________________________________
252 void AliPerformanceMatch::ProcessTPCTRD(AliStack* /*const stack*/, AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
253 {
254   //
255   // Match TPC and TRD min-bias tracks
256   // at radius between detectors. TPC is the reference detector.
257   //
258   if(!esdTrack) return;
259   if(!esdFriendTrack) return;
260   
261   //
262   // Propagate tracks to the radius between TPC-TRD
263   // using B-field and material budget
264   //
265   Double_t radius = fCutsRC->GetTPCTRDMatchingRadius();
266   Double_t mass = esdTrack->GetMass();
267   Double_t step=1.0; // cm
268
269   //
270   // Propagate TPCouter (reference detector)
271   //
272   Bool_t isTPCOK=kFALSE;
273   AliExternalTrackParam *outerTPC=NULL;
274
275   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
276   esdTrack->GetImpactParametersTPC(dca,cov);
277
278   if( (esdTrack->GetNcls(1)>fCutsRC->GetMinNClustersTPC()) && 
279       (TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ()) && 
280       (esdFriendTrack->GetTPCOut()) &&
281       (outerTPC=new AliExternalTrackParam(*(esdFriendTrack->GetTPCOut())))) 
282   {
283      isTPCOK = AliTracker::PropagateTrackToBxByBz(outerTPC,radius,mass,step,kTRUE);
284   }
285   if(!isTPCOK) return;
286
287   //
288   // Propagate TRDinner
289   //
290   Bool_t isTRDOK = kFALSE;
291   AliExternalTrackParam *innerTRD=NULL;
292
293   // get TRD track
294   AliTRDtrackV1 *trdTrack=NULL; //esdFriendTrack = fESDfriend->GetTrack(itrk);
295   TObject *calObject=NULL;
296   Int_t icalib = 0;
297   while((calObject = esdFriendTrack->GetCalibObject(icalib++))) {
298     if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
299     if(!(trdTrack = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
300   }
301
302   //if( (esdTrack->GetNcls(2)>fCutsRC->GetMinNClustersTRD()) &&
303   if( (trdTrack) &&
304       (trdTrack->GetNumberOfTracklets()>fCutsRC->GetMinNTrackletsTRD()) &&
305       (trdTrack->GetTracklet(0)) &&
306       (esdFriendTrack->GetTRDIn()) &&
307       (innerTRD = new AliExternalTrackParam(*(esdFriendTrack->GetTRDIn())))) 
308   {
309     isTRDOK = AliTracker::PropagateTrackToBxByBz(innerTRD,radius,mass,step,kTRUE);
310   }
311
312   //
313   // Fill histograms (TPC reference detector)
314   //
315   if(isTPCOK)
316     FillHistograms(outerTPC,innerTRD,isTRDOK);
317
318   if(outerTPC) delete outerTPC;  
319   if(innerTRD) delete innerTRD;  
320 }
321
322 //_____________________________________________________________________________
323 void AliPerformanceMatch::FillHistograms(AliExternalTrackParam *const refParam, AliExternalTrackParam *const param, Bool_t isRec) 
324 {
325   //
326   // fill performance histograms 
327   // (TPC always as reference)
328   //
329   if(!refParam) return;
330
331   //
332   // Deltas (dy,dz,dphi,dtheta,dpt)
333   //
334   Float_t delta[5] = {0};
335   if(param && isRec) {
336     delta[0] = param->GetY()-refParam->GetY();
337     delta[1] = param->GetZ()-refParam->GetZ();
338     delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
339     delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
340     if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
341   }
342   // 
343   // Pulls (y,z,snp,tanl,1/pt)
344   //
345   Float_t sigma[5] = {0};
346   Float_t pull[5] = {0};
347   if(param && isRec) {
348     sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());   
349     sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
350     sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
351     sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
352     sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
353     if(sigma[0]) pull[0] = delta[0] / sigma[0]; 
354     if(sigma[1]) pull[1] = delta[1] / sigma[1]; 
355     if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2]; 
356     if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3]; 
357     if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4]; 
358   }
359
360   // Fill histograms
361   Double_t vResolHisto[11] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->GetY(),refParam->GetZ(),refParam->Phi(),refParam->Eta(),refParam->Pt(),isRec};
362   fResolHisto->Fill(vResolHisto);
363
364   Double_t vPullHisto[11] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->GetY(),refParam->GetZ(),refParam->GetSnp(),refParam->GetTgl(),refParam->OneOverPt(),isRec};
365   fPullHisto->Fill(vPullHisto);
366 }
367
368 //_____________________________________________________________________________
369 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
370 {
371   // Process comparison information 
372   //
373   if(!esdEvent) 
374   {
375     Error("Exec","esdEvent not available");
376     return;
377   }
378   AliHeader* header = 0;
379   AliGenEventHeader* genHeader = 0;
380   AliStack* stack = 0;
381   TArrayF vtxMC(3);
382   
383   if(bUseMC)
384   {
385     if(!mcEvent) {
386       Error("Exec","mcEvent not available");
387       return;
388     }
389     // get MC event header
390     header = mcEvent->Header();
391     if (!header) {
392       Error("Exec","Header not available");
393       return;
394     }
395     // MC particle stack
396     stack = mcEvent->Stack();
397     if (!stack) {
398       Error("Exec","Stack not available");
399       return;
400     }
401     // get MC vertex
402     genHeader = header->GenEventHeader();
403     if (!genHeader) {
404       Error("Exec","Could not retrieve genHeader from Header");
405       return;
406     }
407     genHeader->PrimaryVertex(vtxMC);
408   } 
409   
410   // use ESD friends
411   if(bUseESDfriend) {
412     if(!esdFriend) {
413       Error("Exec","esdFriend not available");
414       return;
415     }
416   }
417
418   //  Process events
419   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
420   { 
421     AliESDtrack *track = esdEvent->GetTrack(iTrack);
422     if(!track) continue;
423
424     AliESDfriendTrack *friendTrack=0;
425     if(bUseESDfriend) {
426       friendTrack=esdFriend->GetTrack(iTrack);
427       if(!friendTrack) continue;
428     }
429
430     if(GetAnalysisMode() == 0) ProcessTPCITS(stack,track,friendTrack);
431     else if(GetAnalysisMode() == 1) ProcessTPCTRD(stack,track,friendTrack);
432     else {
433       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
434       return;
435     }
436   }
437 }
438
439 //_____________________________________________________________________________
440 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
441   // Create resolution histograms
442  
443   TH1F *hisr, *hism;
444   if (!gPad) new TCanvas;
445   hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
446   if (type) return hism;
447   else 
448     return hisr;
449 }
450
451 //_____________________________________________________________________________
452 void AliPerformanceMatch::Analyse() {
453   // Analyse comparison information and store output histograms
454   // in the folder "folderMatch"
455   //
456   TH1::AddDirectory(kFALSE);
457   TH1F *h=0;
458   TH1F *h2=0;
459   TH2F *h2D=0;
460   TObjArray *aFolderObj = new TObjArray;
461
462   // write results in the folder 
463   TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
464   c->cd();
465
466   char name[256];
467   char title[256];
468
469   fResolHisto->GetAxis(10)->SetRangeUser(1.0,2.0); // only reconstructed
470   fPullHisto->GetAxis(10)->SetRangeUser(1.0,2.0);  // only reconstructed
471   for(Int_t i=0; i<5; i++) 
472   {
473     for(Int_t j=5; j<10; j++) 
474     {
475       //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
476       if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(0.0,0.89); // eta window
477       else fResolHisto->GetAxis(8)->SetRangeUser(-1.5,1.49);
478       fResolHisto->GetAxis(9)->SetRangeUser(0.16,100.); // pt threshold
479
480       h2D = (TH2F*)fResolHisto->Projection(i,j);
481
482       h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
483       sprintf(name,"h_res_%d_vs_%d",i,j);
484       h->SetName(name);
485
486       h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
487       sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(resolution)");
488       h->GetYaxis()->SetTitle(title);
489       sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
490       h->SetTitle(title);
491
492       if(j==9) h->SetBit(TH1::kLogX);    
493       aFolderObj->Add(h);
494
495       h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
496       //h = (TH1F*)arr->At(1);
497       sprintf(name,"h_mean_res_%d_vs_%d",i,j);
498       h->SetName(name);
499
500       h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
501       sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(mean)");
502       h->GetYaxis()->SetTitle(title);
503
504       sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
505       h->SetTitle(title);
506
507       if(j==9) h->SetBit(TH1::kLogX);    
508       aFolderObj->Add(h);
509
510       //
511       //if(j!=8) fPullHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
512       if(j!=8) fPullHisto->GetAxis(8)->SetRangeUser(0.0,0.89); // eta window
513       else  fPullHisto->GetAxis(8)->SetRangeUser(-1.5,1.49); // eta window
514       fPullHisto->GetAxis(9)->SetRangeUser(0.16,100.);  // pt threshold
515
516       h2D = (TH2F*)fPullHisto->Projection(i,j);
517
518       h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
519       sprintf(name,"h_pull_%d_vs_%d",i,j);
520       h->SetName(name);
521
522       h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
523       sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(resolution)");
524       h->GetYaxis()->SetTitle(title);
525       sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
526       h->SetTitle(title);
527
528       //if(j==9) h->SetBit(TH1::kLogX);    
529       aFolderObj->Add(h);
530
531       h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
532       sprintf(name,"h_mean_pull_%d_vs_%d",i,j);
533       h->SetName(name);
534
535       h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
536       sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(mean)");
537       h->GetYaxis()->SetTitle(title);
538       sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
539       h->SetTitle(title);
540
541       //if(j==9) h->SetBit(TH1::kLogX);    
542       aFolderObj->Add(h);
543     }
544   }
545
546   //
547   // Efficiency plots
548   //
549   for(Int_t i=5; i<10; i++) 
550   {
551     if(i!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
552     else fResolHisto->GetAxis(8)->SetRangeUser(-1.5,1.49);
553     fResolHisto->GetAxis(9)->SetRangeUser(0.16,100.); // pt threshold
554
555     fResolHisto->GetAxis(10)->SetRange(1,fResolHisto->GetAxis(10)->GetNbins()); // all 
556     h = (TH1F*)fResolHisto->Projection(i);
557
558     fResolHisto->GetAxis(10)->SetRangeUser(1.0,2.0); // only reconstructed
559     h2 = (TH1F*)fResolHisto->Projection(i);
560
561     TH1F* h2c = (TH1F*)h2->Clone();
562     h2c->Divide(h2c,h,1,1,"B");
563  
564     sprintf(name,"h_eff_%d",i);
565     h2c->SetName(name);
566
567     h2c->GetXaxis()->SetTitle(h2c->GetXaxis()->GetTitle());
568     h2c->GetYaxis()->SetTitle("efficiency");
569     h2c->SetTitle("matching effciency");
570
571     aFolderObj->Add(h2c);
572   }
573   // export objects to analysis folder
574   fAnalysisFolder = ExportToFolder(aFolderObj);
575
576   // delete only TObjArray
577   if(aFolderObj) delete aFolderObj;
578 }
579
580 //_____________________________________________________________________________
581 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array) 
582 {
583   // recreate folder avery time and export objects to new one
584   //
585   AliPerformanceMatch * comp=this;
586   TFolder *folder = comp->GetAnalysisFolder();
587
588   TString name, title;
589   TFolder *newFolder = 0;
590   Int_t i = 0;
591   Int_t size = array->GetSize();
592
593   if(folder) { 
594      // get name and title from old folder
595      name = folder->GetName();  
596      title = folder->GetTitle();  
597
598          // delete old one
599      delete folder;
600
601          // create new one
602      newFolder = CreateFolder(name.Data(),title.Data());
603      newFolder->SetOwner();
604
605          // add objects to folder
606      while(i < size) {
607            newFolder->Add(array->At(i));
608            i++;
609          }
610   }
611
612 return newFolder;
613 }
614
615 //_____________________________________________________________________________
616 Long64_t AliPerformanceMatch::Merge(TCollection* const list) 
617 {
618   // Merge list of objects (needed by PROOF)
619
620   if (!list)
621   return 0;
622
623   if (list->IsEmpty())
624   return 1;
625
626   TIterator* iter = list->MakeIterator();
627   TObject* obj = 0;
628
629   // collection of generated histograms
630   Int_t count=0;
631   while((obj = iter->Next()) != 0) 
632   {
633   AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
634   if (entry == 0) continue; 
635
636   fResolHisto->Add(entry->fResolHisto);
637   fPullHisto->Add(entry->fPullHisto);
638
639   count++;
640   }
641
642 return count;
643 }
644
645 //_____________________________________________________________________________
646 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) { 
647 // create folder for analysed histograms
648 //
649 TFolder *folder = 0;
650   folder = new TFolder(name.Data(),title.Data());
651
652   return folder;
653 }