]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdPt/AlidNdPtCutAnalysis.cxx
c74923839a9b50ab7b35cc26ee32ee71c79645b0
[u/mrichter/AliRoot.git] / PWG0 / dNdPt / AlidNdPtCutAnalysis.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 \r
16 #include <iostream>\r
17 \r
18 #include "TFile.h"\r
19 #include "TCint.h"\r
20 #include "TH1.h"\r
21 #include "TH2.h"\r
22 \r
23 #include "AliHeader.h"  \r
24 #include "AliGenEventHeader.h"  \r
25 #include "AliStack.h"  \r
26 #include "AliESDEvent.h"  \r
27 #include "AliMCEvent.h"  \r
28 #include "AliESDtrackCuts.h"  \r
29 #include "AliLog.h" \r
30 \r
31 #include "AlidNdPtEventCuts.h"\r
32 #include "AlidNdPtAcceptanceCuts.h"\r
33 \r
34 #include "AliPWG0Helper.h"\r
35 #include "AlidNdPtHelper.h"\r
36 #include "AlidNdPtCutAnalysis.h"\r
37 \r
38 using namespace std;\r
39 \r
40 ClassImp(AlidNdPtCutAnalysis)\r
41 \r
42 //_____________________________________________________________________________\r
43   AlidNdPtCutAnalysis::AlidNdPtCutAnalysis(): AlidNdPt(),\r
44   fAnalysisFolder(0),\r
45   fRecEventHist(0),\r
46   fMCEventHist(0),\r
47   fRecMCEventHist(0),\r
48   fRecMCTrackHist(0)\r
49 {\r
50   // default constructor\r
51   Init();\r
52 }\r
53 \r
54 //_____________________________________________________________________________\r
55 AlidNdPtCutAnalysis::AlidNdPtCutAnalysis(Char_t* name, Char_t* title): AlidNdPt(name,title),\r
56   fAnalysisFolder(0),\r
57   fRecEventHist(0),\r
58   fMCEventHist(0),\r
59   fRecMCEventHist(0),\r
60   fRecMCTrackHist(0)\r
61 {\r
62   // constructor\r
63   Init();\r
64 }\r
65 \r
66 //_____________________________________________________________________________\r
67 AlidNdPtCutAnalysis::~AlidNdPtCutAnalysis() {\r
68   // \r
69   if(fRecEventHist) delete fRecEventHist; fRecEventHist=0;\r
70   if(fMCEventHist) delete fMCEventHist; fMCEventHist=0;\r
71   if(fRecMCEventHist) delete fRecMCEventHist; fRecMCEventHist=0;\r
72   if(fRecMCTrackHist) delete fRecMCTrackHist; fRecMCTrackHist=0;\r
73 \r
74   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;\r
75 }\r
76 \r
77 //_____________________________________________________________________________\r
78 void AlidNdPtCutAnalysis::Init(){\r
79   //\r
80   // Init histograms\r
81   //\r
82   const Int_t ptNbins = 56; \r
83   const Double_t ptMin = 0.; \r
84   const Double_t ptMax = 16.; \r
85 \r
86   Double_t binsPt[ptNbins+1] = {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.65,0.7,0.75,0.8,0.85,0.9,0.95,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0,2.2,2.4,2.6,2.8,3.0,3.2,3.4,3.6,3.8,4.0,4.5,5.0,5.5,6.0,6.5,7.0,7.5,8.0,9.0,10.0,11.0,12.0,13.0,14.0,15.0,16.0};\r
87 \r
88   //Xv:Yv:Zv:ResZv:Mult\r
89   Int_t binsRecEventHist[5]={100,100,140,100,150};\r
90   Double_t minRecEventHist[5]={-3.,-3.,-35.,0.,0.}; \r
91   Double_t maxRecEventHist[5]={3.,3.,35.,10.,150.}; \r
92   fRecEventHist = new THnSparseF("fRecEventHist","Xv:Yv:Zv:ResZv:Mult",5,binsRecEventHist,minRecEventHist,maxRecEventHist);\r
93   fRecEventHist->GetAxis(0)->SetTitle("Xv (cm)");\r
94   fRecEventHist->GetAxis(1)->SetTitle("Yv (cm)");\r
95   fRecEventHist->GetAxis(2)->SetTitle("Zv (cm)");\r
96   fRecEventHist->GetAxis(3)->SetTitle("ResZv (cm)");\r
97   fRecEventHist->GetAxis(4)->SetTitle("Mult");\r
98   fRecEventHist->Sumw2();\r
99 \r
100   //Xv:Yv:Zv\r
101   Int_t binsMCEventHist[3]={100,100,140};\r
102   Double_t minMCEventHist[3]={-0.1,-0.1,-35.}; \r
103   Double_t maxMCEventHist[3]={0.1,0.1,35.}; \r
104   fMCEventHist = new THnSparseF("fMCEventHist","mcXv:mcYv:mcZv",3,binsMCEventHist,minMCEventHist,maxMCEventHist);\r
105   fMCEventHist->GetAxis(0)->SetTitle("mcXv (cm)");\r
106   fMCEventHist->GetAxis(1)->SetTitle("mcYv (cm)");\r
107   fMCEventHist->GetAxis(2)->SetTitle("mcZv (cm)");\r
108   fMCEventHist->Sumw2();\r
109 \r
110   //Xv-mcXv:Yv-mcYv:Zv-mcZv:Mult\r
111   Int_t binsRecMCEventHist[4]={100,100,100,150};\r
112   Double_t minRecMCEventHist[4]={-10.0,-10.0,-10.0,0.}; \r
113   Double_t maxRecMCEventHist[4]={10.0,10.0,10.0,150.}; \r
114   fRecMCEventHist = new THnSparseF("fRecMCEventHist","mcXv-Xv:mcYv-Yv:mcZv-Zv:Mult",4,binsRecMCEventHist,minRecMCEventHist,maxRecMCEventHist);\r
115   fRecMCEventHist->GetAxis(0)->SetTitle("mcXv-Xv (cm)");\r
116   fRecMCEventHist->GetAxis(1)->SetTitle("mcYv-Yv (cm)");\r
117   fRecMCEventHist->GetAxis(2)->SetTitle("mcZv-Zv (cm)");\r
118   fRecMCEventHist->GetAxis(3)->SetTitle("Mult");\r
119   fRecMCEventHist->Sumw2();\r
120 \r
121   //\r
122   // THnSparse track histograms\r
123   //\r
124 \r
125   //nClust:chi2PerClust:nClust/nFindableClust:DCAy:DCAz:eta:phi:pt:isKink:isPrim\r
126   Int_t binsRecMCTrackHist[10]={160,100,100,100,100,30,90,ptNbins, 2,2};\r
127   Double_t minRecMCTrackHist[10]={0., 0., 0., -10.,-10.,-1.5, 0., ptMin, 0., 0.};\r
128   Double_t maxRecMCTrackHist[10]={160.,10.,1.2, 10.,10.,1.5, 2.*TMath::Pi(), ptMax, 2.,2.};\r
129 \r
130   fRecMCTrackHist = new THnSparseF("fRecMCTrackHist","nClust:chi2PerClust:nClust/nFindableClust:DCAy:DCAz:eta:phi:pt:isKink:isPrim",10,binsRecMCTrackHist,minRecMCTrackHist,maxRecMCTrackHist);\r
131   fRecMCTrackHist->SetBinEdges(7,binsPt);\r
132 \r
133   fRecMCTrackHist->GetAxis(0)->SetTitle("nClust");\r
134   fRecMCTrackHist->GetAxis(1)->SetTitle("chi2PerClust");\r
135   fRecMCTrackHist->GetAxis(2)->SetTitle("nClust/nFindableClust");\r
136   fRecMCTrackHist->GetAxis(3)->SetTitle("DCAy (cm)");\r
137   fRecMCTrackHist->GetAxis(4)->SetTitle("DCAz (cm)");\r
138   fRecMCTrackHist->GetAxis(5)->SetTitle("#eta");\r
139   fRecMCTrackHist->GetAxis(6)->SetTitle("#phi (rad)");\r
140   fRecMCTrackHist->GetAxis(7)->SetTitle("p_{T} (GeV/c)");\r
141   fRecMCTrackHist->GetAxis(8)->SetTitle("isKink");\r
142   fRecMCTrackHist->GetAxis(9)->SetTitle("isPrim");\r
143   fRecMCTrackHist->Sumw2();\r
144 \r
145   // init output folder\r
146   fAnalysisFolder = CreateFolder("folderdNdPt","Analysis dNdPt Folder");\r
147 \r
148 }\r
149 \r
150 //_____________________________________________________________________________\r
151 void AlidNdPtCutAnalysis::Process(AliESDEvent *const esdEvent, AliMCEvent * const mcEvent)\r
152 {\r
153   //\r
154   // Process real and/or simulated events\r
155   //\r
156   if(!esdEvent) {\r
157     AliDebug(AliLog::kError, "esdEvent not available");\r
158     return;\r
159   }\r
160 \r
161   // get selection cuts\r
162   AlidNdPtEventCuts *evtCuts = GetEventCuts(); \r
163   AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
164   AliESDtrackCuts *esdTrackCuts = GetTrackCuts(); \r
165 \r
166   if(!evtCuts || !accCuts  || !esdTrackCuts) {\r
167     AliDebug(AliLog::kError, "cuts not available");\r
168     return;\r
169   }\r
170 \r
171   // trigger selection\r
172   Bool_t isEventTriggered = kTRUE;\r
173   if(evtCuts->IsTriggerRequired())  {\r
174     static AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;\r
175     isEventTriggered = triggerAnalysis->IsTriggerFired(esdEvent, GetTrigger());\r
176   }\r
177 \r
178   // use MC information\r
179   AliHeader* header = 0;\r
180   AliGenEventHeader* genHeader = 0;\r
181   AliStack* stack = 0;\r
182   TArrayF vtxMC(3);\r
183   //AlidNdPtHelper::MCProcessType evtType = AlidNdPtHelper::kInvalidProcess;\r
184   AliPWG0Helper::MCProcessType evtType = AliPWG0Helper::kInvalidProcess;\r
185 \r
186   if(IsUseMCInfo())\r
187   {\r
188     if(!mcEvent) {\r
189       AliDebug(AliLog::kError, "mcEvent not available");\r
190       return;\r
191     }\r
192 \r
193     // get MC event header\r
194     header = mcEvent->Header();\r
195     if (!header) {\r
196       AliDebug(AliLog::kError, "Header not available");\r
197       return;\r
198     }\r
199 \r
200     // MC particle stack\r
201     stack = mcEvent->Stack();\r
202     if (!stack) {\r
203       AliDebug(AliLog::kError, "Stack not available");\r
204       return;\r
205     }\r
206 \r
207     // get event type (ND=0x1, DD=0x2, SD=0x4)\r
208     evtType = AliPWG0Helper::GetEventProcessType(header);\r
209     AliDebug(AliLog::kDebug+1, Form("Found process type %d", evtType));\r
210 \r
211     // get MC vertex\r
212     genHeader = header->GenEventHeader();\r
213     if (!genHeader) {\r
214       AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");\r
215       return;\r
216     }\r
217     genHeader->PrimaryVertex(vtxMC);\r
218 \r
219     // Fill MC event histogram\r
220     Double_t vMCEventHist[3]={vtxMC[0],vtxMC[1],vtxMC[2]};\r
221     fMCEventHist->Fill(vMCEventHist);\r
222 \r
223   } // end bUseMC\r
224 \r
225   // get reconstructed vertex  \r
226   Bool_t bRedoTPCVertex = evtCuts->IsRedoTPCVertex();\r
227   Bool_t bUseConstraints = evtCuts->IsUseBeamSpotConstraint();\r
228   const AliESDVertex* vtxESD = AlidNdPtHelper::GetVertex(esdEvent,evtCuts,accCuts,esdTrackCuts,GetAnalysisMode(),kFALSE,bRedoTPCVertex,bUseConstraints); \r
229   if(!vtxESD) return; \r
230 \r
231   Bool_t isRecVertex = AlidNdPtHelper::TestRecVertex(vtxESD, GetAnalysisMode(), kFALSE);\r
232   Bool_t isEventOK = evtCuts->AcceptEvent(esdEvent,mcEvent,vtxESD) && isRecVertex;\r
233 \r
234   TObjArray *allChargedTracks=0;\r
235   Int_t multAll=0;\r
236 \r
237   // check event cuts\r
238   if(isEventOK && isEventTriggered)\r
239   {\r
240     // get all charged tracks\r
241     allChargedTracks = AlidNdPtHelper::GetAllChargedTracks(esdEvent,GetAnalysisMode());\r
242     if(!allChargedTracks) return;\r
243 \r
244     Int_t entries = allChargedTracks->GetEntries();\r
245     for(Int_t i=0; i<entries;++i) \r
246     {\r
247       AliESDtrack *track = (AliESDtrack*)allChargedTracks->At(i);\r
248       if(!track) continue;\r
249 \r
250       FillHistograms(track, stack);\r
251       multAll++;\r
252     }\r
253 \r
254   Double_t vRecEventHist[5] = {vtxESD->GetXv(),vtxESD->GetYv(),vtxESD->GetZv(),vtxESD->GetZRes(),multAll};\r
255   fRecEventHist->Fill(vRecEventHist);\r
256 \r
257   if(IsUseMCInfo()) {\r
258     Double_t vRecMCEventHist[5] = {vtxESD->GetXv()-vtxMC[0],vtxESD->GetYv()-vtxMC[1],vtxESD->GetZv()-vtxMC[2],multAll};\r
259     fRecMCEventHist->Fill(vRecMCEventHist);\r
260   }\r
261   }\r
262 \r
263   if(allChargedTracks) delete allChargedTracks; allChargedTracks = 0;\r
264 }\r
265 \r
266 //_____________________________________________________________________________\r
267 void AlidNdPtCutAnalysis::FillHistograms(AliESDtrack *const esdTrack, AliStack *const stack)\r
268 {\r
269   //\r
270   // Fill ESD track and MC histograms \r
271   //\r
272   if(!esdTrack) return;\r
273   if(esdTrack->Charge() == 0.) return;\r
274 \r
275   Float_t pt = esdTrack->Pt();\r
276   Float_t eta = esdTrack->Eta();\r
277   Float_t phi = esdTrack->Phi();\r
278   Int_t nClust = esdTrack->GetTPCclusters(0);\r
279   Int_t nFindableClust = esdTrack->GetTPCNclsF();\r
280 \r
281   Float_t chi2PerCluster = 0.;\r
282   if(nClust>0.) chi2PerCluster = esdTrack->GetTPCchi2()/Float_t(nClust);\r
283 \r
284   Float_t clustPerFindClust = 0.;\r
285   if(nFindableClust>0.) clustPerFindClust = Float_t(nClust)/nFindableClust;\r
286 \r
287   Float_t b[2], bCov[3];\r
288   esdTrack->GetImpactParameters(b,bCov);\r
289 \r
290   // kink daughter\r
291   Bool_t isKink = kFALSE;\r
292   if(esdTrack->GetKinkIndex(0)>0) isKink=kTRUE;\r
293 \r
294   //\r
295   // Fill rec vs MC information\r
296   //\r
297 \r
298   Bool_t isPrim = kTRUE;\r
299 \r
300   if(IsUseMCInfo()) {\r
301     if(!stack) return;\r
302     Int_t label = TMath::Abs(esdTrack->GetLabel()); \r
303     TParticle* particle = stack->Particle(label);\r
304     if(!particle) return;\r
305     if(particle->GetPDG() && particle->GetPDG()->Charge()==0.) return;\r
306     isPrim = stack->IsPhysicalPrimary(label);\r
307   }\r
308 \r
309   // fill histo\r
310   Double_t vRecMCTrackHist[10] = {nClust,chi2PerCluster,clustPerFindClust,b[0],b[1],eta,phi,pt,isKink,isPrim}; \r
311   fRecMCTrackHist->Fill(vRecMCTrackHist);\r
312 }\r
313 \r
314 //_____________________________________________________________________________\r
315 Long64_t AlidNdPtCutAnalysis::Merge(TCollection* list) \r
316 {\r
317   // Merge list of objects (needed by PROOF)\r
318 \r
319   if (!list)\r
320   return 0;\r
321 \r
322   if (list->IsEmpty())\r
323   return 1;\r
324 \r
325   TIterator* iter = list->MakeIterator();\r
326   TObject* obj = 0;\r
327 \r
328   // collection of generated histograms\r
329   Int_t count=0;\r
330   while((obj = iter->Next()) != 0) {\r
331     AlidNdPtCutAnalysis* entry = dynamic_cast<AlidNdPtCutAnalysis*>(obj);\r
332     if (entry == 0) continue; \r
333   \r
334     // event histo\r
335     fRecEventHist->Add(entry->fRecEventHist);\r
336     fRecMCEventHist->Add(entry->fRecMCEventHist);\r
337     fMCEventHist->Add(entry->fMCEventHist);\r
338 \r
339     // track histo\r
340     fRecMCTrackHist->Add(entry->fRecMCTrackHist);\r
341     \r
342   count++;\r
343   }\r
344 \r
345 return count;\r
346 }\r
347 \r
348 //_____________________________________________________________________________\r
349 void AlidNdPtCutAnalysis::Analyse() \r
350 {\r
351   //\r
352   // Analyse histograms\r
353   //\r
354   TH1::AddDirectory(kFALSE);\r
355   TObjArray *aFolderObj = new TObjArray;\r
356   TH1D *h1D = 0; \r
357   TH2D *h2D = 0; \r
358 \r
359 \r
360   //\r
361   // get cuts\r
362   //\r
363   AlidNdPtEventCuts *evtCuts = GetEventCuts(); \r
364   AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
365   AliESDtrackCuts *esdTrackCuts = GetTrackCuts(); \r
366 \r
367   if(!evtCuts || !accCuts || !esdTrackCuts) {\r
368     Error("AlidNdPtCutAnalysis::Analyse()", "cuts not available");\r
369     return;\r
370   }\r
371 \r
372   //\r
373   // set min and max values\r
374   //\r
375   Double_t minPt = accCuts->GetMinPt();\r
376   Double_t maxPt = accCuts->GetMaxPt();\r
377   Double_t minEta = accCuts->GetMinEta();\r
378   Double_t maxEta = accCuts->GetMaxEta()-0.00001;\r
379 \r
380   Double_t maxDCAr = accCuts->GetMaxDCAr();\r
381 \r
382   //\r
383   // Create rec. event histograms\r
384   //\r
385   h2D = (TH2D *)fRecEventHist->Projection(0,1);\r
386   h2D->SetName("rec_xv_vs_yv");\r
387   aFolderObj->Add(h2D);\r
388 \r
389   h2D = (TH2D *)fRecEventHist->Projection(0,2);\r
390   h2D->SetName("rec_xv_vs_zv");\r
391   aFolderObj->Add(h2D);\r
392 \r
393   h2D = (TH2D *)fRecEventHist->Projection(3,4);\r
394   h2D->SetName("rec_resZv_vs_Mult");\r
395   aFolderObj->Add(h2D);\r
396 \r
397 \r
398   //\r
399   // MC available\r
400   //\r
401   if(IsUseMCInfo()) {\r
402 \r
403   //\r
404   // Create mc event histograms\r
405   //\r
406   h2D = (TH2D *)fMCEventHist->Projection(0,1);\r
407   h2D->SetName("mc_xv_vs_yv");\r
408   aFolderObj->Add(h2D);\r
409 \r
410   h2D = (TH2D *)fMCEventHist->Projection(0,2);\r
411   h2D->SetName("mc_xv_vs_zv");\r
412   aFolderObj->Add(h2D);\r
413 \r
414   //\r
415   // Create rec-mc event histograms\r
416   //\r
417   h2D = (TH2D *)fRecMCEventHist->Projection(0,3);\r
418   h2D->SetName("rec_mc_deltaXv_vs_mult");\r
419   aFolderObj->Add(h2D);\r
420 \r
421   h2D = (TH2D *)fRecMCEventHist->Projection(1,3);\r
422   h2D->SetName("rec_mc_deltaYv_vs_mult");\r
423   aFolderObj->Add(h2D);\r
424 \r
425   h2D = (TH2D *)fRecMCEventHist->Projection(2,3);\r
426   h2D->SetName("rec_mc_deltaZv_vs_mult");\r
427   aFolderObj->Add(h2D);\r
428 \r
429   } // end use MC info \r
430 \r
431 \r
432 \r
433   //\r
434   // Create rec-mc track track histograms \r
435   //\r
436 \r
437   // DCA cuts\r
438   fRecMCTrackHist->GetAxis(3)->SetRangeUser(-maxDCAr,maxDCAr);\r
439   fRecMCTrackHist->GetAxis(4)->SetRangeUser(-maxDCAr,maxDCAr);\r
440 \r
441   h2D = (TH2D *)fRecMCTrackHist->Projection(7,5);\r
442   h2D->SetName("pt_vs_eta");\r
443   aFolderObj->Add(h2D);\r
444 \r
445   fRecMCTrackHist->GetAxis(7)->SetRangeUser(minPt,maxPt);  \r
446 \r
447   h2D = (TH2D *)fRecMCTrackHist->Projection(0,5);\r
448   h2D->SetName("nClust_vs_eta");\r
449   aFolderObj->Add(h2D);\r
450 \r
451   h2D = (TH2D *)fRecMCTrackHist->Projection(1,5);\r
452   h2D->SetName("chi2PerClust_vs_eta");\r
453   aFolderObj->Add(h2D);\r
454 \r
455   h2D = (TH2D *)fRecMCTrackHist->Projection(2,5);\r
456   h2D->SetName("ratio_nClust_nFindableClust_vs_eta");\r
457   aFolderObj->Add(h2D);\r
458 \r
459   //\r
460   fRecMCTrackHist->GetAxis(7)->SetRangeUser(minEta,maxEta);  \r
461 \r
462   h2D = (TH2D *)fRecMCTrackHist->Projection(0,6);\r
463   h2D->SetName("nClust_vs_phi");\r
464   aFolderObj->Add(h2D);\r
465 \r
466   h2D = (TH2D *)fRecMCTrackHist->Projection(1,6);\r
467   h2D->SetName("chi2PerClust_vs_phi");\r
468   aFolderObj->Add(h2D);\r
469 \r
470   h2D = (TH2D *)fRecMCTrackHist->Projection(2,6);\r
471   h2D->SetName("ratio_nClust_nFindableClust_vs_phi");\r
472   aFolderObj->Add(h2D);\r
473 \r
474   h2D = (TH2D *)fRecMCTrackHist->Projection(5,6);\r
475   h2D->SetName("eta_vs_phi");\r
476   aFolderObj->Add(h2D);\r
477 \r
478   //\r
479   fRecMCTrackHist->GetAxis(7)->SetRangeUser(0.0,maxPt);  \r
480 \r
481   h2D = (TH2D *)fRecMCTrackHist->Projection(0,7);\r
482   h2D->SetName("nClust_vs_pt");\r
483   aFolderObj->Add(h2D);\r
484 \r
485   h2D = (TH2D *)fRecMCTrackHist->Projection(1,7);\r
486   h2D->SetName("chi2PerClust_vs_pt");\r
487   aFolderObj->Add(h2D);\r
488 \r
489   h2D = (TH2D *)fRecMCTrackHist->Projection(2,7);\r
490   h2D->SetName("ratio_nClust_nFindableClust_vs_pt");\r
491   aFolderObj->Add(h2D);\r
492 \r
493   h2D = (TH2D *)fRecMCTrackHist->Projection(6,7);\r
494   h2D->SetName("phi_vs_pt");\r
495   aFolderObj->Add(h2D);\r
496 \r
497 \r
498   // fiducial volume\r
499   fRecMCTrackHist->GetAxis(5)->SetRangeUser(minEta,maxEta);  \r
500   fRecMCTrackHist->GetAxis(7)->SetRangeUser(minPt,maxPt);  \r
501 \r
502   // DCA cuts\r
503   fRecMCTrackHist->GetAxis(3)->SetRangeUser(-maxDCAr,maxDCAr);\r
504   fRecMCTrackHist->GetAxis(4)->SetRangeUser(-maxDCAr,maxDCAr);\r
505 \r
506   h2D = (TH2D *)fRecMCTrackHist->Projection(0,1);\r
507   h2D->SetName("nClust_vs_chi2PerClust");\r
508   aFolderObj->Add(h2D);\r
509 \r
510   h2D = (TH2D *)fRecMCTrackHist->Projection(0,2);\r
511   h2D->SetName("nClust_vs_ratio_nClust_nFindableClust");\r
512   aFolderObj->Add(h2D);\r
513 \r
514   // DCAy cuts\r
515   fRecMCTrackHist->GetAxis(3)->SetRange(1,fRecMCTrackHist->GetAxis(3)->GetNbins());\r
516   fRecMCTrackHist->GetAxis(4)->SetRangeUser(-1.0,1.0);\r
517 \r
518   // sec\r
519   fRecMCTrackHist->GetAxis(9)->SetRange(1,1);\r
520   h1D = (TH1D *)fRecMCTrackHist->Projection(3);\r
521   h1D->SetName("dcay_sec");\r
522   aFolderObj->Add(h1D);\r
523 \r
524   // prim\r
525   fRecMCTrackHist->GetAxis(9)->SetRange(2,2);\r
526   h1D = (TH1D *)fRecMCTrackHist->Projection(3);\r
527   h1D->SetName("dcay_prim");\r
528   aFolderObj->Add(h1D);\r
529 \r
530   // DCAz cuts\r
531   fRecMCTrackHist->GetAxis(3)->SetRangeUser(-1.0,1.0);\r
532   fRecMCTrackHist->GetAxis(4)->SetRange(1,fRecMCTrackHist->GetAxis(4)->GetNbins());\r
533 \r
534   // sec\r
535   fRecMCTrackHist->GetAxis(9)->SetRange(1,1);\r
536   h1D = (TH1D *)fRecMCTrackHist->Projection(4);\r
537   h1D->SetName("dcaz_sec");\r
538   aFolderObj->Add(h1D);\r
539 \r
540   // prim\r
541   fRecMCTrackHist->GetAxis(9)->SetRange(2,2);\r
542   h1D = (TH1D *)fRecMCTrackHist->Projection(4);\r
543   h1D->SetName("dcaz_prim");\r
544   aFolderObj->Add(h1D);\r
545 \r
546 \r
547   // export objects to analysis folder\r
548   fAnalysisFolder = ExportToFolder(aFolderObj);\r
549 \r
550   // delete only TObjArray\r
551   if(aFolderObj) delete aFolderObj;\r
552 }\r
553 \r
554 //_____________________________________________________________________________\r
555 TFolder* AlidNdPtCutAnalysis::ExportToFolder(TObjArray * array) \r
556 {\r
557   // recreate folder avery time and export objects to new one\r
558   //\r
559   AlidNdPtCutAnalysis * comp=this;\r
560   TFolder *folder = comp->GetAnalysisFolder();\r
561 \r
562   TString name, title;\r
563   TFolder *newFolder = 0;\r
564   Int_t i = 0;\r
565   Int_t size = array->GetSize();\r
566 \r
567   if(folder) { \r
568      // get name and title from old folder\r
569      name = folder->GetName();  \r
570      title = folder->GetTitle();  \r
571 \r
572          // delete old one\r
573      delete folder;\r
574 \r
575          // create new one\r
576      newFolder = CreateFolder(name.Data(),title.Data());\r
577      newFolder->SetOwner();\r
578 \r
579          // add objects to folder\r
580      while(i < size) {\r
581            newFolder->Add(array->At(i));\r
582            i++;\r
583          }\r
584   }\r
585 \r
586 return newFolder;\r
587 }\r
588 \r
589 //_____________________________________________________________________________\r
590 TFolder* AlidNdPtCutAnalysis::CreateFolder(TString name,TString title) { \r
591 // create folder for analysed histograms\r
592 //\r
593 TFolder *folder = 0;\r
594   folder = new TFolder(name.Data(),title.Data());\r
595 \r
596   return folder;\r
597 }\r