]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliPerformanceDEdx.cxx
Update of the alignment-data fileting macro including a fix for the access to the...
[u/mrichter/AliRoot.git] / PWG1 / 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 //------------------------------------------------------------------------------
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   //AliPerformanceDEdx * compObj = (AliPerformanceDEdx*)f.Get("AliPerformanceDEdx");
20   AliPerformanceDEdx * compObj = (AliPerformanceDEdx*)coutput->FindObject("AliPerformanceDEdx");
21
22   // Analyse comparison data
23   compObj->Analyse();
24
25   // the output histograms/graphs will be stored in the folder "folderDEdx" 
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_DEdx.root"."recreate");
31   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
32   fout.Close();
33
34 */
35
36 #include <TDirectory.h>
37 #include <TAxis.h>
38 #include <TCanvas.h>
39 #include <TH1.h>
40 #include <TH2.h>
41 #include <TF1.h>
42
43 #include "AliPerformanceDEdx.h" 
44 #include "AliESDEvent.h"
45 #include "AliMCEvent.h"
46 #include "AliESDtrack.h"
47 #include "AliStack.h"
48 #include "AliLog.h" 
49 #include "AliMCInfoCuts.h" 
50 #include "AliMathBase.h"
51 #include "AliRecInfoCuts.h" 
52 #include "AliTreeDraw.h"
53 #include "AliHeader.h"
54 #include "AliGenEventHeader.h"
55
56 using namespace std;
57
58 ClassImp(AliPerformanceDEdx)
59
60 //_____________________________________________________________________________
61 AliPerformanceDEdx::AliPerformanceDEdx():
62   AliPerformanceObject("AliPerformanceDEdx"),
63
64   // dEdx 
65   fDeDxHisto(0),
66   
67   // Cuts 
68   fCutsRC(0), 
69   fCutsMC(0),
70
71   // histogram folder 
72   fAnalysisFolder(0)
73 {
74   // default constructor        
75   Init();
76 }
77
78 //_____________________________________________________________________________
79 AliPerformanceDEdx::AliPerformanceDEdx(Char_t* name="AliPerformanceDEdx", Char_t* title="AliPerformanceDEdx",Int_t analysisMode=0, Bool_t hptGenerator=kFALSE):
80   AliPerformanceObject(name,title),
81
82   // dEdx 
83   fDeDxHisto(0),
84   
85   // Cuts 
86   fCutsRC(0), 
87   fCutsMC(0),
88
89   // histogram folder 
90   fAnalysisFolder(0)
91 {
92   // named constructor
93
94   SetAnalysisMode(analysisMode);
95   SetHptGenerator(hptGenerator);
96   Init();
97 }
98
99
100 //_____________________________________________________________________________
101 AliPerformanceDEdx::~AliPerformanceDEdx()
102 {
103   // destructor
104   if(fDeDxHisto)  delete fDeDxHisto; fDeDxHisto=0; 
105   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
106 }
107
108 //_____________________________________________________________________________
109 void AliPerformanceDEdx::Init()
110 {
111   // Init histograms
112   
113   // TPC dEdx
114   // set pt bins
115   Int_t nPBins = 50;
116   Double_t pMin = 1.e-2, pMax = 10.;
117
118   Double_t *binsP = 0;
119   if (IsHptGenerator())  { 
120     nPBins = 100; pMax = 100.;
121     binsP = CreateLogAxis(nPBins,pMin,pMax);
122   } else {
123     binsP = CreateLogAxis(nPBins,pMin,pMax);
124   }
125
126
127   /*
128   Int_t nPBins = 31;
129     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.};
130     Double_t pMin = 0., pMax = 10.;
131
132     if(IsHptGenerator() == kTRUE) {
133       nPBins = 100;
134       pMin = 0.; pMax = 100.;
135     }
136    */
137
138    //dedx:alpha:y:z:snp:tgl:ncls:p
139    //dedx:phi:y:z:snp:tgl:ncls:p
140    //Int_t binsQA[8]    = {300, 50, 50,  50, 50, 50, 80, nPBins};
141    //Double_t xminQA[8] = {0, -4,-20,-250, -1, -2, 0, pMin};
142    //Double_t xmaxQA[8] = {300, 4, 20, 250,  1,  2, 160, pMax};
143    Int_t binsQA[8]    = {300, 144, 50,  50, 50, 50, 80, nPBins};
144    Double_t xminQA[8] = {0, -TMath::Pi(),-20,-250, -1, -2, 0, pMin};
145    Double_t xmaxQA[8] = {300, TMath::Pi(), 20, 250,  1,  2, 160, pMax};
146
147    //fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:alpha:y:z:snp:tgl:ncls:momentum",8,binsQA,xminQA,xmaxQA);
148    fDeDxHisto = new THnSparseF("fDeDxHisto","dedx:phi:y:z:snp:tgl:ncls:momentum",8,binsQA,xminQA,xmaxQA);
149    fDeDxHisto->SetBinEdges(7,binsP);
150
151    fDeDxHisto->GetAxis(0)->SetTitle("dedx (a.u.)");
152    fDeDxHisto->GetAxis(1)->SetTitle("#phi (rad)");
153    fDeDxHisto->GetAxis(2)->SetTitle("y (cm)");
154    fDeDxHisto->GetAxis(3)->SetTitle("z (cm)");
155    fDeDxHisto->GetAxis(4)->SetTitle("sin#phi");
156    fDeDxHisto->GetAxis(5)->SetTitle("tan#lambda");
157    fDeDxHisto->GetAxis(6)->SetTitle("ncls");
158    fDeDxHisto->GetAxis(7)->SetTitle("p (GeV/c)");
159    fDeDxHisto->Sumw2();
160
161    // Init cuts
162    if(!fCutsMC) 
163      AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
164    if(!fCutsRC) 
165      AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
166
167    // init folder
168    fAnalysisFolder = CreateFolder("folderDEdx","Analysis de/dx Folder");
169 }
170
171 //_____________________________________________________________________________
172 void AliPerformanceDEdx::ProcessTPC(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
173 {
174   // Fill dE/dx  comparison information
175   AliDebug(AliLog::kWarning, "Warning: Not implemented");
176 }
177
178 //_____________________________________________________________________________
179 void AliPerformanceDEdx::ProcessInnerTPC(AliStack* const stack, AliESDtrack *const esdTrack)
180 {
181   if(!esdTrack) return;
182
183   const AliExternalTrackParam *innerParam =  esdTrack->GetInnerParam();
184   if(!innerParam) return;
185
186   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
187   esdTrack->GetImpactParametersTPC(dca,cov);
188
189   if((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return; // TPC refit
190   if(TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ()) 
191   { 
192     Float_t dedx = esdTrack->GetTPCsignal();
193     Int_t ncls = esdTrack->GetTPCNcls();
194
195     Double_t pt = innerParam->Pt();
196     Double_t lam = TMath::ATan2(innerParam->Pz(),innerParam->Pt());
197     Double_t p = pt/TMath::Cos(lam);
198     //Double_t alpha = innerParam->GetAlpha();
199     Double_t phi = TMath::ATan2(innerParam->Py(),innerParam->Px());
200     //if(phi<0.) phi += 2.*TMath::Phi();
201     Double_t y = innerParam->GetY();
202     Double_t z = innerParam->GetZ();
203     Double_t snp = innerParam->GetSnp();
204     Double_t tgl = innerParam->GetTgl();
205
206     //Double_t vDeDxHisto[8] = {dedx,alpha,y,z,snp,tgl,ncls,p};
207     Double_t vDeDxHisto[8] = {dedx,phi,y,z,snp,tgl,ncls,p};
208     fDeDxHisto->Fill(vDeDxHisto); 
209   }
210
211   if(!stack) return;
212 }
213
214 //_____________________________________________________________________________
215 void AliPerformanceDEdx::ProcessTPCITS(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
216 {
217   // Fill dE/dx  comparison information
218   
219    AliDebug(AliLog::kWarning, "Warning: Not implemented");
220 }
221
222 //_____________________________________________________________________________
223 void AliPerformanceDEdx::ProcessConstrained(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/)
224 {
225   // Fill dE/dx  comparison information
226   
227    AliDebug(AliLog::kWarning, "Warning: Not implemented");
228 }
229
230 //_____________________________________________________________________________
231 Long64_t AliPerformanceDEdx::Merge(TCollection* const list) 
232 {
233   // Merge list of objects (needed by PROOF)
234
235   if (!list)
236   return 0;
237
238   if (list->IsEmpty())
239   return 1;
240
241   TIterator* iter = list->MakeIterator();
242   TObject* obj = 0;
243
244   // collection of generated histograms
245   Int_t count=0;
246   while((obj = iter->Next()) != 0) 
247   {
248     AliPerformanceDEdx* entry = dynamic_cast<AliPerformanceDEdx*>(obj);
249     if (entry == 0) continue;
250
251     fDeDxHisto->Add(entry->fDeDxHisto);
252     count++;
253   }
254
255 return count;
256 }
257
258 //_____________________________________________________________________________
259 void AliPerformanceDEdx::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
260 {
261   // Process comparison information 
262   //
263   if(!esdEvent) 
264   {
265       AliDebug(AliLog::kError, "esdEvent not available");
266       return;
267   }
268   AliHeader* header = 0;
269   AliGenEventHeader* genHeader = 0;
270   AliStack* stack = 0;
271   TArrayF vtxMC(3);
272   
273   if(bUseMC)
274   {
275     if(!mcEvent) {
276       AliDebug(AliLog::kError, "mcEvent not available");
277       return;
278     }
279
280     // get MC event header
281     header = mcEvent->Header();
282     if (!header) {
283       AliDebug(AliLog::kError, "Header not available");
284       return;
285     }
286     // MC particle stack
287     stack = mcEvent->Stack();
288     if (!stack) {
289       AliDebug(AliLog::kError, "Stack not available");
290       return;
291     }
292
293     // get MC vertex
294     genHeader = header->GenEventHeader();
295     if (!genHeader) {
296       AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");
297       return;
298     }
299     genHeader->PrimaryVertex(vtxMC);
300
301   } // end bUseMC
302
303   // use ESD friends
304   if(bUseESDfriend) {
305     if(!esdFriend) {
306       AliDebug(AliLog::kError, "esdFriend not available");
307       return;
308     }
309   }
310
311
312
313
314   //  Process events
315   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
316   { 
317     AliESDtrack *track = esdEvent->GetTrack(iTrack);
318     if(!track) continue;
319
320     if(GetAnalysisMode() == 0) ProcessTPC(stack,track);
321     else if(GetAnalysisMode() == 1) ProcessTPCITS(stack,track);
322     else if(GetAnalysisMode() == 2) ProcessConstrained(stack,track);
323     else if(GetAnalysisMode() == 3) ProcessInnerTPC(stack,track);
324     else {
325       printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
326       return;
327     }
328   }
329 }
330
331 //_____________________________________________________________________________
332 void AliPerformanceDEdx::Analyse()
333 {
334   // Analyze comparison information and store output histograms
335   // in the folder "folderDEdx"
336   //
337   TH1::AddDirectory(kFALSE);
338   TH1F *h1D=0;
339   TH2F *h2D=0;
340   TObjArray *aFolderObj = new TObjArray;
341
342   char name[256];
343   char title[256];
344
345   for(Int_t i=1; i<8; i++) { 
346     //
347     h2D = (TH2F*)fDeDxHisto->Projection(0,i);
348
349     sprintf(name,"h_dedx_%d_vs_%d",0,i);
350     h2D->SetName(name);
351     sprintf(title,"%s vs %s",fDeDxHisto->GetAxis(0)->GetTitle(),fDeDxHisto->GetAxis(i)->GetTitle());
352     h2D->SetTitle(title);
353     h2D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
354     h2D->GetYaxis()->SetTitle(fDeDxHisto->GetAxis(0)->GetTitle());
355
356     if(i==7) h2D->SetBit(TH1::kLogX);
357     aFolderObj->Add(h2D);
358   }
359
360   // resolution histograms for mips
361   fDeDxHisto->GetAxis(2)->SetRangeUser(-15.,14.999);
362   fDeDxHisto->GetAxis(3)->SetRangeUser(-120.,119.999);
363   fDeDxHisto->GetAxis(4)->SetRangeUser(-0.4, 0.399);
364   fDeDxHisto->GetAxis(5)->SetRangeUser(-1.,0.999);
365   fDeDxHisto->GetAxis(6)->SetRangeUser(60.,140.);
366   fDeDxHisto->GetAxis(7)->SetRangeUser(0.4,0.499);
367
368   h1D=(TH1F*)fDeDxHisto->Projection(0);
369   h1D->SetName("dedx_mips");
370
371   h1D->SetTitle("dedx_mips");
372   h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(0)->GetTitle());
373   aFolderObj->Add(h1D);
374
375   //
376   TObjArray *arr[7] = {0};
377   TF1 *f1[7] = {0};
378   
379   for(Int_t i=1; i<8; i++) 
380   { 
381     arr[i] = new TObjArray;
382     f1[i] = new TF1("gaus","gaus");
383     //printf("i %d \n",i);
384
385     h2D = (TH2F*)fDeDxHisto->Projection(0,i);
386
387     f1[i]->SetRange(40,60); // should be pion peak
388     h2D->FitSlicesY(f1[i],0,-1,10,"QNR",arr[i]); // gaus fit of pion peak
389
390     h1D = (TH1F*)arr[i]->At(1);
391     sprintf(name,"mean_dedx_mips_vs_%d",i);
392     h1D->SetName(name);
393     sprintf(title,"%s vs %s","mean_dedx_mips (a.u.)",fDeDxHisto->GetAxis(i)->GetTitle());
394     h1D->SetTitle(title);
395     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
396     h1D->GetYaxis()->SetTitle("mean_dedx_mips (a.u.)");
397     //h1D->SetMinimum(40);
398     //h1D->SetMaximum(60);
399
400     aFolderObj->Add(h1D);
401
402     h1D = (TH1F*)arr[i]->At(2);
403     sprintf(name,"res_dedx_mips_vs_%d",i);
404     h1D->SetName(name);
405     sprintf(title,"%s vs %s","res_dedx_mips (a.u)",fDeDxHisto->GetAxis(i)->GetTitle());
406     h1D->SetTitle(title);
407     h1D->GetXaxis()->SetTitle(fDeDxHisto->GetAxis(i)->GetTitle());
408     h1D->GetYaxis()->SetTitle("res_dedx_mips (a.u.)");
409     //h1D->SetMinimum(0);
410     //h1D->SetMaximum(6);
411
412     aFolderObj->Add(h1D);
413   }
414
415   // export objects to analysis folder
416   fAnalysisFolder = ExportToFolder(aFolderObj);
417
418   // delete only TObjrArray
419   if(aFolderObj) delete aFolderObj;
420 }
421
422 //_____________________________________________________________________________
423 TFolder* AliPerformanceDEdx::ExportToFolder(TObjArray * array) 
424 {
425   // recreate folder avery time and export objects to new one
426   //
427   AliPerformanceDEdx * comp=this;
428   TFolder *folder = comp->GetAnalysisFolder();
429
430   TString name, title;
431   TFolder *newFolder = 0;
432   Int_t i = 0;
433   Int_t size = array->GetSize();
434
435   if(folder) { 
436      // get name and title from old folder
437      name = folder->GetName();  
438      title = folder->GetTitle();  
439
440          // delete old one
441      delete folder;
442
443          // create new one
444      newFolder = CreateFolder(name.Data(),title.Data());
445      newFolder->SetOwner();
446
447          // add objects to folder
448      while(i < size) {
449            newFolder->Add(array->At(i));
450            i++;
451          }
452   }
453
454 return newFolder;
455 }
456
457
458 //_____________________________________________________________________________
459 TFolder* AliPerformanceDEdx::CreateFolder(TString name,TString title) { 
460 // create folder for analysed histograms
461 TFolder *folder = 0;
462   folder = new TFolder(name.Data(),title.Data());
463
464   return folder;
465 }
466