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