]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliComparisonDCA.cxx
1. Using the THnSparse instead of THx and TProfiles
[u/mrichter/AliRoot.git] / PWG1 / AliComparisonDCA.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliComparisonDCA 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
8 // which is a data member of AliComparisonDCA.
9 //  
10 // Author: J.Otwinowski 04/02/2008 
11 //------------------------------------------------------------------------------
12
13 /*
14  
15   // after running comparison task, read the file, and get component
16   gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
17   LoadMyLibs();
18   TFile f("Output.root");
19   //AliComparisonDCA * compObj = (AliComparisonDCA*)f.Get("AliComparisonDCA");
20   AliComparisonDCA * compObj = (AliComparisonDCA*)cOutput->FindObject("AliComparisonDCA");
21
22   // Analyse comparison data
23   compObj->Analyse();
24
25   // the output histograms/graphs will be stored in the folder "folderDCA" 
26   compObj->GetAnalysisFolder()->ls("*");
27  
28   // user can save whole comparison object (or only folder with anlysed histograms) 
29   // in the seperate output file (e.g.)
30   TFile fout("Analysed_DCA.root","recreate");
31   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
32   fout.Close();
33
34 */
35
36 #include "TGraph2D.h"
37 #include "TCanvas.h"
38 #include "TGraph.h"
39
40 #include "AliTracker.h"   
41 #include "AliESDEvent.h"   
42 #include "AliRecInfoCuts.h" 
43 #include "AliMCInfoCuts.h" 
44 #include "AliLog.h" 
45 #include "AliESDVertex.h" 
46 #include "AliMathBase.h"
47
48 #include "AliMCInfo.h" 
49 #include "AliESDRecInfo.h" 
50 #include "AliComparisonDCA.h" 
51
52 using namespace std;
53
54 ClassImp(AliComparisonDCA)
55
56 //_____________________________________________________________________________
57 AliComparisonDCA::AliComparisonDCA():
58   AliComparisonObject("AliComparisonDCA"),
59
60   // DCA histograms
61   fDCAHisto(0),
62   /*
63   fD0TanSPtTPCITS(0),
64   fD1TanSPtTPCITS(0),
65   fD0TanSPt(0),
66   fD1TanSPt(0),
67   fD0TanSPtTPC(0),
68   fD1TanSPtTPC(0),
69   */
70
71   // Cuts 
72   fCutsRC(0), 
73   fCutsMC(0),  
74
75   // histogram folder 
76   fAnalysisFolder(0)
77 {
78   // default constructor        
79 }
80
81 //_____________________________________________________________________________
82 AliComparisonDCA::AliComparisonDCA(Char_t* name="AliComparisonDCA", Char_t* title="AliComparisonDCA",Int_t analysisMode=0, Bool_t hptGenerator=kFALSE):
83   AliComparisonObject(name,title),
84
85   // DCA histograms
86   fDCAHisto(0),
87   /*
88   fD0TanSPtTPCITS(0),
89   fD1TanSPtTPCITS(0),
90   fD0TanSPt(0),
91   fD1TanSPt(0),
92   fD0TanSPtTPC(0),
93   fD1TanSPtTPC(0),
94   */
95
96   // Cuts 
97   fCutsRC(0), 
98   fCutsMC(0),  
99
100   // histogram folder 
101   fAnalysisFolder(0)
102 {
103   // named constructor   
104
105   SetAnalysisMode(analysisMode);
106   SetHptGenerator(hptGenerator);
107   Init();
108 }
109
110
111 //_____________________________________________________________________________
112 AliComparisonDCA::~AliComparisonDCA()
113 {
114   // destructor
115   if(fDCAHisto)  delete fDCAHisto; fDCAHisto=0; 
116   /*
117   if(fD0TanSPtTPCITS) delete fD0TanSPtTPCITS; fD0TanSPtTPCITS=0;
118   if(fD1TanSPtTPCITS) delete fD1TanSPtTPCITS; fD1TanSPtTPCITS=0;
119   if(fD0TanSPt) delete fD0TanSPt; fD0TanSPt=0;
120   if(fD1TanSPt) delete fD1TanSPt; fD1TanSPt=0;
121   if(fD0TanSPtTPC) delete fD0TanSPtTPC; fD0TanSPtTPC=0;
122   if(fD1TanSPtTPC) delete fD1TanSPtTPC; fD1TanSPtTPC=0;
123   */
124   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
125
126 }
127
128 //_____________________________________________________________________________
129 void AliComparisonDCA::Init()
130 {
131   // DCA histograms
132
133  Int_t nPBins = 31;
134     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.};
135     Double_t pMin = 0., pMax = 10.;
136
137     if(IsHptGenerator() == kTRUE) {
138       nPBins = 100;
139       pMin = 0.; pMax = 100.;
140     }
141
142    //dca_r, dca_z, eta, pt
143    Int_t binsQA[4]    = {100,100,20,nPBins};
144    Double_t xminQA[4] = {-10.,-10.,-1., pMin};
145    Double_t xmaxQA[4] = {10.,10.,1., pMax};
146
147    fDCAHisto = new THnSparseF("fDCAHisto","dca_r:dca_z:eta:pt",4,binsQA,xminQA,xmaxQA);
148    if(!IsHptGenerator()) fDCAHisto->SetBinEdges(3,binsP);
149
150    fDCAHisto->GetAxis(0)->SetTitle("dca_r (cm)");
151    fDCAHisto->GetAxis(1)->SetTitle("dca_z (cm)");
152    fDCAHisto->GetAxis(2)->SetTitle("eta");
153    fDCAHisto->GetAxis(3)->SetTitle("pt (GeV/c)");
154    fDCAHisto->Sumw2();
155         
156   /*    
157   fD0TanSPtTPCITS = new TH3F("DCAyTanSPtTPCITS","DCAyTanSPt",40,-2,2, 10,0.3,3, 100,-1,1);
158   fD0TanSPtTPCITS->SetXTitle("tan(#theta)");
159   fD0TanSPtTPCITS->SetYTitle("#sqrt{p_{t}(GeV)}");
160   fD0TanSPtTPCITS->SetZTitle("DCA_{xy}");
161
162   fD1TanSPtTPCITS = new TH3F("DCAzTanSPtTPCITS","DCAzTanSPt",40,-2,2, 10,0.3,3, 100,-1,1);
163   fD1TanSPtTPCITS->SetXTitle("tan(#theta)");
164   fD1TanSPtTPCITS->SetYTitle("#sqrt(p_{t}(GeV))");
165   fD1TanSPtTPCITS->SetZTitle("DCA_{z}");
166
167   fD0TanSPt = new TH3F("DCAyTanSPt","DCAyTanSPt",40,-2,2, 10,0.3,3, 100,-1,1);
168   fD0TanSPt->SetXTitle("tan(#theta)");
169   fD0TanSPt->SetYTitle("#sqrt{p_{t}(GeV)}");
170   fD0TanSPt->SetZTitle("DCA_{xy}");
171
172   fD1TanSPt = new TH3F("DCAzTanSPt","DCAzTanSPt",40,-2,2, 10,0.3,3, 100, -1,1);
173   fD1TanSPt->SetXTitle("tan(#theta)");
174   fD1TanSPt->SetYTitle("#sqrt{p_{t}(GeV)}");
175   fD1TanSPt->SetZTitle("DCA_{z}");
176
177   fD0TanSPtTPC = new TH3F("DCAyTanSPtTPC","DCAyTanSPt",40,-2,2, 10,0.3,3, 100,-1,1);
178   fD0TanSPtTPC->SetXTitle("tan(#theta)");
179   fD0TanSPtTPC->SetYTitle("#sqrt{p_{t}(GeV)}");
180   fD0TanSPtTPC->SetZTitle("DCA_{xy}");
181
182   fD1TanSPtTPC = new TH3F("DCAzTanSPtTPC","DCAzTanSPt",40,-2,2, 10,0.3,3, 100, -1,1);
183   fD1TanSPtTPC->SetXTitle("tan(#theta)");
184   fD1TanSPtTPC->SetYTitle("#sqrt{p_{t}(GeV)}");
185   fD1TanSPtTPC->SetZTitle("DCA_{z}");
186   */
187
188   // init cuts
189   if(!fCutsMC) 
190     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
191   if(!fCutsRC) 
192     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
193  
194   // init folder
195   fAnalysisFolder = CreateFolder("folderDCA","Analysis DCA Folder");
196 }
197
198 //_____________________________________________________________________________
199 void AliComparisonDCA::ProcessTPC(AliMCInfo* const infoMC, AliESDRecInfo * const infoRC)
200 {
201   // Fill DCA comparison information
202   AliExternalTrackParam *track = 0;
203   Double_t field      = AliTracker::GetBz(); // nominal Bz field [kG]
204   Double_t kMaxD      = 123456.0; // max distance
205
206   Double_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
207
208   Float_t mcpt = infoMC->GetParticle().Pt();
209   Float_t mceta = infoMC->GetParticle().Eta();
210   //Float_t spt = TMath::Sqrt(mcpt);
211
212   // distance to Prim. vertex 
213   const Double_t* dv = infoMC->GetVDist(); 
214
215   Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
216
217   // Check selection cuts
218   if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return; 
219   if (!isPrim) return;
220   if (infoRC->GetStatus(1)!=3) return;
221   if (!infoRC->GetESDtrack()) return;  
222   if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
223   //if (!infoRC->GetESDtrack()->GetConstrainedParam()) return;
224
225   // calculate and set prim. vertex
226   AliESDVertex vertexMC;
227   vertexMC.SetXv( infoMC->GetParticle().Vx() - dv[0] );
228   vertexMC.SetYv( infoMC->GetParticle().Vy() - dv[1] );
229   vertexMC.SetZv( infoMC->GetParticle().Vz() - dv[2] );
230
231   // calculate track parameters at vertex
232   if (infoRC->GetESDtrack()->GetTPCInnerParam())
233   {
234     if ((track = new AliExternalTrackParam(*infoRC->GetESDtrack()->GetTPCInnerParam())) != 0 )
235     {
236       Bool_t bDCAStatus = track->PropagateToDCA(&vertexMC,field,kMaxD,dca,cov);
237
238       if(bDCAStatus) {
239          Double_t vDCAHisto[4]={dca[0],dca[1],mceta,mcpt};
240          fDCAHisto->Fill(vDCAHisto);
241       }
242     delete track;
243     }
244   }
245 }
246
247 //_____________________________________________________________________________
248 void AliComparisonDCA::ProcessTPCITS(AliMCInfo* const infoMC, AliESDRecInfo * const infoRC)
249 {
250   // Fill DCA comparison information
251   Int_t clusterITS[200];
252   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
253
254   Float_t mcpt = infoMC->GetParticle().Pt();
255   Float_t mceta = infoMC->GetParticle().Eta();
256   //Float_t spt = TMath::Sqrt(mcpt);
257
258   // distance to Prim. vertex 
259   const Double_t* dv = infoMC->GetVDist(); 
260   Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
261
262   // Check selection cuts
263   if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return; 
264   if (!isPrim) return;
265   if (infoRC->GetStatus(1)!=3) return;
266   if (!infoRC->GetESDtrack()) return;  
267   if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
268   //if (!infoRC->GetESDtrack()->GetConstrainedParam()) return;
269
270   infoRC->GetESDtrack()->GetImpactParameters(dca,cov);
271
272   // ITS + TPC
273   if(infoRC->GetESDtrack()->GetITSclusters(clusterITS)>fCutsRC->GetMinNClustersITS())
274   {
275     Double_t vDCAHisto[4]={dca[0],dca[1],mceta,mcpt};
276     fDCAHisto->Fill(vDCAHisto);
277   }
278 }
279
280 void AliComparisonDCA::ProcessConstrained(AliMCInfo* const infoMC, AliESDRecInfo * const infoRC)
281 {
282   // Fill DCA comparison information
283   
284   AliDebug(AliLog::kWarning, "Warning: Not implemented");
285 }
286
287 //_____________________________________________________________________________
288 Long64_t AliComparisonDCA::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   TIterator* iter = list->MakeIterator();
299   TObject* obj = 0;
300
301   // collection of generated histograms
302   Int_t count=0;
303   while((obj = iter->Next()) != 0) 
304   {
305     AliComparisonDCA* entry = dynamic_cast<AliComparisonDCA*>(obj);
306     if (entry == 0) continue; 
307
308     fDCAHisto->Add(entry->fDCAHisto);
309     /*
310     fD0TanSPtTPCITS->Add(entry->fD0TanSPtTPCITS);
311     fD1TanSPtTPCITS->Add(entry->fD1TanSPtTPCITS);
312     fD0TanSPt->Add(entry->fD0TanSPt);
313     fD1TanSPt->Add(entry->fD1TanSPt);
314     fD0TanSPtTPC->Add(entry->fD0TanSPtTPC);
315     fD1TanSPtTPC->Add(entry->fD1TanSPtTPC);
316     */
317
318     count++;
319   }
320
321 return count;
322 }
323
324 //_____________________________________________________________________________
325 void AliComparisonDCA::Exec(AliMCInfo* const infoMC, AliESDRecInfo * const infoRC)
326 {
327   // Process comparison information
328   if(GetAnalysisMode() == 0) ProcessTPC(infoMC,infoRC);
329   else if(GetAnalysisMode() == 1) ProcessTPCITS(infoMC,infoRC);
330   else if(GetAnalysisMode() == 2) ProcessConstrained(infoMC,infoRC);
331   else {
332     printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
333     return;
334   }
335 }
336
337 //_____________________________________________________________________________
338 void AliComparisonDCA::Analyse()
339 {
340   //
341   // Analyse comparison information and store output histograms
342   // in the analysis folder "folderDCA" 
343   //
344   
345   TH1::AddDirectory(kFALSE);
346   TObjArray *aFolderObj = new TObjArray;
347
348   /*
349   TGraph * gr[8]= { 0,0,0,0,0,0,0,0 };
350   TGraph2D *gr2[8]= { 0,0,0,0,0,0,0,0};
351   AliComparisonDCA * comp=this;
352
353   // write results in the folder 
354   // Canvas to draw analysed histograms
355   TCanvas * c = new TCanvas("canDCA","DCA resolution");
356   c->Divide(2,4);
357   //
358   // DCA resolution
359   //
360   c->cd(1);
361   gr[0] = AliMathBase::MakeStat1D(comp->fD0TanSPtTPC,2,5);
362   gr[0]->GetXaxis()->SetTitle("Tan(#theta)");
363   gr[0]->GetYaxis()->SetTitle("#sigmaDCA_xy (cm)");
364   gr[0]->SetName("DCAXYResolTanTPC");
365   gr[0]->SetTitle("resol. DCA_xy (TPC only)");
366   gr[0]->Draw("Al*");
367
368   aFolderObj->Add(gr[0]);
369
370   c->cd(2);
371   gr[1] = AliMathBase::MakeStat1D(comp->fD1TanSPtTPC,2,5);
372   gr[1]->GetXaxis()->SetTitle("Tan(#theta)");
373   gr[1]->GetYaxis()->SetTitle("#sigmaDCA_z (cm)");
374   gr[1]->SetName("DCAZResolTanTPC");
375   gr[1]->SetTitle("resol. DCA_z (TPC only)");
376   gr[1]->Draw("Al*");
377
378   aFolderObj->Add(gr[1]);
379
380   c->cd(3);
381   gr[2] = AliMathBase::MakeStat1D(comp->fD0TanSPtTPCITS,2,5);
382   gr[2]->GetXaxis()->SetTitle("Tan(#theta)");
383   gr[2]->GetYaxis()->SetTitle("#sigmaDCA_xy (cm)");
384   gr[2]->SetName("DCAXYResolTanTPCITS");
385   gr[2]->SetTitle("resol. DCA_xy (TPC+ITS)");
386   gr[2]->Draw("Al*");
387
388   aFolderObj->Add(gr[2]);
389
390   c->cd(4);
391   gr[3] = AliMathBase::MakeStat1D(comp->fD1TanSPtTPCITS,2,5);
392   gr[3]->GetXaxis()->SetTitle("Tan(#theta)");
393   gr[3]->GetYaxis()->SetTitle("#sigmaDCA_z (cm)");
394   gr[3]->SetName("DCAZResolTanTPCITS");
395   gr[3]->SetTitle("resol. DCA_z (TPC+ITS)");
396   gr[3]->Draw("Al*");
397
398   aFolderObj->Add(gr[3]);
399
400   //
401   // DCA mean value
402   //
403   c->cd(5);
404   gr[4] = AliMathBase::MakeStat1D(comp->fD0TanSPtTPC,2,4);
405   gr[4]->GetXaxis()->SetTitle("Tan(#theta)");
406   gr[4]->GetYaxis()->SetTitle("mean DCA_xy (cm)");
407   gr[4]->SetName("DCAXYMeanTanTPC");
408   gr[4]->SetTitle("mean DCA_xy (TPC only)");
409   gr[4]->Draw("Al*");
410
411   aFolderObj->Add(gr[4]);
412
413   c->cd(6);
414   gr[5] = AliMathBase::MakeStat1D(comp->fD1TanSPtTPC,2,4);
415   gr[5]->GetXaxis()->SetTitle("Tan(#theta)");
416   gr[5]->GetYaxis()->SetTitle("mean DCA_z (cm)");
417   gr[5]->SetName("DCAZMeanTanTPC");
418   gr[5]->SetTitle("mean DCA_z (TPC only)");
419   gr[5]->Draw("Al*");
420
421   aFolderObj->Add(gr[5]);
422
423   c->cd(7);
424   gr[6] = AliMathBase::MakeStat1D(comp->fD0TanSPtTPCITS,2,4);
425   gr[6]->GetXaxis()->SetTitle("Tan(#theta)");
426   gr[6]->GetYaxis()->SetTitle("mean DCA_xy (cm)");
427   gr[6]->SetName("DCAXYMeanTanTPCITS");
428   gr[6]->SetTitle("mean DCA_xy (TPC+ITS)");
429   gr[6]->Draw("Al*");
430
431   aFolderObj->Add(gr[6]);
432
433   c->cd(8);
434   gr[7] = AliMathBase::MakeStat1D(comp->fD1TanSPtTPCITS,2,4);
435   gr[7]->GetXaxis()->SetTitle("Tan(#theta)");
436   gr[7]->GetYaxis()->SetTitle("mean DCA_z (cm)");
437   gr[7]->SetName("DCAZMeanTanTPCITS");
438   gr[7]->SetTitle("mean DCA_z (TPC+ITS)");
439   gr[7]->Draw("Al*");
440
441   aFolderObj->Add(gr[7]);
442
443   // 2D DCA resolution 
444   TCanvas * c1 = new TCanvas("canDCA1","2D DCA resolution");
445   c1->Divide(2,4);
446
447   // TPC only
448   c1->cd(1);
449   gr2[0] = AliMathBase::MakeStat2D(comp->fD0TanSPtTPC,4,2,5); 
450   gr2[0]->GetXaxis()->SetTitle("Tan(#theta)");
451   gr2[0]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
452   gr2[0]->GetZaxis()->SetTitle("#sigmaDCA_xy (cm)");
453   gr2[0]->SetName("DCAXYResolSPTTanTPC");
454   gr2[0]->SetTitle("#sigma DCA_xy (TPC only)");
455   gr2[0]->GetHistogram()->Draw("colz");
456
457   gr2[0]->GetHistogram()->SetName("DCAXYResolSPTTanTPC");
458   aFolderObj->Add(gr2[0]->GetHistogram());
459
460   c1->cd(2);
461   gr2[1] = AliMathBase::MakeStat2D(comp->fD1TanSPtTPC,4,2,5); 
462   gr2[1]->GetXaxis()->SetTitle("Tan(#theta)");
463   gr2[1]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
464   gr2[1]->GetZaxis()->SetTitle("#sigmaDCA_z (cm)");
465   gr2[1]->SetName("DCAZResolSPTTanTPC");
466   gr2[1]->SetTitle("#sigma DCA_z (TPC only)");
467   gr2[1]->GetHistogram()->Draw("colz");
468
469   gr2[1]->GetHistogram()->SetName("DCAZResolSPTTanTPC");
470   aFolderObj->Add(gr2[1]->GetHistogram());
471
472   // TPC+ITS
473   c1->cd(3);
474   gr2[2] = AliMathBase::MakeStat2D(comp->fD0TanSPtTPCITS,4,2,5); 
475   gr2[2]->GetXaxis()->SetTitle("Tan(#theta)");
476   gr2[2]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
477   gr2[2]->GetZaxis()->SetTitle("#sigmaDCA_xy (cm)");
478   gr2[2]->SetName("DCAXYResolSPTTanTPCITS");
479   gr2[2]->SetTitle("#sigma DCA_xy (TPC+ITS)");
480   gr2[2]->GetHistogram()->Draw("colz");
481
482   gr2[2]->GetHistogram()->SetName("DCAXYResolSPTTanTPCITS");
483   aFolderObj->Add(gr2[2]->GetHistogram());
484
485   c1->cd(4);
486   gr2[3] = AliMathBase::MakeStat2D(comp->fD1TanSPtTPCITS,4,2,5); 
487   gr2[3]->GetXaxis()->SetTitle("Tan(#theta)");
488   gr2[3]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
489   gr2[3]->GetZaxis()->SetTitle("#sigmaDCA_z (cm)");
490   gr2[3]->SetName("DCAZResolSPTTanTPCITS");
491   gr2[3]->SetTitle("#sigma DCA_z (TPC+ITS)");
492   gr2[3]->GetHistogram()->Draw("colz");
493
494   gr2[3]->GetHistogram()->SetName("DCAZResolSPTTanTPCITS");
495   aFolderObj->Add(gr2[3]->GetHistogram());
496
497   // 2D DCA mean value  
498   c1->cd(5);
499   gr2[4] = AliMathBase::MakeStat2D(comp->fD0TanSPtTPC,4,2,4); 
500   gr2[4]->GetXaxis()->SetTitle("Tan(#theta)");
501   gr2[4]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
502   gr2[4]->GetZaxis()->SetTitle("mean DCA_xy (cm)");
503   gr2[4]->SetName("DCAXYMeanSPTTanTPC");
504   gr2[4]->SetTitle("mean DCA_xy (TPC only)");
505   gr2[4]->GetHistogram()->Draw("colz");
506
507   gr2[4]->GetHistogram()->SetName("DCAXYMeanSPTTanTPC");
508   aFolderObj->Add(gr2[4]->GetHistogram());
509
510   c1->cd(6);
511   gr2[5] = AliMathBase::MakeStat2D(comp->fD1TanSPtTPC,4,2,4); 
512   gr2[5]->GetXaxis()->SetTitle("Tan(#theta)");
513   gr2[5]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
514   gr2[5]->GetZaxis()->SetTitle("mean DCA_z (cm)");
515   gr2[5]->SetName("DCAZMeanSPTTanTPC");
516   gr2[5]->SetTitle("mean DCA_z (TPC only)");
517   gr2[5]->GetHistogram()->Draw("colz");
518
519   gr2[5]->GetHistogram()->SetName("DCAZMeanSPTTanTPC");
520   aFolderObj->Add(gr2[5]->GetHistogram());
521
522   c1->cd(7);
523   gr2[6] = AliMathBase::MakeStat2D(comp->fD0TanSPtTPCITS,4,2,4); 
524   gr2[6]->GetXaxis()->SetTitle("Tan(#theta)");
525   gr2[6]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
526   gr2[6]->GetZaxis()->SetTitle("mean DCA_xy (cm)");
527   gr2[6]->SetName("DCAXYMeanSPTTanTPCITS");
528   gr2[6]->SetTitle("mean DCA_xy (TPC+ITS)");
529   gr2[6]->GetHistogram()->Draw("colz");
530
531   gr2[6]->GetHistogram()->SetName("DCAXYMeanSPTTanTPCITS");
532   aFolderObj->Add(gr2[6]->GetHistogram());
533
534   c1->cd(8);
535   gr2[7] = AliMathBase::MakeStat2D(comp->fD1TanSPtTPCITS,4,2,4); 
536   gr2[7]->GetXaxis()->SetTitle("Tan(#theta)");
537   gr2[7]->GetYaxis()->SetTitle("#sqrt{p_{t}(GeV)}");
538   gr2[7]->GetZaxis()->SetTitle("mean DCA_z (cm)");
539   gr2[7]->SetName("DCAZMeanSPTTanTPCITS");
540   gr2[7]->SetTitle("mean DCA_z (TPC+ITS)");
541   gr2[7]->GetHistogram()->Draw("colz");
542
543   gr2[7]->GetHistogram()->SetName("DCAZMeanSPTTanTPCITS");
544   aFolderObj->Add(gr2[7]->GetHistogram());
545
546   */
547   // export objects to analysis folder
548   fAnalysisFolder = ExportToFolder(aFolderObj);
549
550   // delete only TObjArray
551   if(aFolderObj) delete aFolderObj;
552 }
553
554 //_____________________________________________________________________________
555 TFolder* AliComparisonDCA::ExportToFolder(TObjArray * array) 
556 {
557   // recreate folder avery time and export objects to new one
558   //
559   AliComparisonDCA * comp=this;
560   TFolder *folder = comp->GetAnalysisFolder();
561
562   TString name, title;
563   TFolder *newFolder = 0;
564   Int_t i = 0;
565   Int_t size = array->GetSize();
566
567   if(folder) { 
568      // get name and title from old folder
569      name = folder->GetName();  
570      title = folder->GetTitle();  
571
572          // delete old one
573      delete folder;
574
575          // create new one
576      newFolder = CreateFolder(name.Data(),title.Data());
577      newFolder->SetOwner();
578
579          // add objects to folder
580      while(i < size) {
581            newFolder->Add(array->At(i));
582            i++;
583          }
584   }
585
586 return newFolder;
587 }
588
589
590 //_____________________________________________________________________________
591 TFolder* AliComparisonDCA::CreateFolder(TString name,TString title) { 
592 // create folder for analysed histograms
593 TFolder *folder = 0;
594   folder = new TFolder(name.Data(),title.Data());
595
596   return folder;
597 }