]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliPerformanceDEdx.cxx
coverity fix
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliPerformanceDEdx.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliPerformanceDEdx class. It keeps information from 
3 // comparison of reconstructed and MC particle tracks. In addtion, 
4 // it keeps selection cuts used during comparison. The comparison 
5 // information is stored in the ROOT histograms. Analysis of these 
6 // histograms can be done by using Analyse() class function. The result of 
7 // the analysis (histograms/graphs) are stored in the folder which is 
8 // a data of AliPerformanceDEdx.
9 //  
10 // Author: J.Otwinowski 04/02/2008 
11 // Changes by M.Knichel 15/10/2010
12 //------------------------------------------------------------------------------
13
14 /*
15  
16   // after running comparison task, read the file, and get component
17   gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
18   LoadMyLibs();
19   TFile f("Output.root");
20   //AliPerformanceDEdx * compObj = (AliPerformanceDEdx*)f.Get("AliPerformanceDEdx");
21   AliPerformanceDEdx * compObj = (AliPerformanceDEdx*)coutput->FindObject("AliPerformanceDEdx");
22
23   // Analyse comparison data
24   compObj->Analyse();
25
26   // the output histograms/graphs will be stored in the folder "folderDEdx" 
27   compObj->GetAnalysisFolder()->ls("*");
28
29   // user can save whole comparison object (or only folder with anlysed histograms) 
30   // in the seperate output file (e.g.)
31   TFile fout("Analysed_DEdx.root","recreate");
32   //compObj->Write(); // compObj->GetAnalysisFolder()->Write();
33   compObj->GetAnalysisFolder()->Write();
34   fout.Close();
35
36 */
37
38 #include "TDirectory.h"
39 #include "TAxis.h"
40 #include "TCanvas.h"
41 #include "TH1.h"
42 #include "TH2.h"
43 #include "TF1.h"
44 #include "TSystem.h"
45 #include "TChain.h"
46
47 #include "AliPerformanceDEdx.h"
48 #include "AliPerformanceTPC.h"
49 #include "AliTPCPerformanceSummary.h"
50 #include "AliESDEvent.h"
51 #include "AliTracker.h"
52 #include "AliMCEvent.h"
53 #include "AliESDtrack.h"
54 #include "AliStack.h"
55 #include "AliLog.h" 
56 #include "AliMCInfoCuts.h" 
57 #include "AliMathBase.h"
58 #include "AliRecInfoCuts.h" 
59 #include "AliTreeDraw.h"
60 #include "AliHeader.h"
61 #include "AliGenEventHeader.h"
62
63 using namespace std;
64
65 ClassImp(AliPerformanceDEdx)
66
67 Bool_t AliPerformanceDEdx::fgMergeTHnSparse = kFALSE;
68 Bool_t AliPerformanceDEdx::fgUseMergeTHnSparse = kFALSE;
69
70 //_____________________________________________________________________________
71 AliPerformanceDEdx::AliPerformanceDEdx():
72   AliPerformanceObject("AliPerformanceDEdx"),
73
74   // dEdx 
75   fDeDxHisto(0),
76   fFolderObj(0),
77   
78   // Cuts 
79   fCutsRC(0), 
80   fCutsMC(0),
81
82   // histogram folder 
83   fAnalysisFolder(0)
84 {
85   // default constructor        
86   Init();
87 }
88
89 //_____________________________________________________________________________
90 AliPerformanceDEdx::AliPerformanceDEdx(Char_t* name="AliPerformanceDEdx", Char_t* title="AliPerformanceDEdx",Int_t analysisMode=0, Bool_t hptGenerator=kFALSE):
91   AliPerformanceObject(name,title),
92
93   // dEdx 
94   fDeDxHisto(0),
95   fFolderObj(0),
96   
97   // Cuts 
98   fCutsRC(0), 
99   fCutsMC(0),
100
101   // histogram folder 
102   fAnalysisFolder(0)
103 {
104   // named constructor
105
106   SetAnalysisMode(analysisMode);
107   SetHptGenerator(hptGenerator);
108   Init();
109 }
110
111
112 //_____________________________________________________________________________
113 AliPerformanceDEdx::~AliPerformanceDEdx()
114 {
115   // destructor
116   if(fDeDxHisto)  delete fDeDxHisto; fDeDxHisto=0; 
117   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
118 }
119
120 //_____________________________________________________________________________
121 void AliPerformanceDEdx::Init()
122 {
123   // Init histograms
124   
125   // TPC dEdx
126   // set pt bins
127   Int_t nPBins = 50;
128   Double_t pMin = 1.e-2, pMax = 10.;
129
130   Double_t *binsP = 0;
131   if (IsHptGenerator())  { 
132     nPBins = 100; pMax = 100.;
133     binsP = CreateLogAxis(nPBins,pMin,pMax);
134   } else {
135     binsP = CreateLogAxis(nPBins,pMin,pMax);
136   }
137
138
139   /*
140   Int_t nPBins = 31;
141     Double_t binsP[32] = {0.,0.05,0.1,0.15,0.2,0.25,0.3,0.35,0.4,0.45,0.5,0.55,0.6,0.7,0.8,0.9,1.0,1.2,1.4,1.6,1.8,2.0,2.25,2.5,2.75,3.,3.5,4.,5.,6.,8.,10.};
142     Double_t pMin = 0., pMax = 10.;
143
144     if(IsHptGenerator() == kTRUE) {
145       nPBins = 100;
146       pMin = 0.; pMax = 100.;
147     }
148    */
149
150    //dedx:alpha:y:z:snp:tgl:ncls:p
151    //dedx:phi:y:z:snp:tgl:ncls:p
152    //Int_t binsQA[8]    = {300, 50, 50,  50, 50, 50, 80, nPBins};
153    //Double_t xminQA[8] = {0, -4,-20,-250, -1, -2, 0, pMin};
154    //Double_t xmaxQA[8] = {300, 4, 20, 250,  1,  2, 160, pMax};
155    Int_t binsQA[9]    = {300, 144, 50,  50, 50, 50, 80, nPBins, 160};
156    Double_t xminQA[9] = {0, -TMath::Pi(),-20,-250, -1, -2, 0, pMin, 0.};
157    Double_t xmaxQA[9] = {300, TMath::Pi(), 20, 250,  1,  2, 160, pMax ,160.};
158
159    //fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:alpha:y:z:snp:tgl:ncls:momentum",8,binsQA,xminQA,xmaxQA);
160    fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:phi:y:z:snp:tgl:ncls:momentum:TPCSignalN",9,binsQA,xminQA,xmaxQA);
161    fDeDxHisto->SetBinEdges(7,binsP);
162
163    fDeDxHisto->GetAxis(0)->SetTitle("dedx (a.u.)");
164    fDeDxHisto->GetAxis(1)->SetTitle("#phi (rad)");
165    fDeDxHisto->GetAxis(2)->SetTitle("y (cm)");
166    fDeDxHisto->GetAxis(3)->SetTitle("z (cm)");
167    fDeDxHisto->GetAxis(4)->SetTitle("sin#phi");
168    fDeDxHisto->GetAxis(5)->SetTitle("tan#lambda");
169    fDeDxHisto->GetAxis(6)->SetTitle("ncls");
170    fDeDxHisto->GetAxis(7)->SetTitle("p (GeV/c)");
171    fDeDxHisto->GetAxis(8)->SetTitle("number of cls used for dEdx");
172    //fDeDxHisto->Sumw2();
173
174    // Init cuts
175    if(!fCutsMC) {
176      AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
177    }
178    if(!fCutsRC) {
179      AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
180    }
181
182    // init folder
183    fAnalysisFolder = CreateFolder("folderDEdx","Analysis de/dx Folder");
184    
185    // save merge status in object
186    fMergeTHnSparseObj = fgMergeTHnSparse;
187
188 }
189
190 //_____________________________________________________________________________
191 void AliPerformanceDEdx::ProcessTPC(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
192 {
193   // Fill dE/dx  comparison information
194   AliDebug(AliLog::kWarning, "Warning: Not implemented");
195 }
196
197 //_____________________________________________________________________________
198 void AliPerformanceDEdx::ProcessInnerTPC(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent* const esdEvent)
199 {
200  //
201  // Fill TPC track information at inner TPC wall
202  // 
203   if(!esdEvent) return;
204   if(!esdTrack) return;
205
206   if( IsUseTrackVertex() ) 
207   { 
208     // Relate TPC inner params to prim. vertex
209     const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
210     Double_t x[3]; esdTrack->GetXYZ(x);
211     Double_t b[3]; AliTracker::GetBxByBz(x,b);
212     Bool_t isOK = kFALSE;
213     if(fabs(b[2])>0.000001)
214       isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
215     //    Bool_t isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
216     if(!isOK) return;
217
218     /*
219       // JMT -- recaluclate DCA for HLT if not present
220       if ( dca[0] == 0. && dca[1] == 0. ) {
221         track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
222       }
223     */
224   }
225
226   // get external param. at inner TPC wall
227   const AliExternalTrackParam *innerParam =  esdTrack->GetInnerParam();
228   if(!innerParam) return;
229
230   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
231   esdTrack->GetImpactParametersTPC(dca,cov);
232
233   if((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return; // TPC refit
234
235   //
236   // select primaries
237   //
238   Double_t dcaToVertex = -1;
239   if( fCutsRC->GetDCAToVertex2D() ) 
240   {
241       dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY()                    + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
242   }
243   if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
244   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
245   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
246
247   Float_t dedx = esdTrack->GetTPCsignal();
248   Int_t ncls = esdTrack->GetTPCNcls();
249   Int_t TPCSignalN = esdTrack->GetTPCsignalN();
250
251
252   Double_t pt = innerParam->Pt();
253   Double_t lam = TMath::ATan2(innerParam->Pz(),innerParam->Pt());
254   Double_t p = pt/TMath::Cos(lam);
255   //Double_t alpha = innerParam->GetAlpha();
256   Double_t phi = TMath::ATan2(innerParam->Py(),innerParam->Px());
257   //if(phi<0.) phi += 2.*TMath::Phi();
258   Double_t y = innerParam->GetY();
259   Double_t z = innerParam->GetZ();
260   Double_t snp = innerParam->GetSnp();
261   Double_t tgl = innerParam->GetTgl();
262
263   //Double_t vDeDxHisto[8] = {dedx,alpha,y,z,snp,tgl,ncls,p};
264   Double_t vDeDxHisto[9] = {dedx,phi,y,z,snp,tgl,ncls,p,TPCSignalN};
265   fDeDxHisto->Fill(vDeDxHisto); 
266
267   if(!stack) return;
268 }
269
270 //_____________________________________________________________________________
271 void AliPerformanceDEdx::ProcessTPCITS(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
272 {
273   // Fill dE/dx  comparison information
274   
275    AliDebug(AliLog::kWarning, "Warning: Not implemented");
276 }
277
278 //_____________________________________________________________________________
279 void AliPerformanceDEdx::ProcessConstrained(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
280 {
281   // Fill dE/dx  comparison information
282   
283    AliDebug(AliLog::kWarning, "Warning: Not implemented");
284 }
285
286 //_____________________________________________________________________________
287 Long64_t AliPerformanceDEdx::Merge(TCollection* const list) 
288 {
289   // Merge list of objects (needed by PROOF)
290
291   if (!list)
292   return 0;
293
294   if (list->IsEmpty())
295   return 1;
296   
297   Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
298
299   TIterator* iter = list->MakeIterator();
300   TObject* obj = 0;
301   TObjArray* objArrayList = 0;
302   objArrayList = new TObjArray();
303
304   // collection of generated histograms
305   Int_t count=0;
306   while((obj = iter->Next()) != 0) 
307   {
308     AliPerformanceDEdx* entry = dynamic_cast<AliPerformanceDEdx*>(obj);
309     if (entry == 0) continue; 
310     if (merge) {
311         if ((fDeDxHisto) && (entry->fDeDxHisto)) { fDeDxHisto->Add(entry->fDeDxHisto); }        
312     }
313     // the analysisfolder is only merged if present
314     if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
315
316     count++;
317   }
318   if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
319   // to signal that track histos were not merged: reset
320   if (!merge) { fDeDxHisto->Reset(); }
321   // delete
322   if (objArrayList)  delete objArrayList;  objArrayList=0;  
323
324 return count;
325 }
326
327 //_____________________________________________________________________________
328 void AliPerformanceDEdx::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
329 {
330   // Process comparison information 
331   //
332   if(!esdEvent) 
333   {
334       AliDebug(AliLog::kError, "esdEvent not available");
335       return;
336   }
337   AliHeader* header = 0;
338   AliGenEventHeader* genHeader = 0;
339   AliStack* stack = 0;
340   TArrayF vtxMC(3);
341   
342   if(bUseMC)
343   {
344     if(!mcEvent) {
345       AliDebug(AliLog::kError, "mcEvent not available");
346       return;
347     }
348
349     // get MC event header
350     header = mcEvent->Header();
351     if (!header) {
352       AliDebug(AliLog::kError, "Header not available");
353       return;
354     }
355     // MC particle stack
356     stack = mcEvent->Stack();
357     if (!stack) {
358       AliDebug(AliLog::kError, "Stack not available");
359       return;
360     }
361
362     // get MC vertex
363     genHeader = header->GenEventHeader();
364     if (!genHeader) {
365       AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");
366       return;
367     }
368     genHeader->PrimaryVertex(vtxMC);
369
370   } // end bUseMC
371
372   // use ESD friends
373   if(bUseESDfriend) {
374     if(!esdFriend) {
375       AliDebug(AliLog::kError, "esdFriend not available");
376       return;
377     }
378   }
379
380   // trigger
381   if(!bUseMC && GetTriggerClass()) {
382     Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
383     if(!isEventTriggered) return; 
384   }
385
386   // get event vertex
387   const AliESDVertex *vtxESD = NULL;
388   if( IsUseTrackVertex() ) 
389   { 
390     // track vertex
391     vtxESD = esdEvent->GetPrimaryVertexTracks();
392   }
393   else {
394     // TPC track vertex
395     vtxESD = esdEvent->GetPrimaryVertexTPC();
396   }
397   if(vtxESD && (vtxESD->GetStatus()<=0)) return;
398
399   //  Process events
400   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
401   { 
402     AliESDtrack *track = esdEvent->GetTrack(iTrack);
403     if(!track) continue;
404
405     if(GetAnalysisMode() == 0) ProcessTPC(stack,track);
406     else if(GetAnalysisMode() == 1) ProcessTPCITS(stack,track);
407     else if(GetAnalysisMode() == 2) ProcessConstrained(stack,track);
408     else if(GetAnalysisMode() == 3) ProcessInnerTPC(stack,track,esdEvent);
409     else {
410       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
411       return;
412     }
413   }
414 }
415
416 //_____________________________________________________________________________
417 void AliPerformanceDEdx::Analyse()
418 {
419   // Analyze comparison information and store output histograms
420   // in the folder "folderDEdx"
421   //
422   TH1::AddDirectory(kFALSE);
423   TH1::SetDefaultSumw2(kFALSE);
424   TH1F *h1D=0;
425   TH2F *h2D=0;
426   TObjArray *aFolderObj = new TObjArray;
427   TString selString;
428
429   char name[256];
430   char title[256];
431
432   for(Int_t i=1; i<9; i++) { 
433     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, i);
434   }
435
436   // resolution histograms for mips
437   //dedx:phi:y:z:snp:tgl:ncls:p
438   fDeDxHisto->GetAxis(2)->SetRangeUser(-15.,14.999);
439   fDeDxHisto->GetAxis(3)->SetRangeUser(-120.,119.999);
440   fDeDxHisto->GetAxis(4)->SetRangeUser(-0.4, 0.399);
441   fDeDxHisto->GetAxis(5)->SetRangeUser(-0.9,0.89);
442   fDeDxHisto->GetAxis(6)->SetRangeUser(60.,160.);
443   fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.499);
444   fDeDxHisto->GetAxis(8)->SetRangeUser(60.,160.);
445   
446   
447   selString = "mipsres";
448   AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
449
450   //
451   TObjArray *arr[8] = {0};
452   TF1 *f1[8] = {0};
453   
454   for(Int_t i=0; i<8; i++) 
455   { 
456     arr[i] = new TObjArray;
457     f1[i] = new TF1("gaus","gaus");
458     //printf("i %d \n",i);
459
460     h2D = (TH2F*)fDeDxHisto->Projection(0,i+1);
461
462     f1[i]->SetRange(40,60); // should be pion peak
463     h2D->FitSlicesY(f1[i],0,-1,10,"QNR",arr[i]); // gaus fit of pion peak
464
465     h1D = (TH1F*)arr[i]->At(1);
466     snprintf(name,256,"mean_dedx_mips_vs_%d",i);
467     h1D->SetName(name);
468     snprintf(title,256,"%s vs %s","mean_dedx_mips (a.u.)",fDeDxHisto->GetAxis(i)->GetTitle());
469     h1D->SetTitle(title);
470     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
471     h1D->GetYaxis()->SetTitle("mean_dedx_mips (a.u.)");
472     //h1D->SetMinimum(40);
473     //h1D->SetMaximum(60);
474
475     aFolderObj->Add(h1D);
476
477     h1D = (TH1F*)arr[i]->At(2);
478     snprintf(name,256,"res_dedx_mips_vs_%d",i);
479     h1D->SetName(name);
480     snprintf(title,256,"%s vs %s","res_dedx_mips (a.u)",fDeDxHisto->GetAxis(i)->GetTitle());
481     h1D->SetTitle(title);
482     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
483     h1D->GetYaxis()->SetTitle("res_dedx_mips (a.u.)");
484     //h1D->SetMinimum(0);
485     //h1D->SetMaximum(6);
486
487     aFolderObj->Add(h1D);
488   }
489
490     // select MIPs (version from AliTPCPerfomanceSummary)
491     fDeDxHisto->GetAxis(0)->SetRangeUser(35,60);
492     fDeDxHisto->GetAxis(2)->SetRangeUser(-20,20);
493     fDeDxHisto->GetAxis(3)->SetRangeUser(-250,250);
494     fDeDxHisto->GetAxis(4)->SetRangeUser(-1, 1);
495     fDeDxHisto->GetAxis(5)->SetRangeUser(-1,1);
496     fDeDxHisto->GetAxis(6)->SetRangeUser(80,160);
497     fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.55);
498     fDeDxHisto->GetAxis(8)->SetRangeUser(80,160);
499
500     selString = "mips";
501     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
502     
503     selString = "mips_C";
504     fDeDxHisto->GetAxis(5)->SetRangeUser(-3,0);
505     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);
506     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
507     
508     selString = "mips_A";
509     fDeDxHisto->GetAxis(5)->SetRangeUser(0,3);
510     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);    
511     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
512     
513     //restore cuts
514     fDeDxHisto->GetAxis(0)->SetRangeUser(0,300);
515     fDeDxHisto->GetAxis(2)->SetRangeUser(-20,20);
516     fDeDxHisto->GetAxis(3)->SetRangeUser(-250,250);
517     fDeDxHisto->GetAxis(4)->SetRangeUser(-1, 1);
518     fDeDxHisto->GetAxis(5)->SetRangeUser(-3,3);
519     fDeDxHisto->GetAxis(6)->SetRangeUser(0,160);
520     fDeDxHisto->GetAxis(7)->SetRangeUser(1e-2,10);    
521     fDeDxHisto->GetAxis(8)->SetRangeUser(0,160);
522
523     printf("exportToFolder\n");
524     // export objects to analysis folder
525     fAnalysisFolder = ExportToFolder(aFolderObj);
526     if (fFolderObj) delete fFolderObj;
527     fFolderObj = aFolderObj;
528     aFolderObj=0;
529
530
531   for(Int_t i=0;i<7;i++) { 
532     if(f1[i]) delete f1[i]; f1[i]=0;
533   }
534
535 }
536
537 //_____________________________________________________________________________
538 TFolder* AliPerformanceDEdx::ExportToFolder(TObjArray * array) 
539 {
540   // recreate folder avery time and export objects to new one
541   //
542   AliPerformanceDEdx * comp=this;
543   TFolder *folder = comp->GetAnalysisFolder();
544
545   TString name, title;
546   TFolder *newFolder = 0;
547   Int_t i = 0;
548   Int_t size = array->GetSize();
549
550   if(folder) { 
551      // get name and title from old folder
552      name = folder->GetName();  
553      title = folder->GetTitle();  
554
555          // delete old one
556      delete folder;
557
558          // create new one
559      newFolder = CreateFolder(name.Data(),title.Data());
560      newFolder->SetOwner();
561
562          // add objects to folder
563      while(i < size) {
564            newFolder->Add(array->At(i));
565            i++;
566          }
567   }
568
569 return newFolder;
570 }
571
572
573 //_____________________________________________________________________________
574 TFolder* AliPerformanceDEdx::CreateFolder(TString name,TString title) { 
575 // create folder for analysed histograms
576 TFolder *folder = 0;
577   folder = new TFolder(name.Data(),title.Data());
578
579   return folder;
580 }
581
582 //_____________________________________________________________________________
583 TTree* AliPerformanceDEdx::CreateSummary()
584 {
585     // implementaion removed, switched back to use AliPerformanceSummary (now called in AliPerformanceTask)
586     return 0;
587 }
588