]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceDEdx.cxx
Removing some (but not all!) annoying warnings
[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    //dedx:alpha:y:z:snp:tgl:ncls:p
150    //dedx:phi:y:z:snp:tgl:ncls:p
151    //Int_t binsQA[8]    = {300, 50, 50,  50, 50, 50, 80, nPBins};
152    //Double_t xminQA[8] = {0, -4,-20,-250, -1, -2, 0, pMin};
153    //Double_t xmaxQA[8] = {300, 4, 20, 250,  1,  2, 160, pMax};
154   Int_t binsQA[10]    = {300, 144, 50,  50, 50, 50, 80, nPBins, 80, 80};
155   Double_t xminQA[10] = {0, -TMath::Pi(),-20,-250, -1, -2, 0, pMin, 0., 0.};
156   Double_t xmaxQA[10] = {300, TMath::Pi(), 20, 250,  1,  2, 160, pMax ,160., 160.};
157
158    //fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:alpha:y:z:snp:tgl:ncls:momentum",8,binsQA,xminQA,xmaxQA);
159    fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:phi:y:z:snp:tgl:ncls:momentum:TPCSignalN:CrossedRows",10,binsQA,xminQA,xmaxQA);
160    fDeDxHisto->SetBinEdges(7,binsP);
161
162    fDeDxHisto->GetAxis(0)->SetTitle("dedx (a.u.)");
163    fDeDxHisto->GetAxis(1)->SetTitle("#phi (rad)");
164    fDeDxHisto->GetAxis(2)->SetTitle("y (cm)");
165    fDeDxHisto->GetAxis(3)->SetTitle("z (cm)");
166    fDeDxHisto->GetAxis(4)->SetTitle("sin#phi");
167    fDeDxHisto->GetAxis(5)->SetTitle("tan#lambda");
168    fDeDxHisto->GetAxis(6)->SetTitle("ncls");
169    fDeDxHisto->GetAxis(7)->SetTitle("p (GeV/c)");
170    fDeDxHisto->GetAxis(8)->SetTitle("number of cls used for dEdx");
171    fDeDxHisto->GetAxis(9)->SetTitle("Crossed rows in the TPC");
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   Float_t nCrossedRows = esdTrack->GetTPCClusterInfo(2,1);
251
252
253   Double_t pt = innerParam->Pt();
254   Double_t lam = TMath::ATan2(innerParam->Pz(),innerParam->Pt());
255   Double_t p = pt/TMath::Cos(lam);
256   //Double_t alpha = innerParam->GetAlpha();
257   Double_t phi = TMath::ATan2(innerParam->Py(),innerParam->Px());
258   //if(phi<0.) phi += 2.*TMath::Phi();
259   Double_t y = innerParam->GetY();
260   Double_t z = innerParam->GetZ();
261   Double_t snp = innerParam->GetSnp();
262   Double_t tgl = innerParam->GetTgl();
263
264   //Double_t vDeDxHisto[8] = {dedx,alpha,y,z,snp,tgl,ncls,p};
265   Double_t vDeDxHisto[10] = {dedx,phi,y,z,snp,tgl,ncls,p,TPCSignalN,nCrossedRows};
266   fDeDxHisto->Fill(vDeDxHisto); 
267
268   if(!stack) return;
269 }
270
271 //_____________________________________________________________________________
272 void AliPerformanceDEdx::ProcessTPCITS(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
273 {
274   // Fill dE/dx  comparison information
275   
276    AliDebug(AliLog::kWarning, "Warning: Not implemented");
277 }
278
279 //_____________________________________________________________________________
280 void AliPerformanceDEdx::ProcessConstrained(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
281 {
282   // Fill dE/dx  comparison information
283   
284    AliDebug(AliLog::kWarning, "Warning: Not implemented");
285 }
286
287 //_____________________________________________________________________________
288 Long64_t AliPerformanceDEdx::Merge(TCollection* const list) 
289 {
290   // Merge list of objects (needed by PROOF)
291
292   if (!list)
293   return 0;
294
295   if (list->IsEmpty())
296   return 1;
297   
298   Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
299
300   TIterator* iter = list->MakeIterator();
301   TObject* obj = 0;
302   TObjArray* objArrayList = 0;
303   objArrayList = new TObjArray();
304
305   // collection of generated histograms
306   Int_t count=0;
307   while((obj = iter->Next()) != 0) 
308   {
309     AliPerformanceDEdx* entry = dynamic_cast<AliPerformanceDEdx*>(obj);
310     if (entry == 0) continue; 
311     if (merge) {
312         if ((fDeDxHisto) && (entry->fDeDxHisto)) { fDeDxHisto->Add(entry->fDeDxHisto); }        
313     }
314     // the analysisfolder is only merged if present
315     if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
316
317     count++;
318   }
319   if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
320   // to signal that track histos were not merged: reset
321   if (!merge) { fDeDxHisto->Reset(); }
322   // delete
323   if (objArrayList)  delete objArrayList;  objArrayList=0;  
324
325 return count;
326 }
327
328 //_____________________________________________________________________________
329 void AliPerformanceDEdx::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
330 {
331   // Process comparison information 
332   //
333   if(!esdEvent) 
334   {
335       AliDebug(AliLog::kError, "esdEvent not available");
336       return;
337   }
338   AliHeader* header = 0;
339   AliGenEventHeader* genHeader = 0;
340   AliStack* stack = 0;
341   TArrayF vtxMC(3);
342   
343   if(bUseMC)
344   {
345     if(!mcEvent) {
346       AliDebug(AliLog::kError, "mcEvent not available");
347       return;
348     }
349
350     // get MC event header
351     header = mcEvent->Header();
352     if (!header) {
353       AliDebug(AliLog::kError, "Header not available");
354       return;
355     }
356     // MC particle stack
357     stack = mcEvent->Stack();
358     if (!stack) {
359       AliDebug(AliLog::kError, "Stack not available");
360       return;
361     }
362
363     // get MC vertex
364     genHeader = header->GenEventHeader();
365     if (!genHeader) {
366       AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");
367       return;
368     }
369     genHeader->PrimaryVertex(vtxMC);
370
371   } // end bUseMC
372
373   // use ESD friends
374   if(bUseESDfriend) {
375     if(!esdFriend) {
376       AliDebug(AliLog::kError, "esdFriend not available");
377       return;
378     }
379   }
380
381   // trigger
382   if(!bUseMC && GetTriggerClass()) {
383     Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
384     if(!isEventTriggered) return; 
385   }
386
387   // get event vertex
388   const AliESDVertex *vtxESD = NULL;
389   if( IsUseTrackVertex() ) 
390   { 
391     // track vertex
392     vtxESD = esdEvent->GetPrimaryVertexTracks();
393   }
394   else {
395     // TPC track vertex
396     vtxESD = esdEvent->GetPrimaryVertexTPC();
397   }
398   if(vtxESD && (vtxESD->GetStatus()<=0)) return;
399
400   //  Process events
401   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
402   { 
403     AliESDtrack *track = esdEvent->GetTrack(iTrack);
404     if(!track) continue;
405
406     if(GetAnalysisMode() == 0) ProcessTPC(stack,track);
407     else if(GetAnalysisMode() == 1) ProcessTPCITS(stack,track);
408     else if(GetAnalysisMode() == 2) ProcessConstrained(stack,track);
409     else if(GetAnalysisMode() == 3) ProcessInnerTPC(stack,track,esdEvent);
410     else {
411       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
412       return;
413     }
414   }
415 }
416
417 //_____________________________________________________________________________
418 void AliPerformanceDEdx::Analyse()
419 {
420   // Analyze comparison information and store output histograms
421   // in the folder "folderDEdx"
422   //
423   TH1::AddDirectory(kFALSE);
424   TH1::SetDefaultSumw2(kFALSE);
425   TH1F *h1D=0;
426   TH2F *h2D=0;
427   TObjArray *aFolderObj = new TObjArray;
428   TString selString;
429
430   char name[256];
431   char title[256];
432
433   for(Int_t i=1; i<10; i++) { 
434     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, i);
435   }
436
437     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 6, 7);
438     AddProjection(aFolderObj, "dedx", fDeDxHisto, 7, 8, 9);
439     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 8, 9);
440     AddProjection(aFolderObj, "dedx", fDeDxHisto, 6, 8, 9);
441
442   // resolution histograms for mips
443   //dedx:phi:y:z:snp:tgl:ncls:p
444   fDeDxHisto->GetAxis(2)->SetRangeUser(-15.,14.999);
445   fDeDxHisto->GetAxis(3)->SetRangeUser(-120.,119.999);
446   fDeDxHisto->GetAxis(4)->SetRangeUser(-0.4, 0.399);
447   fDeDxHisto->GetAxis(5)->SetRangeUser(-0.9,0.89);
448   fDeDxHisto->GetAxis(6)->SetRangeUser(60.,160.);
449   fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.499);
450   fDeDxHisto->GetAxis(8)->SetRangeUser(60.,160.);
451   
452   
453   selString = "mipsres";
454   AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
455
456   //
457   TObjArray *arr[9] = {0};
458   TF1 *f1[9] = {0};
459   
460   for(Int_t i=0; i<9; i++) 
461   { 
462     arr[i] = new TObjArray;
463     f1[i] = new TF1("gaus","gaus");
464     //printf("i %d \n",i);
465
466     h2D = (TH2F*)fDeDxHisto->Projection(0,i+1);
467
468     f1[i]->SetRange(40,60); // should be pion peak
469     h2D->FitSlicesY(f1[i],0,-1,10,"QNR",arr[i]); // gaus fit of pion peak
470
471     h1D = (TH1F*)arr[i]->At(1);
472     snprintf(name,256,"mean_dedx_mips_vs_%d",i);
473     h1D->SetName(name);
474     snprintf(title,256,"%s vs %s","mean_dedx_mips (a.u.)",fDeDxHisto->GetAxis(i)->GetTitle());
475     h1D->SetTitle(title);
476     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
477     h1D->GetYaxis()->SetTitle("mean_dedx_mips (a.u.)");
478     //h1D->SetMinimum(40);
479     //h1D->SetMaximum(60);
480
481     aFolderObj->Add(h1D);
482
483     h1D = (TH1F*)arr[i]->At(2);
484     snprintf(name,256,"res_dedx_mips_vs_%d",i);
485     h1D->SetName(name);
486     snprintf(title,256,"%s vs %s","res_dedx_mips (a.u)",fDeDxHisto->GetAxis(i)->GetTitle());
487     h1D->SetTitle(title);
488     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
489     h1D->GetYaxis()->SetTitle("res_dedx_mips (a.u.)");
490     //h1D->SetMinimum(0);
491     //h1D->SetMaximum(6);
492
493     aFolderObj->Add(h1D);
494   }
495
496     // select MIPs (version from AliTPCPerfomanceSummary)
497     fDeDxHisto->GetAxis(0)->SetRangeUser(35,60);
498     fDeDxHisto->GetAxis(2)->SetRangeUser(-20,20);
499     fDeDxHisto->GetAxis(3)->SetRangeUser(-250,250);
500     fDeDxHisto->GetAxis(4)->SetRangeUser(-1, 1);
501     fDeDxHisto->GetAxis(5)->SetRangeUser(-1,1);
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     selString = "mips";
508     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, &selString);
509     
510     selString = "mips_C";
511     fDeDxHisto->GetAxis(5)->SetRangeUser(-3,0);
512     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);
513     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
514     
515     selString = "mips_A";
516     fDeDxHisto->GetAxis(5)->SetRangeUser(0,3);
517     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 5, &selString);    
518     AddProjection(aFolderObj, "dedx", fDeDxHisto, 0, 1, &selString);
519     
520     //restore cuts
521     fDeDxHisto->GetAxis(0)->SetRangeUser(0,300);
522     fDeDxHisto->GetAxis(2)->SetRangeUser(-20,20);
523     fDeDxHisto->GetAxis(3)->SetRangeUser(-250,250);
524     fDeDxHisto->GetAxis(4)->SetRangeUser(-1, 1);
525     fDeDxHisto->GetAxis(5)->SetRangeUser(-3,3);
526     fDeDxHisto->GetAxis(6)->SetRangeUser(0,160);
527     fDeDxHisto->GetAxis(7)->SetRangeUser(1e-2,10);    
528     fDeDxHisto->GetAxis(8)->SetRangeUser(0,160);
529     fDeDxHisto->GetAxis(9)->SetRangeUser(0.,1.);
530
531
532
533
534     printf("exportToFolder\n");
535     // export objects to analysis folder
536     fAnalysisFolder = ExportToFolder(aFolderObj);
537     if (fFolderObj) delete fFolderObj;
538     fFolderObj = aFolderObj;
539     aFolderObj=0;
540
541
542   for(Int_t i=0;i<9;i++) { 
543     if(f1[i]) delete f1[i]; f1[i]=0;
544   }
545
546 }
547
548 //_____________________________________________________________________________
549 TFolder* AliPerformanceDEdx::ExportToFolder(TObjArray * array) 
550 {
551   // recreate folder avery time and export objects to new one
552   //
553   AliPerformanceDEdx * comp=this;
554   TFolder *folder = comp->GetAnalysisFolder();
555
556   TString name, title;
557   TFolder *newFolder = 0;
558   Int_t i = 0;
559   Int_t size = array->GetSize();
560
561   if(folder) { 
562      // get name and title from old folder
563      name = folder->GetName();  
564      title = folder->GetTitle();  
565
566          // delete old one
567      delete folder;
568
569          // create new one
570      newFolder = CreateFolder(name.Data(),title.Data());
571      newFolder->SetOwner();
572
573          // add objects to folder
574      while(i < size) {
575            newFolder->Add(array->At(i));
576            i++;
577          }
578   }
579
580 return newFolder;
581 }
582
583
584 //_____________________________________________________________________________
585 TFolder* AliPerformanceDEdx::CreateFolder(TString name,TString title) { 
586 // create folder for analysed histograms
587 TFolder *folder = 0;
588   folder = new TFolder(name.Data(),title.Data());
589
590   return folder;
591 }
592
593 //_____________________________________________________________________________
594 TTree* AliPerformanceDEdx::CreateSummary()
595 {
596     // implementaion removed, switched back to use AliPerformanceSummary (now called in AliPerformanceTask)
597     return 0;
598 }
599