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