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