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