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