1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // Produces the data needed to calculate the quality assurance. //
21 // All data must be mergeable objects. //
24 // Sylwester Radomski (radomski@physi.uni-heidelberg.de) //
26 ////////////////////////////////////////////////////////////////////////////
28 // --- ROOT system ---
29 #include <TClonesArray.h>
38 // --- AliRoot header files ---
39 #include "AliESDEvent.h"
41 #include "AliTRDcluster.h"
42 #include "AliTRDQADataMakerRec.h"
43 #include "AliTRDgeometry.h"
44 #include "AliTRDdataArrayI.h"
45 #include "AliTRDrawStreamTB.h"
47 #include "AliQAChecker.h"
49 ClassImp(AliTRDQADataMakerRec)
51 //____________________________________________________________________________
52 AliTRDQADataMakerRec::AliTRDQADataMakerRec() :
53 AliQADataMakerRec(AliQA::GetDetName(AliQA::kTRD), "TRD Quality Assurance Data Maker")
56 // Default constructor
59 //____________________________________________________________________________
60 AliTRDQADataMakerRec::AliTRDQADataMakerRec(const AliTRDQADataMakerRec& qadm) :
67 SetName((const char*)qadm.GetName()) ;
68 SetTitle((const char*)qadm.GetTitle());
72 //__________________________________________________________________
73 AliTRDQADataMakerRec& AliTRDQADataMakerRec::operator=(const AliTRDQADataMakerRec& qadm)
79 this->~AliTRDQADataMakerRec();
80 new(this) AliTRDQADataMakerRec(qadm);
85 //____________________________________________________________________________
86 void AliTRDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
89 // Detector specific actions at end of cycle
94 //AliInfo(Form("Fitting RecPoints %d", task))
96 if (task == AliQA::kRECPOINTS) {
98 TH1D *hist = new TH1D("fitHist", "", 200, -0.5, 199.5);
101 // fill detector map;
102 for(int i=0; i<540; i++) {
103 Double_t v = ((TH1D*)list->At(0))->GetBinContent(i+1);
107 TH2D *detMap = (TH2D*)list->At(51);
108 Int_t bin = detMap->FindBin(sm, det);
109 detMap->SetBinContent(bin, v);
113 // Rec points full chambers
114 for (Int_t i=0; i<540; i++) {
116 //AliInfo(Form("I = %d", i));
118 //TH1D *h = ((TH2D*)list->At(1))->ProjectionY(Form("qaTRD_recPoints_amp_%d",i), i+1, i+1);
120 for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
121 Double_t xvalue = hist->GetBinCenter(b);
122 Int_t bin = ((TH2D*)list->At(1))->FindBin(i,xvalue);
123 Double_t value = ((TH2D*)list->At(1))->GetBinContent(bin);
124 hist->SetBinContent(b, value);
127 //printf("Sum = %d %f\n", i, hist->GetSum());
128 if (hist->GetSum() < 100) continue; // chamber not present
130 hist->Fit("landau", "q0", "goff", 10, 180);
131 TF1 *fit = hist->GetFunction("landau");
132 ((TH1D*)list->At(12))->Fill(fit->GetParameter(1));
133 ((TH1D*)list->At(13))->Fill(fit->GetParameter(2));
136 // time-bin by time-bin sm by sm
137 for(Int_t i=0; i<18; i++) { // loop over super-modules
139 for(Int_t j=0; j<35; j++) { // loop over time bins
141 //TH1D *h = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, j+1, j+1);
143 for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
144 Double_t xvalue = hist->GetBinCenter(b);
147 for(Int_t det=i*30; det<(i+1)*30; det++) { // loop over detectors
148 Int_t bin = ((TH3D*)list->At(10))->FindBin(det,j,xvalue);
149 Double_t value = ((TH3D*)list->At(10))->GetBinContent(bin);
152 //printf("v = %f\n", value);
153 hist->SetBinContent(b, svalue);
156 if (hist->GetSum() < 100) continue;
157 //printf("fitting %d %d %f\n", i, j, hist->GetSum());
159 hist->Fit("landau", "q0", "goff", 10, 180);
160 TF1 *fit = hist->GetFunction("landau");
162 TH1D *h1 = (TH1D*)list->At(14+18+i);
163 Int_t bin = h1->FindBin(j);
164 // printf("%d %d %d\n", det, j, bin);
165 h1->SetBinContent(bin, TMath::Abs(fit->GetParameter(1)));
170 // time-bin by time-bin chamber by chamber
171 for (Int_t i=0; i<540; i++) {
173 //TH1D *test = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, 0, 35);
174 //if (test->GetSum() < 100) continue;
176 //AliInfo(Form("fitting det = %d", i));
178 for(Int_t j=0; j<35; j++) {
180 //TH1D *h = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, j+1, j+1);
182 for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
183 Double_t xvalue = hist->GetBinCenter(b);
184 Int_t bin = ((TH3D*)list->At(10))->FindBin(i,j,xvalue);
185 Double_t value = ((TH3D*)list->At(10))->GetBinContent(bin);
186 //printf("v = %f\n", value);
187 hist->SetBinContent(b, value);
190 if (hist->GetSum() < 100) continue;
191 //printf("fitting %d %d %f\n", i, j, hist->GetSum());
193 hist->Fit("landau", "q0", "goff", 10, 180);
194 TF1 *fit = hist->GetFunction("landau");
198 TH2D *h2 = (TH2D*)list->At(14+sm);
199 Int_t bin = h2->FindBin(det,j);
200 // printf("%d %d %d\n", det, j, bin);
201 h2->SetBinContent(bin, TMath::Abs(fit->GetParameter(1)));
202 h2->SetBinError(bin,fit->GetParError(1));
206 if (hist) delete hist;
209 //////////////////////////
210 // const Int_t knbits = 6;
211 // const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
212 //const char *sufRatio[4] = {"TRDrTRDo", "TRDoTPCo", "TRDrTPCo", "TRDzTPCo"};
214 if (task == AliQA::kESDS) {
216 const Int_t knRatio = 4;
217 const Int_t kN[knRatio] = {4,3,4,5};
218 const Int_t kD[knRatio] = {3,1,1,3};
221 for(Int_t type=0; type<2; type++) {
222 for(Int_t i=0; i<knRatio; i++) {
224 TH1D *ratio = (TH1D*)list->At(19 + 2*i + type);
225 TH1D *histN = (TH1D*)list->At(3 + 2*kN[i] + type);
226 TH1D *histD = (TH1D*)list->At(3 + 2*kD[i] + type);
228 BuildRatio(ratio, histN, histD);
231 //ratio->Divide(histD);
238 AliQAChecker::Instance()->Run(AliQA::kTRD, task, list) ;
241 //____________________________________________________________________________
242 void AliTRDQADataMakerRec::InitESDs()
245 // Create ESDs histograms in ESDs subdir
248 const Int_t kNhist = 27;
250 Int_t histoCounter = -1 ;
252 hist[++histoCounter] = new TH1D("qaTRD_esd_ntracks", ";Number of tracks", 300, -0.5, 299.5);
253 hist[++histoCounter] = new TH1D("qaTRD_esd_sector", ";Sector", 18, -0.5, 17.7);
254 hist[++histoCounter] = new TH1D("qaTRD_esd_bits", ";Bits", 64, -0.5, 63.5);
256 const Int_t knbits = 6;
257 const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
260 for(Int_t i=0; i<knbits; i++) {
261 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_pt%s",suf[i]), ";p_{T} (GeV/c);", 100, 0, 10);
262 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_trdz%s", suf[i]), ";z (cm)", 200, -400, 400);
265 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDo", "TRDo;number of clusters", 130, -0.5, 129.5);;
266 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDr", "TRDr;number of clusters", 130, -0.5, 129.5);;
267 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDz", "TRDz;number of clusters", 130, -0.5, 129.5);;
268 //hist[++histoCounter] = new TH1D("qaTRD_esd_clsRatio", ";cluster ratio", 100, 0., 1.3);;
270 hist[++histoCounter] = new TH2D("qaTRD_esd_sigMom", ";momentum (GeV/c);signal", 100, 0, 5, 200, 0, 1e3);
272 // end of cycle plots (hist 19)
273 const char *sufRatio[4] = {"TRDrTRDo", "TRDoTPCo", "TRDrTPCo", "TRDzTPCo"};
275 for(int i=0; i<4; i++) {
276 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_pt%s",sufRatio[i]),
277 Form("Efficiency in Pt %s;p_{T};eff", sufRatio[i]),
280 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_trdz%s",sufRatio[i]),
281 Form("Efficiency in Z %s;z (cm);eff", sufRatio[i]),
285 for(Int_t i=0; i<=histoCounter; i++) {
287 Add2ESDsList(hist[i], i);
292 //____________________________________________________________________________
293 void AliTRDQADataMakerRec::InitRecPoints()
296 // Create Reconstructed Points histograms in RecPoints subdir
299 const Int_t kNhist = 14 + 18 + 18 + 2;
302 hist[0] = new TH1D("qaTRD_recPoints_det", ";Detector ID of the cluster", 540, -0.5, 539.5);
303 hist[1] = new TH2D("qaTRD_recPoints_amp", ";Amplitude", 540, -0.5, 539, 200, -0.5, 199.5);
304 hist[2] = new TH1D("qaTRD_recPoints_npad", ";Number of Pads", 12, -0.5, 11.5);
306 hist[3] = new TH1D("qaTRD_recPoints_dist2", ";residuals [2pad]", 100, -1, 1);
307 hist[4] = new TH1D("qaTRD_recPoints_dist3", ";residuals [3pad]", 100, -1, 1);
308 hist[5] = new TH1D("qaTRD_recPoints_dist4", ";residuals [4pad]", 100, -1, 1);
309 hist[6] = new TH1D("qaTRD_recPoints_dist5", ";residuals [5pad]", 100, -1, 1);
311 hist[7] = new TH2D("qaTRD_recPoints_rowCol", ";row;col", 16, -0.5, 15.5, 145, -0.5, 144.5);
312 hist[8] = new TH1D("qaTRD_recPoints_time", ";time bin", 35, -0.5, 34.5);
313 hist[9] = new TH1D("qaTRD_recPoints_nCls", ";number of clusters", 500, -0.5, 499.5);
315 hist[10] = new TH3D("qaTRD_recPoints_sigTime", ";chamber;time bin;signal",
316 540, -0.5, 539.5, 35, -0.5, 34.5, 100, 0, 200);
317 hist[11] = new TProfile("qaTRD_recPoints_prf", ";distance;center of gravity"
318 , 120, -0.6, 0.6, -1.2, 1.2, "");
320 hist[12] = new TH1D("qaTRD_recPoints_ampMPV", ";amplitude MPV", 200, 0, 200);
321 hist[13] = new TH1D("qaTRD_recPoints_ampSigma", ";amplitude Sigma", 200, 0, 200);
323 // chamber bu chamber
324 for(Int_t i=0; i<18; i++) {
325 hist[14+i] = new TH2D(Form("qaTRD_recPoints_sigTime_sm%d",i), Form("sm%d;det;time bin"),
326 30, -0.5, 29.5, 35, -0.5, 34.5);
327 hist[14+i]->SetMinimum(20);
328 hist[14+i]->SetMaximum(40);
331 // time bin by time bin sm-by-sm
332 for(Int_t i=0; i<18; i++) {
333 hist[14+18+i] = new TH1D(Form("qaTRD_recPoints_sigTimeShape_sm%d", i),
334 Form("sm%d;time bin;signal"),
337 hist[14+18+i]->SetMaximum(120);
340 hist[50] = new TH1D("qaTRD_recPoints_signal", ";amplitude", 200, -0.5, 199.5);
341 hist[51] = new TH2D("qaTRD_recPoints_detMap", ";sm;chamber", 18, -0.5, 17.5, 30, -0.5, 29.5);
344 for(Int_t i=0; i<kNhist; i++) {
346 Add2RecPointsList(hist[i], i);
351 //____________________________________________________________________________
352 void AliTRDQADataMakerRec::InitRaws()
355 // create Raws histograms in Raws subdir
358 const Int_t kSM = 18;
359 //const Int_t kNCh = 540;
360 const Int_t kNhist = 4+kSM;
363 // four histograms to be published
364 hist[0] = new TH1D("qaTRD_raws_det", ";detector", 540, -0.5, 539.5);
365 hist[1] = new TH1D("qaTRD_raws_sig", ";signal", 100, -0.5, 99.5);
366 hist[2] = new TH1D("qaTRD_raws_timeBin", ";time bin", 40, -0.5, 39.5);
367 hist[3] = new TH1D("qaTRD_raws_smId", ";supermodule", 18, -0.5, 17.5);
370 // one double per MCM (not published)
371 const Int_t kNMCM = 30 * 8 * 16;
372 for(Int_t i=0; i<kSM; i++)
373 hist[4+i] = new TH1D(Form("qaTRD_raws_sm%d",i),"",kNMCM, -0.5, kNMCM-0.5);
376 for(Int_t i=0; i<kNhist; i++) {
378 Add2RawsList(hist[i], i);
383 //____________________________________________________________________________
384 void AliTRDQADataMakerRec::MakeESDs(AliESDEvent * esd)
387 // Make QA data from ESDs
390 Int_t nTracks = esd->GetNumberOfTracks();
391 GetESDsData(0)->Fill(nTracks);
394 for (Int_t i=0; i<nTracks; i++) {
396 AliESDtrack *track = esd->GetTrack(i);
397 const AliExternalTrackParam *paramOut = track->GetOuterParam();
398 const AliExternalTrackParam *paramIn = track->GetInnerParam();
401 if (!paramIn) continue;
402 if (!paramOut) continue;
405 if (track->GetKinkIndex(0) > 0) continue;
407 Double_t extZ = GetExtZ(paramIn);
408 if (TMath::Abs(extZ) > 320) continue; // acceptance cut
410 // .. in the acceptance
411 Int_t sector = GetSector(paramOut->GetAlpha());
412 GetESDsData(1)->Fill(sector);
415 UInt_t status = track->GetStatus();
416 for(Int_t bit=0; bit<32; bit++)
417 if (u<<bit & status) GetESDsData(2)->Fill(bit);
419 const Int_t knbits = 6;
420 Int_t bit[6] = {0,0,0,0,0,0};
421 bit[0] = status & AliESDtrack::kTPCin;
422 bit[1] = status & AliESDtrack::kTPCout;
423 bit[2] = (status & AliESDtrack::kTPCout) && !(status & AliESDtrack::kTRDout);
424 bit[3] = status & AliESDtrack::kTRDout;
425 bit[4] = status & AliESDtrack::kTRDrefit;
426 bit[5] = (status & AliESDtrack::kTRDout) && !(status & AliESDtrack::kTRDrefit);
428 // transverse momentum
429 //const Double_t *val = paramOut->GetParameter(); // parameters at the Outer plane
430 Double_t pt = paramOut->Pt(); //1./TMath::Abs(val[4]);
432 for(Int_t b=0; b<knbits; b++) {
434 GetESDsData(2*b+3)->Fill(pt);
435 GetESDsData(2*b+4)->Fill(extZ);
440 for(Int_t b=0; b<3; b++)
441 if (bit[3+b]) GetESDsData(b+15)->Fill(track->GetTRDncls());
444 if (!bit[4]) continue;
446 //fQuality->Fill(track->GetTRDQuality());
447 //fBudget->Fill(track->GetTRDBudget());
448 //fSignal->Fill(track->GetTRDsignal());
450 GetESDsData(18)->Fill(track->GetP(), track->GetTRDsignal());
454 if (status & AliESDtrack::kTRDpid) {
456 for(Int_t l=0; l<6; l++) fTime->Fill(track->GetTRDTimBin(l));
458 // fill pid histograms
459 Double_t trdr0 = 0; //, tpcr0 = 0;
460 Int_t trdBestPid = 5; //, tpcBestPid = 5; // charged
461 const Double_t kminPidValue = 0.9;
464 //track->GetTPCpid(pp); // ESD inconsequence
466 for(Int_t pid=0; pid<5; pid++) {
468 trdr0 += track->GetTRDpid(pid);
471 fTrdPID[pid]->Fill(track->GetTRDpid(pid));
472 //fTpcPID[pid]->Fill(pp[pid]);
474 if (track->GetTRDpid(pid) > kminPidValue) trdBestPid = pid;
475 //if (pp[pid] > kminPidValue) tpcBestPid = pid;
478 fTrdPID[5]->Fill(trdr0); // check unitarity
479 fTrdSigMomPID[trdBestPid]->Fill(track->GetP(), track->GetTRDsignal());
481 //fTpcPID[5]->Fill(tpcr0); // check unitarity
482 //fTpcSigMomPID[tpcBestPid]->Fill(track->GetP(), track->GetTPCsignal());
490 //______________________________________________________________________________
491 Int_t AliTRDQADataMakerRec::GetSector(const Double_t alpha) const
494 // Gets the sector number
497 Double_t size = TMath::DegToRad() * 20.; // shall use TRDgeo
498 Int_t sector = (Int_t)((alpha + TMath::Pi())/size);
503 //______________________________________________________________________________
504 Double_t AliTRDQADataMakerRec::GetExtZ(const AliExternalTrackParam *in) const
507 // Returns the Z position at the entry to TRD
508 // using parameters from the TPC in
511 const Double_t kX0 = 300;
513 Double_t x = in->GetX();
514 const Double_t *par = in->GetParameter();
515 Double_t theta = par[3];
516 Double_t z = in->GetZ();
518 Double_t zz = z + (kX0-x) * TMath::Tan(theta);
523 //____________________________________________________________________________
524 void AliTRDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
527 // Makes QA data from raw data
533 //const Int_t kSM = 18;
534 //const Int_t kROC = 30;
535 const Int_t kROB = 8;
536 //const Int_t kLayer = 6;
537 //const Int_t kStack = 5;
538 const Int_t kMCM = 16;
539 // const Int_t kADC = 22;
541 AliTRDrawStreamTB *raw = new AliTRDrawStreamTB(rawReader);
543 raw->SetRawVersion(3);
546 while (raw->Next()) {
548 GetRawsData(0)->Fill(raw->GetDet());
550 // possibly needs changes with the new reader !!
551 Int_t *sig = raw->GetSignals();
552 for(Int_t i=0; i<3; i++) GetRawsData(1)->Fill(sig[i]);
555 GetRawsData(2)->Fill(raw->GetTimeBin());
557 // calculate the index;
558 Int_t sm = raw->GetSM();
559 Int_t roc = raw->GetROC();
560 Int_t rob = raw->GetROB();
561 Int_t mcm = raw->GetMCM();
562 //Int_t adc = raw->GetADC();
564 //Int_t index = roc * (kROB*kMCM*kADC) + rob * (kMCM*kADC) + mcm * kADC + adc;
565 Int_t index = roc * (kROB*kMCM) + rob * kMCM + mcm;
566 GetRawsData(3)->Fill(sm);
567 GetRawsData(4+sm)->Fill(index);
571 //____________________________________________________________________________
572 void AliTRDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
575 // Makes data from RecPoints
578 // Info("MakeRecPoints", "making");
580 Int_t nsize = Int_t(clustersTree->GetTotBytes() / (sizeof(AliTRDcluster)));
581 TObjArray *clusterArray = new TObjArray(nsize+1000);
583 TBranch *branch = clustersTree->GetBranch("TRDcluster");
585 AliError("Can't get the branch !");
588 branch->SetAddress(&clusterArray);
590 // Loop through all entries in the tree
591 Int_t nEntries = (Int_t) clustersTree->GetEntries();
593 AliTRDcluster *c = 0;
595 for (Int_t i=0; i<540; i++) nDet[i] = 0;
597 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
600 nbytes += clustersTree->GetEvent(iEntry);
602 // Get the number of points in the detector
603 Int_t nCluster = clusterArray->GetEntriesFast();
605 // Loop through all TRD digits
606 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
607 c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster);
609 Int_t iDet = c->GetDetector();
611 GetRecPointsData(0)->Fill(iDet);
612 GetRecPointsData(50)->Fill(c->GetQ());
613 GetRecPointsData(1)->Fill(iDet, c->GetQ());
614 GetRecPointsData(2)->Fill(c->GetNPads());
615 if (c->GetNPads() < 6)
616 GetRecPointsData(1+c->GetNPads())->Fill(c->GetCenter());
618 //if (c->GetPadTime() < 5)
619 ((TH2D*)GetRecPointsData(7))->Fill(c->GetPadRow(), c->GetPadCol());
620 GetRecPointsData(8)->Fill(c->GetPadTime());
622 ((TH3D*)GetRecPointsData(10))->Fill(iDet, c->GetPadTime(), c->GetQ());
625 //if (c->GetNPads() == 2) {
626 Short_t *sig = c->GetSignals();
629 if (sig[0] == 0 && sig[1] == 0 && sig[2] == 0 && sig[5] == 0 && sig[6] == 0)
630 frac = 1. * sig[4] / (sig[3] + sig[4]);
632 if (sig[0] == 0 && sig[1] == 0 && sig[4] == 0 && sig[5] == 0 && sig[6] == 0)
633 frac = -1. * sig[2] / (sig[2] + sig[3]);
635 if (frac > -10) ((TProfile*)GetRecPointsData(11))->Fill(c->GetCenter(), frac);
641 for(Int_t i=0; i<540; i++)
642 if (nDet[i] > 0) GetRecPointsData(9)->Fill(nDet[i]);
648 //____________________________________________________________________________
649 void AliTRDQADataMakerRec::StartOfDetectorCycle()
652 // Detector specific actions at start of cycle
657 //__________________________________________________________________________
658 Int_t AliTRDQADataMakerRec::CheckPointer(TObject *obj, const char *name)
661 // Checks initialization of pointers
664 if (!obj) AliWarning(Form("null pointer: %s", name));
667 //__________________________________________________________________________
668 void AliTRDQADataMakerRec::BuildRatio(TH1D *ratio, TH1D *histN, TH1D*histD) {
670 // Calculate the ratio of two histograms
671 // error are calculated assuming the histos have the same counts
676 Int_t nbins = histN->GetXaxis()->GetNbins();
677 for(Int_t i=1; i<nbins+2; i++) {
679 Double_t valueN = histN->GetBinContent(i);
680 Double_t valueD = histD->GetBinContent(i);
683 ratio->SetBinContent(i, 0);
684 ratio->SetBinError(i, 0);
688 Double_t eps = (valueN < valueD-valueN)? valueN : valueD-valueN;
690 ratio->SetBinContent(i, valueN/valueD);
691 ratio->SetBinError(i, TMath::Sqrt(eps)/valueD);
695 ratio->SetMinimum(-0.1);
696 ratio->SetMaximum(1.1);
697 ratio->SetMarkerStyle(20);
699 //__________________________________________________________________________