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