]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceDEdx.cxx
- come back to number of find over findable clusters
[u/mrichter/AliRoot.git] / PWGPP / 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/PWGPP/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(const Char_t* name="AliPerformanceDEdx", const 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 p bins
127   Int_t nPBins = 50;
128   Double_t pMin = 1.e-2, pMax = 20.;
129
130   Double_t *binsP = 0;
131
132   if (IsHptGenerator())  { 
133         pMax = 100.;
134   } 
135    binsP = CreateLogAxis(nPBins,pMin,pMax);
136
137
138   /*
139   Int_t nPBins = 31;
140     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.};
141     Double_t pMin = 0., pMax = 10.;
142
143     if(IsHptGenerator() == kTRUE) {
144       nPBins = 100;
145       pMin = 0.; pMax = 100.;
146     }
147    */
148
149    //Int_t binsQA[8]    = {300, 50, 50,  50, 50, 50, 80, nPBins};
150    //Double_t xminQA[8] = {0, -4,-20,-250, -1, -2, 0, pMin};
151    //Double_t xmaxQA[8] = {300, 4, 20, 250,  1,  2, 160, pMax};
152  // signal:phi:y:z:snp:tgl:ncls:p:nclsDEdx:nclsF
153   Int_t binsQA[10]    = {300, 144, 50,  50, 50, 50, 80, nPBins, 160, 80};
154   Double_t xminQA[10] = {0, -TMath::Pi(),-20,-250, -1, -2, 0, pMin, 0., 0.};
155   Double_t xmaxQA[10] = {300, TMath::Pi(), 20, 250,  1,  2, 160, pMax ,160., 1.};
156
157    fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:phi:y:z:snp:tgl:ncls:momentum:TPCSignalN:clsF",10,binsQA,xminQA,xmaxQA);
158    fDeDxHisto->SetBinEdges(7,binsP);
159
160    fDeDxHisto->GetAxis(0)->SetTitle("dedx (a.u.)");
161    fDeDxHisto->GetAxis(1)->SetTitle("#phi (rad)");
162    fDeDxHisto->GetAxis(2)->SetTitle("y (cm)");
163    fDeDxHisto->GetAxis(3)->SetTitle("z (cm)");
164    fDeDxHisto->GetAxis(4)->SetTitle("sin#phi");
165    fDeDxHisto->GetAxis(5)->SetTitle("tan#lambda");
166    fDeDxHisto->GetAxis(6)->SetTitle("ncls");
167    fDeDxHisto->GetAxis(7)->SetTitle("p (GeV/c)");
168    fDeDxHisto->GetAxis(8)->SetTitle("number of cls used for dEdx");
169    fDeDxHisto->GetAxis(9)->SetTitle("number of cls found over findable");
170    //fDeDxHisto->Sumw2();
171
172    // Init cuts
173    if(!fCutsMC) {
174      AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
175    }
176    if(!fCutsRC) {
177      AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
178    }
179
180    // init folder
181    fAnalysisFolder = CreateFolder("folderDEdx","Analysis de/dx Folder");
182    
183    // save merge status in object
184    fMergeTHnSparseObj = fgMergeTHnSparse;
185
186 }
187
188 //_____________________________________________________________________________
189 void AliPerformanceDEdx::ProcessTPC(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
190 {
191   // Fill dE/dx  comparison information
192   AliDebug(AliLog::kWarning, "Warning: Not implemented");
193 }
194
195 //_____________________________________________________________________________
196 void AliPerformanceDEdx::ProcessInnerTPC(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent* const esdEvent)
197 {
198  //
199  // Fill TPC track information at inner TPC wall
200  // 
201   if(!esdEvent) return;
202   if(!esdTrack) return;
203
204   if( IsUseTrackVertex() ) 
205   { 
206     // Relate TPC inner params to prim. vertex
207     const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
208     Double_t x[3]; esdTrack->GetXYZ(x);
209     Double_t b[3]; AliTracker::GetBxByBz(x,b);
210     Bool_t isOK = kFALSE;
211     if(fabs(b[2])>0.000001)
212       isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
213     //    Bool_t isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
214     if(!isOK) return;
215
216     /*
217       // JMT -- recaluclate DCA for HLT if not present
218       if ( dca[0] == 0. && dca[1] == 0. ) {
219         track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
220       }
221     */
222   }
223
224   // get external param. at inner TPC wall
225   const AliExternalTrackParam *innerParam =  esdTrack->GetInnerParam();
226   if(!innerParam) return;
227
228   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
229   esdTrack->GetImpactParametersTPC(dca,cov);
230
231   if((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return; // TPC refit
232
233   //
234   // select primaries
235   //
236   Double_t dcaToVertex = -1;
237   if( fCutsRC->GetDCAToVertex2D() ) 
238   {
239       dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY()                    + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
240   }
241   if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
242   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
243   if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
244
245   Float_t dedx = esdTrack->GetTPCsignal();
246   Int_t ncls = esdTrack->GetTPCNcls();
247   Int_t TPCSignalN = esdTrack->GetTPCsignalN();
248   //Float_t nCrossedRows = esdTrack->GetTPCClusterInfo(2,1);
249   Float_t nClsF = esdTrack->GetTPCClusterInfo(2,0);
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[10] = {dedx,phi,y,z,snp,tgl,ncls,p,TPCSignalN,nCrossedRows};
264   Double_t vDeDxHisto[10] = {dedx,phi,y,z,snp,tgl,ncls,p,TPCSignalN,nClsF};
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<10; i++) { 
433     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, i);
434   }
435
436     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 6, 7);
437     AddProjection(aFolderObj, "dedx", fDeDxHisto, 7, 8, 9);
438     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 8, 9);
439     AddProjection(aFolderObj, "dedx", fDeDxHisto, 6, 8, 9);
440
441   // resolution histograms for mips
442   //-> signal:phi:y:z:snp:tgl:ncls:p:nclsDEdx:nclsF
443   fDeDxHisto->GetAxis(2)->SetRangeUser(-15.,14.999);
444   fDeDxHisto->GetAxis(3)->SetRangeUser(-120.,119.999);
445   fDeDxHisto->GetAxis(4)->SetRangeUser(-0.4, 0.399);
446   fDeDxHisto->GetAxis(5)->SetRangeUser(-0.9,0.89);
447   fDeDxHisto->GetAxis(6)->SetRangeUser(60.,160.);
448   fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.499);
449   fDeDxHisto->GetAxis(8)->SetRangeUser(60.,160.);
450   
451   
452   selString = "mipsres";
453   AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
454
455   //
456   TObjArray *arr[10] = {0};
457   TF1 *f1[10] = {0};
458   
459   for(Int_t i=1; i<10; i++) 
460   { 
461     arr[i] = new TObjArray;
462     f1[i] = new TF1("gaus","gaus");
463     //printf("i %d \n",i);
464
465     h2D = (TH2F*)fDeDxHisto->Projection(0,i);
466
467     f1[i]->SetRange(40,60); // should be pion peak
468     h2D->FitSlicesY(f1[i],0,-1,10,"QNR",arr[i]); // gaus fit of pion peak
469
470     h1D = (TH1F*)arr[i]->At(1);
471     snprintf(name,256,"mean_dedx_mips_vs_%d",i);
472     h1D->SetName(name);
473     snprintf(title,256,"%s vs %s","mean_dedx_mips (a.u.)",fDeDxHisto->GetAxis(i)->GetTitle());
474     h1D->SetTitle(title);
475     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
476     h1D->GetYaxis()->SetTitle("mean_dedx_mips (a.u.)");
477     //h1D->SetMinimum(40);
478     //h1D->SetMaximum(60);
479
480     aFolderObj->Add(h1D);
481
482     h1D = (TH1F*)arr[i]->At(2);
483     snprintf(name,256,"res_dedx_mips_vs_%d",i);
484     h1D->SetName(name);
485     snprintf(title,256,"%s vs %s","res_dedx_mips (a.u)",fDeDxHisto->GetAxis(i)->GetTitle());
486     h1D->SetTitle(title);
487     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
488     h1D->GetYaxis()->SetTitle("res_dedx_mips (a.u.)");
489     //h1D->SetMinimum(0);
490     //h1D->SetMaximum(6);
491
492     aFolderObj->Add(h1D);
493   }
494
495
496     // select MIPs (version from AliTPCPerfomanceSummary)
497     fDeDxHisto->GetAxis(0)->SetRangeUser(35,60);
498     fDeDxHisto->GetAxis(2)->SetRangeUser(-20,19.999);
499     fDeDxHisto->GetAxis(3)->SetRangeUser(-250,249.999);
500     fDeDxHisto->GetAxis(4)->SetRangeUser(-1, 0.99);
501     fDeDxHisto->GetAxis(5)->SetRangeUser(-1,0.99);
502     fDeDxHisto->GetAxis(6)->SetRangeUser(80,160);
503     fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.55);
504     fDeDxHisto->GetAxis(8)->SetRangeUser(80,160);
505     fDeDxHisto->GetAxis(9)->SetRangeUser(0.5,1.);
506
507
508
509     selString = "mips";
510     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
511     
512     selString = "mips_C";
513     fDeDxHisto->GetAxis(5)->SetRangeUser(-3,0);
514     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);
515     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
516     
517     selString = "mips_A";
518     fDeDxHisto->GetAxis(5)->SetRangeUser(0,3);
519     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);    
520     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
521     
522     //restore cuts
523     for (Int_t i=0; i<fDeDxHisto->GetNdimensions(); i++) {
524       fDeDxHisto->GetAxis(i)->SetRange(1,fDeDxHisto->GetAxis(i)->GetNbins());
525     }
526
527     printf("exportToFolder\n");
528     // export objects to analysis folder
529     fAnalysisFolder = ExportToFolder(aFolderObj);
530     if (fFolderObj) delete fFolderObj;
531     fFolderObj = aFolderObj;
532     aFolderObj=0;
533
534
535   for(Int_t i=0;i<10;i++) { 
536     if(f1[i]) delete f1[i]; f1[i]=0;
537   }
538
539 }
540
541 //_____________________________________________________________________________
542 TFolder* AliPerformanceDEdx::ExportToFolder(TObjArray * array) 
543 {
544   // recreate folder avery time and export objects to new one
545   //
546   AliPerformanceDEdx * comp=this;
547   TFolder *folder = comp->GetAnalysisFolder();
548
549   TString name, title;
550   TFolder *newFolder = 0;
551   Int_t i = 0;
552   Int_t size = array->GetSize();
553
554   if(folder) { 
555      // get name and title from old folder
556      name = folder->GetName();  
557      title = folder->GetTitle();  
558
559          // delete old one
560      delete folder;
561
562          // create new one
563      newFolder = CreateFolder(name.Data(),title.Data());
564      newFolder->SetOwner();
565
566          // add objects to folder
567      while(i < size) {
568            newFolder->Add(array->At(i));
569            i++;
570          }
571   }
572
573 return newFolder;
574 }
575
576
577 //_____________________________________________________________________________
578 TFolder* AliPerformanceDEdx::CreateFolder(TString name,TString title) { 
579 // create folder for analysed histograms
580 TFolder *folder = 0;
581   folder = new TFolder(name.Data(),title.Data());
582
583   return folder;
584 }
585
586 //_____________________________________________________________________________
587 TTree* AliPerformanceDEdx::CreateSummary()
588 {
589     // implementaion removed, switched back to use AliPerformanceSummary (now called in AliPerformanceTask)
590     return 0;
591 }
592