]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliTRDQATask.cxx
Reconstruction QA by Sylwester
[u/mrichter/AliRoot.git] / ESDCheck / AliTRDQATask.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // An analysis task to check the TRD data in simulated data
20 //
21 //*-- Sylwester Radomski
22 //////////////////////////////////////////////////////////////////////////////
23 // track type codding
24 //
25 // tpci = kTPCin
26 // tpco = kTPCout
27 // tpcz = kTPCout && !kTRDout
28 // trdo = kTRDout
29 // trdr = kTRDref
30 // trdz = kTRDout && !kTRDref
31 // 
32
33 #include <TCanvas.h>
34 #include <TChain.h>
35 #include <TFile.h>
36 #include <TGaxis.h>
37 #include <TH1D.h>
38 #include <TH2D.h>
39 #include <TROOT.h>
40 #include <TStyle.h>
41
42 #include "AliTRDQATask.h"
43 #include "AliESD.h"
44 #include "AliLog.h"
45
46 //______________________________________________________________________________
47 AliTRDQATask::AliTRDQATask(const char *name) : 
48   AliAnalysisTask(name,""),  
49   fChain(0),
50   fESD(0)
51 {
52   // Constructor.
53   // Input slot #0 works with an Ntuple
54   DefineInput(0, TChain::Class());
55   // Output slot #0 writes into a TH1 container
56   DefineOutput(0,  TObjArray::Class()) ; 
57 }
58
59 //______________________________________________________________________________
60 void AliTRDQATask::ConnectInputData(const Option_t *)
61 {
62   // Initialisation of branch container and histograms 
63
64   AliInfo(Form("*** Initialization of %s", GetName())) ; 
65   
66   // Get input data
67   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
68   if (!fChain) {
69     AliError(Form("Input 0 for %s not found\n", GetName()));
70     return ;
71   }
72
73   // One should first check if the branch address was taken by some other task
74   char ** address = (char **)GetBranchAddress(0, "ESD");
75   if (address) {
76     fESD = (AliESD*)(*address);
77   } else {
78     fESD = new AliESD();
79     SetBranchAddress(0, "ESD", &fESD);
80   }
81 }
82
83 //________________________________________________________________________
84 void AliTRDQATask::CreateOutputObjects()
85 {
86   // create histograms 
87   fNTracks     = new TH1D("ntracks", ";number of all tracks", 500, -0.5, 499.5); 
88   fEventSize   = new TH1D("evSize", ";event size (MB)", 100, 0, 5);
89
90   fTrackStatus = new TH1D("trackStatus", ";status bit", 32, -0.5, 31.5);
91   fKinkIndex   = new TH1D("kinkIndex", ";kink index", 41, -20.5, 20.5);
92   
93   fParIn  = new TH1D("parIn", "Inner Plane", 2, -0.5, 1.5);
94   fParOut = new TH1D("parOut", "Outer Plane", 2, -0.5, 1.5);
95
96   fXIn  = new TH1D("xIn", ";X at the inner plane (cm)", 200, 50, 250);
97   fXOut = new TH1D("xOut", ";X at the outer plane (cm)", 300, 50, 400);
98   
99   const int knNameAlpha = 4;
100   const char *namesAlpha[knNameAlpha] = {"alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr"};
101   //TH1D *fAlpha[4];
102   for(int i=0; i<knNameAlpha; i++) {
103     fAlpha[i] = new TH1D(namesAlpha[i], "alpha", 100, -4, 4);
104   }
105   fSectorTRD = new TH1D ("sectorTRD", ";sector TRD", 20, -0.5, 19.5);
106
107
108   // track parameters
109   const int knbits = 6;
110   const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
111   for(int i=0; i<knbits; i++) {
112     fPt[i]      = new TH1D(Form("pt%s",suf[i]), ";p_{T} (GeV/c);entries TPC", 50, 0, 10);
113     fTheta[i]   = new TH1D(Form("theta%s", suf[i]), "theta (rad)", 100, -4, 4); 
114     fSigmaY[i]  = new TH1D(Form("sigmaY%s",suf[i]),  ";sigma Y (cm)", 200, 0, 1);
115     fChi2[i]    = new TH1D(Form("Chi2%s", suf[i]), ";#chi2 / ndf", 100, 0, 10);
116     fPlaneYZ[i] = new TH2D(Form("planeYZ%s", suf[i]), Form("%sy (cm);z (cm)", suf[i]), 
117                            100, -60, 60, 500, -500, 500);
118   }
119   
120   // efficiency
121   fEffPt[0] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[0], suf[1]));
122   fEffPt[1] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[1], suf[3]));
123   fEffPt[2] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[3], suf[4]));
124   fEffPt[3] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[1], suf[4]));
125
126   for(int i=0; i<4; i++) {
127     fEffPt[i]->Sumw2();
128     fEffPt[i]->SetMarkerStyle(20);
129     fEffPt[i]->SetMinimum(0);
130     fEffPt[i]->SetMaximum(1.1);
131   }
132
133   // track features
134   fClustersTRD[0] = new TH1D("clsTRDo", "TRDo;number of clusters", 130, -0.5, 129.5);;
135   fClustersTRD[1] = new TH1D("clsTRDr", "TRDr;number of clusters", 130, -0.5, 129.5);;
136   fClustersTRD[2] = new TH1D("clsTRDz", "TRDz;number of clusters", 130, -0.5, 129.5);;
137
138   // for good refitted tracks only
139   fTime    = new TH1D("time", ";time bin", 25, -0.5, 24.5);
140   fBudget  = new TH1D("budget", ";material budget", 100, 0, 100);
141   fQuality = new TH1D("quality", ";track quality", 100, 0, 1.1);
142   fSignal  = new TH1D("signal", ";signal", 100, 0, 1e3);  
143   
144   // dEdX and PID
145   fTrdSigMom = new TH2D("trdSigMom", ";momentum (GeV/c);signal", 100, 0, 3, 100, 0, 1e3);
146   fTpcSigMom = new TH2D("tpcSigMom", ";momentum (GeV/c);signal", 100, 0, 3, 100, 0, 200);
147   
148   const char *pidName[6] = {"El", "Mu", "Pi", "K", "P", "Ch"};
149   for(int i=0; i<6; i++) {
150     
151     // TPC
152     fTpcPID[i] = new TH1D(Form("tpcPid%s",pidName[i]), pidName[i], 100, 0, 1.5);
153     fTpcPID[i]->GetXaxis()->SetTitle("probability");
154     
155     fTpcSigMomPID[i] = new TH2D(Form("tpcSigMom%s",pidName[i]), "", 100, 0, 3, 100, 0, 200);
156     fTpcSigMomPID[i]->SetTitle(Form("%s;momentum (GeV/c);signal",pidName[i])); 
157     
158     // TRD
159     fTrdPID[i] = new TH1D(Form("trdPid%s",pidName[i]), pidName[i], 100, 0, 1.5);
160     fTrdPID[i]->GetXaxis()->SetTitle("probability");
161      
162     fTrdSigMomPID[i] = new TH2D(Form("trdSigMom%s",pidName[i]), "", 100, 0, 3, 100, 0, 1e3);
163     fTrdSigMomPID[i]->SetTitle(Form("%s;momentum (GeV/c);signal",pidName[i]));  
164   }
165
166   
167   // create output container
168   fOutputContainer = new TObjArray(150); 
169   
170   // register histograms to the container  
171   int counter = 0;
172   
173   fOutputContainer->AddAt(fNTracks,     counter++);
174   fOutputContainer->AddAt(fEventSize,   counter++);
175   fOutputContainer->AddAt(fTrackStatus, counter++);
176   fOutputContainer->AddAt(fKinkIndex,   counter++);
177   fOutputContainer->AddAt(fParIn,       counter++);
178   fOutputContainer->AddAt(fParOut,      counter++);
179   fOutputContainer->AddAt(fXIn,         counter++);
180   fOutputContainer->AddAt(fXOut,        counter++);
181   fOutputContainer->AddAt(fAlpha[0],    counter++);
182   fOutputContainer->AddAt(fAlpha[1],    counter++);
183   fOutputContainer->AddAt(fAlpha[2],    counter++);
184   fOutputContainer->AddAt(fAlpha[3],    counter++);
185
186   fOutputContainer->AddAt(fSectorTRD,   counter++);
187   for(int i=0; i<knbits; i++) {
188      fOutputContainer->AddAt(fPt[i],      counter++);
189      fOutputContainer->AddAt(fTheta[i],   counter++);
190      fOutputContainer->AddAt(fSigmaY[i],  counter++);
191      fOutputContainer->AddAt(fChi2[i],    counter++);
192      fOutputContainer->AddAt(fPlaneYZ[i], counter++);
193   }   
194   fOutputContainer->AddAt(fEffPt[0], counter++);
195   fOutputContainer->AddAt(fEffPt[1], counter++);
196   fOutputContainer->AddAt(fEffPt[2], counter++);
197   fOutputContainer->AddAt(fEffPt[3], counter++);
198
199   fOutputContainer->AddAt(fClustersTRD[0], counter++);
200   fOutputContainer->AddAt(fClustersTRD[1], counter++);
201   fOutputContainer->AddAt(fClustersTRD[2], counter++);
202   fOutputContainer->AddAt(fTime,      counter++);
203   fOutputContainer->AddAt(fBudget,    counter++);
204   fOutputContainer->AddAt(fQuality,   counter++);
205   fOutputContainer->AddAt(fSignal,    counter++);
206   fOutputContainer->AddAt(fTrdSigMom, counter++);
207   fOutputContainer->AddAt(fTpcSigMom, counter++);
208   for(int i=0; i<6; i++) {
209      fOutputContainer->AddAt(fTpcPID[i],       counter++);
210      fOutputContainer->AddAt(fTpcSigMomPID[i], counter++);
211      fOutputContainer->AddAt(fTrdPID[i],       counter++);
212      fOutputContainer->AddAt(fTrdSigMomPID[i], counter++);
213   }
214
215   //AliInfo(Form("Number of histograms = %d", counter));
216
217  }
218
219 //______________________________________________________________________________
220 void AliTRDQATask::Exec(Option_t *) 
221 {
222   // Process one event
223   
224    Long64_t entry = fChain->GetReadEntry() ;
225   
226   // Processing of one event 
227    
228   if (!fESD) {
229     AliError("fESD is not connected to the input!") ; 
230     return ; 
231   }
232   
233   if ( !((entry-1)%100) ) 
234     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
235
236   int nTracks = fESD->GetNumberOfTracks();
237   fNTracks->Fill(nTracks); 
238
239   // track loop
240   for(int i=0; i<nTracks; i++) {
241     
242     AliESDtrack *track = fESD->GetTrack(i);
243     const AliExternalTrackParam *paramOut = track->GetOuterParam();
244     const AliExternalTrackParam *paramIn = track->GetInnerParam();
245
246     fParIn->Fill(!!paramIn);
247     if (!paramIn) continue;
248     fXIn->Fill(paramIn->GetX());
249
250     fParOut->Fill(!!paramOut);
251     if (!paramOut) continue;
252     fXOut->Fill(paramOut->GetX());
253  
254     int sector = GetSector(paramOut->GetAlpha());
255     if (!CheckSector(sector)) continue;
256     fSectorTRD->Fill(sector);
257
258     fKinkIndex->Fill(track->GetKinkIndex(0));
259     if (track->GetKinkIndex(0)) continue;    
260
261     UInt_t u = 1;
262     UInt_t status = track->GetStatus();
263     for(int bit=0; bit<32; bit++) 
264       if (u<<bit & status) fTrackStatus->Fill(bit);
265
266     const int knbits = 6; 
267     int bit[6] = {0,0,0,0,0,0};    
268     bit[0] = status & AliESDtrack::kTPCin;
269     bit[1] = status & AliESDtrack::kTPCout;
270     bit[2] = (status & AliESDtrack::kTPCout) && !(status & AliESDtrack::kTRDout);
271     bit[3] = status & AliESDtrack::kTRDout;
272     bit[4] = status & AliESDtrack::kTRDrefit;
273     bit[5] = (status & AliESDtrack::kTRDout) && !(status & AliESDtrack::kTRDrefit);
274
275     
276     // transverse momentum
277     const double *val = track->GetParameter(); // parameters at the vertex
278     double pt = 1./TMath::Abs(val[4]);
279
280     for(int b=0; b<knbits; b++) {
281       if (bit[b]) {
282         fPt[b]->Fill(pt); 
283         fTheta[b]->Fill(val[3]);
284         fSigmaY[b]->Fill(TMath::Sqrt(paramOut->GetSigmaY2()));
285         fChi2[b]->Fill(track->GetTRDchi2()/track->GetTRDncls());    
286         fPlaneYZ[b]->Fill(paramOut->GetY(), paramOut->GetZ()); 
287       }
288     }
289
290     // sectors
291     if (bit[1]) {
292       fAlpha[0]->Fill(paramIn->GetAlpha());
293       fAlpha[1]->Fill(paramOut->GetAlpha());
294     }
295     
296     if (bit[3]) fAlpha[2]->Fill(paramOut->GetAlpha());
297     if (bit[4]) fAlpha[3]->Fill(paramOut->GetAlpha());
298
299     // clusters
300     for(int b=0; b<3; b++) 
301       if (bit[3+b]) fClustersTRD[b]->Fill(track->GetTRDncls());
302
303     // refitted only
304     if (!bit[4]) continue;
305
306     fQuality->Fill(track->GetTRDQuality());
307     fBudget->Fill(track->GetTRDBudget());
308     fSignal->Fill(track->GetTRDsignal());
309         
310     fTrdSigMom->Fill(track->GetP(), track->GetTRDsignal());
311     fTpcSigMom->Fill(track->GetP(), track->GetTPCsignal());
312
313     // PID only
314     if (status & AliESDtrack::kTRDpid) {
315       
316       for(int l=0; l<6; l++) fTime->Fill(track->GetTRDTimBin(l));
317
318       // fill pid histograms
319       double trdr0 = 0, tpcr0 = 0;
320       int trdBestPid = 5, tpcBestPid = 5;  // charged
321       const double kminPidValue =  0.9;
322
323       double pp[5];
324       track->GetTPCpid(pp); // ESD inconsequence
325
326       for(int pid=0; pid<5; pid++) {
327         
328         trdr0 += track->GetTRDpid(pid);
329         tpcr0 += pp[pid];
330         
331         fTrdPID[pid]->Fill(track->GetTRDpid(pid));
332         fTpcPID[pid]->Fill(pp[pid]);
333         
334         if (track->GetTRDpid(pid) > kminPidValue) trdBestPid = pid;
335         if (pp[pid] > kminPidValue) tpcBestPid = pid;
336       }
337       
338       fTrdPID[5]->Fill(trdr0); // check unitarity
339       fTrdSigMomPID[trdBestPid]->Fill(track->GetP(), track->GetTRDsignal());
340       
341       fTpcPID[5]->Fill(tpcr0); // check unitarity
342       fTpcSigMomPID[tpcBestPid]->Fill(track->GetP(), track->GetTPCsignal());
343     }
344     
345   }
346
347   CalculateEff();
348   PostData(0, fOutputContainer);
349 }
350
351 //______________________________________________________________________________
352 void AliTRDQATask::Terminate(Option_t *)
353 {
354   // Processing when the event loop is ended
355   fOutputContainer = (TObjArray*)GetOutputData(0);
356   int counter = 0;
357   fNTracks     = (TH1D*)fOutputContainer->At(counter++);
358   fEventSize   = (TH1D*)fOutputContainer->At(counter++);
359   fTrackStatus = (TH1D*)fOutputContainer->At(counter++);
360   fKinkIndex   = (TH1D*)fOutputContainer->At(counter++);
361   fParIn       = (TH1D*)fOutputContainer->At(counter++);
362   fParOut      = (TH1D*)fOutputContainer->At(counter++);
363   fXIn         = (TH1D*)fOutputContainer->At(counter++);
364   fXOut        = (TH1D*)fOutputContainer->At(counter++);
365   fAlpha[0]    = (TH1D*)fOutputContainer->At(counter++);
366   fAlpha[1]    = (TH1D*)fOutputContainer->At(counter++);
367   fAlpha[2]    = (TH1D*)fOutputContainer->At(counter++);
368   fAlpha[3]    = (TH1D*)fOutputContainer->At(counter++);
369
370   fSectorTRD   = (TH1D*)fOutputContainer->At(counter++);
371   const int knbits = 6;
372   for(int i=0; i<knbits; i++) {
373      fPt[i]      = (TH1D*)fOutputContainer->At(counter++);
374      fTheta[i]   = (TH1D*)fOutputContainer->At(counter++);
375      fSigmaY[i]  = (TH1D*)fOutputContainer->At(counter++);
376      fChi2[i]    = (TH1D*)fOutputContainer->At(counter++);
377      fPlaneYZ[i] = (TH2D*)fOutputContainer->At(counter++);
378   }   
379   fEffPt[0] = (TH1D*)fOutputContainer->At(counter++);
380   fEffPt[1] = (TH1D*)fOutputContainer->At(counter++);
381   fEffPt[2] = (TH1D*)fOutputContainer->At(counter++);
382   fEffPt[3] = (TH1D*)fOutputContainer->At(counter++);
383
384   fClustersTRD[0] = (TH1D*)fOutputContainer->At(counter++);
385   fClustersTRD[1] = (TH1D*)fOutputContainer->At(counter++);
386   fClustersTRD[2] = (TH1D*)fOutputContainer->At(counter++);
387   fTime      = (TH1D*)fOutputContainer->At(counter++);
388   fBudget    = (TH1D*)fOutputContainer->At(counter++);
389   fQuality   = (TH1D*)fOutputContainer->At(counter++);
390   fSignal    = (TH1D*)fOutputContainer->At(counter++);
391   fTrdSigMom = (TH2D*)fOutputContainer->At(counter++);
392   fTpcSigMom = (TH2D*)fOutputContainer->At(counter++);
393   for(int i=0; i<6; i++) {
394      fTpcPID[i]       = (TH1D*)fOutputContainer->At(counter++);
395      fTpcSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++);
396      fTrdPID[i]       = (TH1D*)fOutputContainer->At(counter++);
397      fTrdSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++);
398   }
399
400   // create efficiency histograms
401   Bool_t problem = kFALSE ; 
402   AliInfo(Form(" *** %s Report:", GetName())) ; 
403   
404   CalculateEff();
405
406   DrawESD() ; 
407   DrawGeoESD() ; 
408   //DrawConvESD() ; 
409   DrawPidESD() ; 
410
411   char line[1024] ; 
412   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
413   gROOT->ProcessLine(line);
414   
415   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
416
417   char * report ; 
418   if(problem)
419     report="Problems found, please check!!!";  
420   else 
421     report="OK";
422
423   AliInfo(Form("*** %s Summary Report: %s\n",GetName(), report)) ; 
424
425 }
426
427 //______________________________________________________________________________
428 const int AliTRDQATask::GetSector(const double alpha) const
429 {
430   // Gets the sector number 
431
432   double size = TMath::DegToRad() * 20.;
433   int sector = (int)((alpha + TMath::Pi())/size);
434   return sector;
435 }
436
437 //______________________________________________________________________________
438 const int AliTRDQATask::CheckSector(const int sector) const  
439 {  
440   // Checks the sector number
441   const int knSec = 8;
442   int sec[] = {2,3,5,6,11,12,13,15};
443   
444   for(int i=0; i<knSec; i++) 
445     if (sector == sec[i]) return 1;
446   
447   return 0;
448 }
449
450 //______________________________________________________________________________
451 void AliTRDQATask::CalculateEff() 
452 {
453   // calculates the efficiency
454   
455   for(int i=0; i<4; i++) fEffPt[i]->Reset();
456   
457   fEffPt[0]->Add(fPt[1]);
458   fEffPt[0]->Divide(fPt[0]);
459   
460   fEffPt[1]->Add(fPt[3]);
461   fEffPt[1]->Divide(fPt[1]);
462   
463   fEffPt[2]->Add(fPt[4]);
464   fEffPt[2]->Divide(fPt[3]);
465   
466   fEffPt[3]->Add(fPt[4]);
467   fEffPt[3]->Divide(fPt[1]);
468 }
469
470 //______________________________________________________________________________
471 void AliTRDQATask::DrawESD() const 
472 {
473   // Makes a few plots
474
475   TCanvas * cTRD = new TCanvas("cTRD", "TRD ESD Test", 400, 10, 600, 700) ;
476   cTRD->Divide(6,3) ;
477
478   gROOT->SetStyle("Plain");
479   gStyle->SetPalette(1);
480   gStyle->SetOptStat(0);
481   
482   TGaxis::SetMaxDigits(3);
483   
484   gStyle->SetLabelFont(52, "XYZ");
485   gStyle->SetTitleFont(62, "XYZ");
486   gStyle->SetPadRightMargin(0.02);
487
488   // draw all 
489   
490   const int knplots = 18;
491   const int knover[knplots] = {1,1,1,4,1,1,1,1,1,1,2,1,1,3,1,1,1,1};
492   const int knnames = 24;
493   const char *names[knnames] = {
494     "ntracks", "kinkIndex", "trackStatus", 
495     "ptTPCi", "ptTPCo", "ptTRDo", "ptTRDr",  "ptTPCz", "ptTRDz",
496     "eff_TPCi_TPCo",  "eff_TPCo_TRDo", "eff_TRDo_TRDr",  "eff_TPCo_TRDr",
497     "clsTRDo", "clsTRDr", "clsTRDz", 
498     "alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr", "sectorTRD",
499     "time", "budget", "signal"
500   };
501   
502   const int klogy[knnames] = {
503     1,1,1,
504     1,1,1,
505     0,0,0,0,
506     1,1,
507     0,0,0,0,0,
508     0,1,1
509   };
510
511   int nhist=0;
512   for(int i=0; i<knplots; i++) {
513   cTRD->cd(i+1) ;
514   
515   //  new TCanvas(names[i], names[nhist], 500, 300);
516       
517     for(int j=0; j<knover[i]; j++) {
518       TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++]));
519       if (!hist) continue;
520       if (hist->GetMaximum() > 0 ) 
521         gPad->SetLogy(klogy[i]);
522       if (strstr(hist->GetName(), "eff")) {
523         hist->SetMarkerStyle(20);
524         hist->SetMinimum(0);
525         hist->SetMaximum(1.2);
526       }
527
528       if (!j) hist->Draw();
529       else hist->Draw("SAME");
530     }
531   }
532   cTRD->Print("TRD_ESD.eps");
533 }
534
535 //______________________________________________________________________________
536 void AliTRDQATask::DrawGeoESD() const
537 {
538   // Makes a few plots
539
540   TCanvas * cTRDGeo = new TCanvas("cTRDGeo", "TRD ESDGeo Test", 400, 10, 600, 700) ;
541   cTRDGeo->Divide(4,2) ;
542
543   gStyle->SetOptStat(0);
544   TGaxis::SetMaxDigits(3);
545   
546   gStyle->SetLabelFont(52, "XYZ");
547   gStyle->SetTitleFont(62, "XYZ");
548   
549   gStyle->SetPadTopMargin(0.06);
550   gStyle->SetTitleBorderSize(0);
551   
552   // draw all   
553   const int knnames = 7;
554   const char *names[knnames] = {
555     "xIn", "xOut",
556     "planeYZTPCo", "planeYZTPCz", "planeYZTRDo", "planeYZTRDr", "planeYZTRDz",
557   };
558   
559   const char *opt[knnames] = {
560     "", "",
561     "colz","colz", "colz", "colz", "colz"
562   };
563   
564   const int klogy[knnames] = {
565     1,1,
566     0,0,0,0,0
567   };
568   
569   for(int i=0; i<knnames; i++) {
570   cTRDGeo->cd(i+1) ;
571     TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i]));
572     if (!hist) continue;
573     
574     //if (i<2) new TCanvas(names[i], names[i], 500, 300);
575     //else new TCanvas(names[i], names[i], 300, 900);
576    
577     if (hist->GetMaximum() > 0 ) 
578       gPad->SetLogy(klogy[i]);
579     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
580     
581     hist->Draw(opt[i]);    
582     AliInfo(Form("%s\t%d", names[i], hist->GetEntries()));
583   }
584   
585   cTRDGeo->Print("TRD_Geo.eps");
586 }
587
588 //______________________________________________________________________________
589 void AliTRDQATask::DrawConvESD() const
590 {
591   // Makes a few plots
592
593   AliInfo("Plotting....") ; 
594   TCanvas * cTRDConv = new TCanvas("cTRDConv", "TRD ESDConv Test", 400, 10, 600, 700) ;
595   cTRDConv->Divide(3,2) ;
596
597   gROOT->SetStyle("Plain");
598   gROOT->ForceStyle();
599   gStyle->SetPalette(1);
600   
601   TGaxis::SetMaxDigits(3);
602   
603   gStyle->SetLabelFont(52, "XYZ");
604   gStyle->SetTitleFont(62, "XYZ");
605   gStyle->SetPadRightMargin(0.02);
606
607   const int knnames = 9;
608   const int knplots = 5;
609   const int knover[knplots] = {3,1,1,3,1}; 
610   
611   const char *names[knnames] = {
612     "sigmaYTPCo","sigmaYTRDo", "sigmaYTRDr", "sigmaYTPCz", "sigmaYTRDz",
613     "Chi2TPCo", "Chi2TRDo", "Chi2TRDr", "Chi2TRDz"
614   };
615   
616   const char *opt[knplots] = {
617     "", "", "","","",
618   };
619   
620   const int klogy[knplots] = {
621     0,0,0,1,1
622   };
623
624   int nhist = 0;
625   for(int i=0; i<knplots; i++) {
626     cTRDConv->cd(i+1) ;
627     //new TCanvas(names[i], names[i], 500, 300);
628     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
629    
630     for(int j=0; j<knover[i]; j++) {
631       TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++]));
632       if ( hist->GetMaximum() > 0 ) 
633         gPad->SetLogy(klogy[i]);
634       if (!j) hist->Draw(opt[i]);
635       else hist->Draw("same");
636     }
637
638   }
639     cTRDConv->Print("TRD_Conv.eps");
640 }
641
642 //______________________________________________________________________________
643 void AliTRDQATask::DrawPidESD() const
644 {
645   // Makes a few plots
646
647   TCanvas * cTRDPid = new TCanvas("cTRDPid", "TRD ESDPid Test", 400, 10, 600, 700) ;
648   cTRDPid->Divide(9,3) ;
649
650   gROOT->SetStyle("Plain");
651   gROOT->ForceStyle();
652   gStyle->SetPalette(1);
653   gStyle->SetOptStat(0);
654   
655   TGaxis::SetMaxDigits(3);
656   
657   gStyle->SetLabelFont(52, "XYZ");
658   gStyle->SetTitleFont(62, "XYZ");
659
660   gStyle->SetPadTopMargin(0.05);
661   gStyle->SetPadRightMargin(0.02);
662
663   // draw all 
664  
665   const int knnames = 27;
666   const char *names[knnames] = {
667
668     "signal", "trdSigMom", "tpcSigMom",
669
670     "trdPidEl", "trdPidMu", "trdPidPi", "trdPidK", "trdPidP", "trdPidCh",
671     "trdSigMomEl", "trdSigMomMu", "trdSigMomPi", "trdSigMomK", "trdSigMomP", "trdSigMomCh",
672     
673     "tpcPidEl", "tpcPidMu", "tpcPidPi", "tpcPidK", "tpcPidP", "tpcPidCh",
674     "tpcSigMomEl", "tpcSigMomMu", "tpcSigMomPi", "tpcSigMomK", "tpcSigMomP", "tpcSigMomCh"
675
676   };
677   
678   const char *opt[knnames] = {
679
680     "", "colz", "colz",
681
682     "", "", "", "", "", "" ,
683     "colz", "colz", "colz", "colz", "colz", "colz",
684
685     "", "", "", "", "", "" ,
686     "colz", "colz", "colz", "colz", "colz", "colz" 
687   };
688   
689   const int klogy[knnames] = {
690
691     0,0,0,
692
693     1,1,1,1,1,1,
694     0,0,0,0,0,0,
695
696     1,1,1,1,1,1,
697     0,0,0,0,0,0    
698   };
699
700   for(int i=0; i<knnames; i++) {
701     cTRDPid->cd(i+1) ;
702
703     TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i]));
704     if (!hist) continue;
705     
706     //new TCanvas(names[i], names[i], 500, 300);
707     if ( hist->GetMaximum() > 0  ) 
708       gPad->SetLogy(klogy[i]);
709     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
710     
711     if (strstr(names[i],"sigMom")) gPad->SetLogz(1);
712     if (strstr(names[i],"SigMom")) gPad->SetLogz(1);
713
714     hist->Draw(opt[i]);    
715     AliInfo(Form("%s\t%d", names[i], hist->GetEntries()));
716   }
717    cTRDPid->Print("TRD_Pid.eps");
718 }