]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliTRDQATask.cxx
PROOF-aware version of the analysis framework (Andrei)
[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 nNameAlpha = 4;
100   const char *namesAlpha[nNameAlpha] = {"alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr"};
101   //TH1D *fAlpha[4];
102   for(int i=0; i<nNameAlpha; 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 nbits = 6;
110   const char *suf[nbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
111   for(int i=0; i<nbits; 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<nbits; 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 nbits = 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<nbits; 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 minPidValue =  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) > minPidValue) trdBestPid = pid;
335         if (pp[pid] > minPidValue) 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   AliInfo("TRD QA module");
356   fOutputContainer = (TObjArray*)GetOutputData(0);
357   int counter = 0;
358   fNTracks     = (TH1D*)fOutputContainer->At(counter++);
359   fEventSize   = (TH1D*)fOutputContainer->At(counter++);
360   fTrackStatus = (TH1D*)fOutputContainer->At(counter++);
361   fKinkIndex   = (TH1D*)fOutputContainer->At(counter++);
362   fParIn       = (TH1D*)fOutputContainer->At(counter++);
363   fParOut      = (TH1D*)fOutputContainer->At(counter++);
364   fXIn         = (TH1D*)fOutputContainer->At(counter++);
365   fXOut        = (TH1D*)fOutputContainer->At(counter++);
366   fAlpha[0]    = (TH1D*)fOutputContainer->At(counter++);
367   fAlpha[1]    = (TH1D*)fOutputContainer->At(counter++);
368   fAlpha[2]    = (TH1D*)fOutputContainer->At(counter++);
369   fAlpha[3]    = (TH1D*)fOutputContainer->At(counter++);
370
371   fSectorTRD   = (TH1D*)fOutputContainer->At(counter++);
372   const int nbits = 6;
373   for(int i=0; i<nbits; i++) {
374      fPt[i]      = (TH1D*)fOutputContainer->At(counter++);
375      fTheta[i]   = (TH1D*)fOutputContainer->At(counter++);
376      fSigmaY[i]  = (TH1D*)fOutputContainer->At(counter++);
377      fChi2[i]    = (TH1D*)fOutputContainer->At(counter++);
378      fPlaneYZ[i] = (TH2D*)fOutputContainer->At(counter++);
379   }   
380   fEffPt[0] = (TH1D*)fOutputContainer->At(counter++);
381   fEffPt[1] = (TH1D*)fOutputContainer->At(counter++);
382   fEffPt[2] = (TH1D*)fOutputContainer->At(counter++);
383   fEffPt[3] = (TH1D*)fOutputContainer->At(counter++);
384
385   fClustersTRD[0] = (TH1D*)fOutputContainer->At(counter++);
386   fClustersTRD[1] = (TH1D*)fOutputContainer->At(counter++);
387   fClustersTRD[2] = (TH1D*)fOutputContainer->At(counter++);
388   fTime      = (TH1D*)fOutputContainer->At(counter++);
389   fBudget    = (TH1D*)fOutputContainer->At(counter++);
390   fQuality   = (TH1D*)fOutputContainer->At(counter++);
391   fSignal    = (TH1D*)fOutputContainer->At(counter++);
392   fTrdSigMom = (TH2D*)fOutputContainer->At(counter++);
393   fTpcSigMom = (TH2D*)fOutputContainer->At(counter++);
394   for(int i=0; i<6; i++) {
395      fTpcPID[i]       = (TH1D*)fOutputContainer->At(counter++);
396      fTpcSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++);
397      fTrdPID[i]       = (TH1D*)fOutputContainer->At(counter++);
398      fTrdSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++);
399   }
400
401   // create efficiency histograms
402   
403   CalculateEff();
404 //  PostData(0, fOutputContainer);
405
406   DrawESD() ; 
407   DrawGeoESD() ; 
408   //DrawConvESD() ; 
409   DrawPidESD() ; 
410 }
411
412 //______________________________________________________________________________
413 int AliTRDQATask::GetSector(double alpha) 
414 {
415   // Gets the sector number 
416
417   double size = TMath::DegToRad() * 20.;
418   int sector = (int)((alpha + TMath::Pi())/size);
419   return sector;
420 }
421
422 //______________________________________________________________________________
423 int AliTRDQATask::CheckSector(int sector) 
424 {  
425   // Checks the sector number
426   const int nSec = 8;
427   int sec[] = {2,3,5,6,11,12,13,15};
428   
429   for(int i=0; i<nSec; i++) 
430     if (sector == sec[i]) return 1;
431   
432   return 0;
433 }
434
435 //______________________________________________________________________________
436 void AliTRDQATask::CalculateEff() 
437 {
438   // calculates the efficiency
439   
440   for(int i=0; i<4; i++) fEffPt[i]->Reset();
441   
442   fEffPt[0]->Add(fPt[1]);
443   fEffPt[0]->Divide(fPt[0]);
444   
445   fEffPt[1]->Add(fPt[3]);
446   fEffPt[1]->Divide(fPt[1]);
447   
448   fEffPt[2]->Add(fPt[4]);
449   fEffPt[2]->Divide(fPt[3]);
450   
451   fEffPt[3]->Add(fPt[4]);
452   fEffPt[3]->Divide(fPt[1]);
453 }
454
455 //______________________________________________________________________________
456 void AliTRDQATask::DrawESD() 
457 {
458   // Makes a few plots
459
460   AliInfo("Plotting....") ; 
461   
462   TCanvas * cTRD = new TCanvas("cTRD", "TRD ESD Test", 400, 10, 600, 700) ;
463   cTRD->Divide(6,3) ;
464
465   gROOT->SetStyle("Plain");
466   gStyle->SetPalette(1);
467   gStyle->SetOptStat(0);
468   
469   TGaxis::SetMaxDigits(3);
470   
471   gStyle->SetLabelFont(52, "XYZ");
472   gStyle->SetTitleFont(62, "XYZ");
473   gStyle->SetPadRightMargin(0.02);
474
475   // draw all 
476   
477   const int nplots = 18;
478   const int nover[nplots] = {1,1,1,4,1,1,1,1,1,1,2,1,1,3,1,1,1,1};
479   const int nnames = 24;
480   const char *names[nnames] = {
481     "ntracks", "kinkIndex", "trackStatus", 
482     "ptTPCi", "ptTPCo", "ptTRDo", "ptTRDr",  "ptTPCz", "ptTRDz",
483     "eff_TPCi_TPCo",  "eff_TPCo_TRDo", "eff_TRDo_TRDr",  "eff_TPCo_TRDr",
484     "clsTRDo", "clsTRDr", "clsTRDz", 
485     "alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr", "sectorTRD",
486     "time", "budget", "signal"
487   };
488   
489   const int logy[nnames] = {
490     1,1,1,
491     1,1,1,
492     0,0,0,0,
493     1,1,
494     0,0,0,0,0,
495     0,1,1
496   };
497
498   int nhist=0;
499   for(int i=0; i<nplots; i++) {
500   cTRD->cd(i+1) ;
501   
502   //  new TCanvas(names[i], names[nhist], 500, 300);
503     gPad->SetLogy(logy[i]);
504       
505     for(int j=0; j<nover[i]; j++) {
506       TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++]));
507       if (!hist) continue;
508       
509       if (strstr(hist->GetName(), "eff")) {
510         hist->SetMarkerStyle(20);
511         hist->SetMinimum(0);
512         hist->SetMaximum(1.2);
513       }
514
515       if (!j) hist->Draw();
516       else hist->Draw("SAME");
517     }
518   }
519   cTRD->Print("TRD_ESD.gif");
520 }
521
522 //______________________________________________________________________________
523 void AliTRDQATask::DrawGeoESD() 
524 {
525   // Makes a few plots
526
527   AliInfo("Plotting....") ; 
528
529   TCanvas * cTRDGeo = new TCanvas("cTRDGeo", "TRD ESDGeo Test", 400, 10, 600, 700) ;
530   cTRDGeo->Divide(4,2) ;
531
532   gStyle->SetOptStat(0);
533   TGaxis::SetMaxDigits(3);
534   
535   gStyle->SetLabelFont(52, "XYZ");
536   gStyle->SetTitleFont(62, "XYZ");
537   
538   gStyle->SetPadTopMargin(0.06);
539   gStyle->SetTitleBorderSize(0);
540   
541   // draw all   
542   const int nnames = 7;
543   const char *names[nnames] = {
544     "xIn", "xOut",
545     "planeYZTPCo", "planeYZTPCz", "planeYZTRDo", "planeYZTRDr", "planeYZTRDz",
546   };
547   
548   const char *opt[nnames] = {
549     "", "",
550     "colz","colz", "colz", "colz", "colz"
551   };
552   
553   const int logy[nnames] = {
554     1,1,
555     0,0,0,0,0
556   };
557   
558   for(int i=0; i<nnames; i++) {
559   cTRDGeo->cd(i+1) ;
560     TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i]));
561     if (!hist) continue;
562     
563     //if (i<2) new TCanvas(names[i], names[i], 500, 300);
564     //else new TCanvas(names[i], names[i], 300, 900);
565    
566     gPad->SetLogy(logy[i]);
567     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
568     
569     hist->Draw(opt[i]);    
570     AliInfo(Form("%s\t%d", names[i], hist->GetEntries()));
571   }
572   
573   cTRDGeo->Print("TRD_Geo.gif");
574 }
575
576 //______________________________________________________________________________
577 void AliTRDQATask::DrawConvESD() 
578 {
579   // Makes a few plots
580
581   AliInfo("Plotting....") ; 
582   TCanvas * cTRDConv = new TCanvas("cTRDConv", "TRD ESDConv Test", 400, 10, 600, 700) ;
583   cTRDConv->Divide(3,2) ;
584
585   gROOT->SetStyle("Plain");
586   gROOT->ForceStyle();
587   gStyle->SetPalette(1);
588   
589   TGaxis::SetMaxDigits(3);
590   
591   gStyle->SetLabelFont(52, "XYZ");
592   gStyle->SetTitleFont(62, "XYZ");
593   gStyle->SetPadRightMargin(0.02);
594
595   const int nnames = 9;
596   const int nplots = 5;
597   const int nover[nplots] = {3,1,1,3,1}; 
598   
599   const char *names[nnames] = {
600     "sigmaYTPCo","sigmaYTRDo", "sigmaYTRDr", "sigmaYTPCz", "sigmaYTRDz",
601     "Chi2TPCo", "Chi2TRDo", "Chi2TRDr", "Chi2TRDz"
602   };
603   
604   const char *opt[nplots] = {
605     "", "", "","","",
606   };
607   
608   const int logy[nplots] = {
609     0,0,0,1,1
610   };
611
612   int nhist = 0;
613   for(int i=0; i<nplots; i++) {
614     cTRDConv->cd(i+1) ;
615     //new TCanvas(names[i], names[i], 500, 300);
616     gPad->SetLogy(logy[i]);
617     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
618    
619     for(int j=0; j<nover[i]; j++) {
620       TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++]));
621       if (!j) hist->Draw(opt[i]);
622       else hist->Draw("same");
623     }
624
625   }
626     cTRDConv->Print("TRD_Conv.eps");
627 }
628
629 //______________________________________________________________________________
630 void AliTRDQATask::DrawPidESD() 
631 {
632   // Makes a few plots
633
634   AliInfo("Plotting....") ; 
635   TCanvas * cTRDPid = new TCanvas("cTRDPid", "TRD ESDPid Test", 400, 10, 600, 700) ;
636   cTRDPid->Divide(9,3) ;
637
638   gROOT->SetStyle("Plain");
639   gROOT->ForceStyle();
640   gStyle->SetPalette(1);
641   gStyle->SetOptStat(0);
642   
643   TGaxis::SetMaxDigits(3);
644   
645   gStyle->SetLabelFont(52, "XYZ");
646   gStyle->SetTitleFont(62, "XYZ");
647
648   gStyle->SetPadTopMargin(0.05);
649   gStyle->SetPadRightMargin(0.02);
650
651   // draw all 
652  
653   const int nnames = 27;
654   const char *names[nnames] = {
655
656     "signal", "trdSigMom", "tpcSigMom",
657
658     "trdPidEl", "trdPidMu", "trdPidPi", "trdPidK", "trdPidP", "trdPidCh",
659     "trdSigMomEl", "trdSigMomMu", "trdSigMomPi", "trdSigMomK", "trdSigMomP", "trdSigMomCh",
660     
661     "tpcPidEl", "tpcPidMu", "tpcPidPi", "tpcPidK", "tpcPidP", "tpcPidCh",
662     "tpcSigMomEl", "tpcSigMomMu", "tpcSigMomPi", "tpcSigMomK", "tpcSigMomP", "tpcSigMomCh"
663
664   };
665   
666   const char *opt[nnames] = {
667
668     "", "colz", "colz",
669
670     "", "", "", "", "", "" ,
671     "colz", "colz", "colz", "colz", "colz", "colz",
672
673     "", "", "", "", "", "" ,
674     "colz", "colz", "colz", "colz", "colz", "colz" 
675   };
676   
677   const int logy[nnames] = {
678
679     0,0,0,
680
681     1,1,1,1,1,1,
682     0,0,0,0,0,0,
683
684     1,1,1,1,1,1,
685     0,0,0,0,0,0    
686   };
687
688   for(int i=0; i<nnames; i++) {
689     cTRDPid->cd(i+1) ;
690
691     TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i]));
692     if (!hist) continue;
693     
694     //new TCanvas(names[i], names[i], 500, 300);
695     gPad->SetLogy(logy[i]);
696     if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1);
697     
698     if (strstr(names[i],"sigMom")) gPad->SetLogz(1);
699     if (strstr(names[i],"SigMom")) gPad->SetLogz(1);
700
701     hist->Draw(opt[i]);    
702     AliInfo(Form("%s\t%d", names[i], hist->GetEntries()));
703   }
704    cTRDPid->Print("TRD_Pid.gif");
705 }