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