]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibQAChecker.cxx
Remove for the moment the QA analysis from the PartCorr wagon and move it to its...
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibQAChecker.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 ///////////////////////////////////////////////////////////////////////////////
17 //                                                                           //
18 //  Class to process a tree and create alarms based on thresholds            //
19 //  origin: jens wiechula: jens.wiechula@cern.ch                             //
20 //                                                                           //
21 ///////////////////////////////////////////////////////////////////////////////
22
23 #include <iostream>
24 #include <TObjArray.h>
25 #include <TString.h>
26 #include <TObjString.h>
27 #include <TStyle.h>
28 #include <TMarker.h>
29 #include <TAxis.h>
30 #include <TLine.h>
31 #include <TList.h>
32 #include <TTree.h>
33 #include <TMath.h>
34 #include <TGraph.h>
35 #include <TFrame.h>
36 #include <TIterator.h>
37 #include <TPad.h>
38 #include <TCanvas.h>
39 #include <TH1.h>
40 #include <TH2.h>
41 #include <TStopwatch.h>
42
43 #include <AliLog.h>
44
45 #include "AliTPCCalibQAChecker.h"
46
47 using namespace std;
48
49 AliTPCCalibQAChecker::AliTPCCalibQAChecker() :
50   TNamed("AliTPCCalibQAChecker","AliTPCCalibQAChecker"),
51   fTreePtr(0x0),
52   fHistPtr(0x0),
53   fGraphPtr(0x0),
54   fNumberPtr(0x0),
55   fHist(0x0),
56   fIterSubCheckers(0x0),
57   fArrSubCheckers(0x0),
58   fArrAlarmDescriptions(0x0),
59   fStrDrawRep(""),
60   fStrDrawRepOpt(""),
61   fStrDraw(""),
62   fStrDrawOpt(""),
63   fStrCuts(""),
64   fAlarmType(kMean),
65   fQualityLevel(kINFO),
66   fHistRep(0x0)
67 {
68   //
69   // Default ctor
70   //
71   ResetAlarmThresholds();
72 }
73 //_________________________________________________________________________
74 AliTPCCalibQAChecker::AliTPCCalibQAChecker(const char* name, const char *title) :
75   TNamed(name,title),
76   fTreePtr(0x0),
77   fHistPtr(0x0),
78   fGraphPtr(0x0),
79   fNumberPtr(0x0),
80   fHist(0x0),
81   fIterSubCheckers(0x0),
82   fArrSubCheckers(0x0),
83   fArrAlarmDescriptions(0x0),
84   fStrDrawRep(""),
85   fStrDrawRepOpt(""),
86   fStrDraw(""),
87   fStrDrawOpt(""),
88   fStrCuts(""),
89   fAlarmType(kMean),
90   fQualityLevel(kINFO),
91   fHistRep(0x0)
92 {
93   //
94   // TNamed ctor
95   //
96   ResetAlarmThresholds();
97 }
98 //_________________________________________________________________________
99 AliTPCCalibQAChecker::~AliTPCCalibQAChecker()
100 {
101   //
102   // Default ctor
103   //
104   if (fHistRep) delete fHistRep;
105   if (fIterSubCheckers) delete fIterSubCheckers;
106   if (fArrAlarmDescriptions) delete fArrAlarmDescriptions;
107 }
108 //_________________________________________________________________________
109 void AliTPCCalibQAChecker::AddSubChecker(AliTPCCalibQAChecker *alarm)
110 {
111   //
112   // add a sub checker to this checker
113   //
114   if (!alarm) return;
115   if (!fArrSubCheckers) {
116     fArrSubCheckers=new TObjArray;
117     fArrSubCheckers->SetOwner();
118   }
119   fArrSubCheckers->Add(alarm);
120 }
121 //_________________________________________________________________________
122 void AliTPCCalibQAChecker::Process()
123 {
124   //
125   // Process the alarm thresholds, decide the alarm level, create the representation histogram
126   //
127
128   //reset quality level
129   fQualityLevel=kINFO;
130
131   TStopwatch s;
132   s.Start();
133   //decide which type of checker to use
134   if (fArrSubCheckers && fArrSubCheckers->GetEntries()>0) ProcessSub();
135   else if (fTreePtr && *fTreePtr) ProcessTree();
136   else if (fHistPtr && *fHistPtr) ProcessHist();
137   else if (fGraphPtr && *fGraphPtr) ProcessGraph();
138   else if (fNumberPtr) ProcessNumber();
139   s.Stop();
140   AliInfo(Form("Processing Time (%s): %fs",GetName(),s.RealTime()));
141 }
142 //_________________________________________________________________________
143 void AliTPCCalibQAChecker::ProcessSub()
144 {
145   //
146   // sub checker type checker
147   //
148   QualityFlag_t quality=kINFO;
149   if (fArrSubCheckers && fArrSubCheckers->GetEntries()>0){
150     TIter next(fArrSubCheckers);
151     TObject *o=0x0;
152     while ( (o=next()) ) {
153       AliTPCCalibQAChecker *al=(AliTPCCalibQAChecker*)o;
154       al->Process();
155       QualityFlag_t subQuality=al->GetQuality();
156       if (subQuality>quality) quality=subQuality;
157     }
158   }
159   fQualityLevel=quality;
160 }
161 //_________________________________________________________________________
162 void AliTPCCalibQAChecker::ProcessTree()
163 {
164   //
165   // process tree type checker
166   //
167
168   //Create Representation Histogram
169   CreateRepresentationHist();
170   //
171 //   if (!fTree) return;
172   //chek for the quality
173
174   switch (fAlarmType){
175   case kNentries:
176     ProcessEntries();
177     break;
178   case kMean:
179   case kBinAny:
180   case kBinAll:
181     CreateAlarmHist();
182     ProcessHist();
183     ResetAlarmHist();
184     break;
185   }
186   
187 }
188 //_________________________________________________________________________
189 void AliTPCCalibQAChecker::ProcessHist()
190 {
191   //
192   // process histogram type checker
193   //
194   
195   if (!(fHistPtr && *fHistPtr)) return;
196     
197   switch (fAlarmType){
198   case kNentries:
199   case kMean:
200     ProcessMean();
201     break;
202   case kBinAny:
203   case kBinAll:
204     ProcessBin();
205     break;
206   }
207
208   //create representation histogram if we are not in tree mode
209   if (!fTreePtr){
210     if (!fHistRep) {
211       fHistRep=(*fHistPtr)->Clone();
212       TH1* h=(TH1*)fHistRep;
213       h->SetDirectory(0);
214       h->SetNameTitle(Form("h%s",GetName()),GetTitle());
215     }
216     TH1* h=(TH1*)fHistRep;
217     h->Reset();
218     h->Add(*fHistPtr);
219     if (!fHistRep->InheritsFrom(TH2::Class())) AddQualityLines(h);
220   }
221 }
222 //_________________________________________________________________________
223 void AliTPCCalibQAChecker::ProcessGraph()
224 {
225   //
226   // process graph type checker
227   //
228   if (!(fGraphPtr && *fGraphPtr)) return;
229   Int_t npoints=(*fGraphPtr)->GetN();
230   fQualityLevel=GetQuality(npoints,(*fGraphPtr)->GetY());
231 }
232 //_________________________________________________________________________
233 void AliTPCCalibQAChecker::ProcessNumber()
234 {
235   //
236   // process number type checker
237   //
238   if (!fNumberPtr) return;
239   fQualityLevel=GetQuality(*fNumberPtr);
240   //create representation histogram
241   if (!fHistRep){
242     TH1* h=new TH1F(Form("h%s",GetName()),GetTitle(),1,0,1);
243     h->GetXaxis()->SetBinLabel(1,"Value");
244     AddQualityLines(h);
245     h->SetDirectory(0);
246     fHistRep=h;
247   }
248   TH1 *h=(TH1*)fHistRep;
249   TMarker *marker=(TMarker*)h->GetListOfFunctions()->FindObject("TMarker");
250   if (!marker) {
251     marker=new TMarker(.5,0,20);
252     h->GetListOfFunctions()->Add(marker);
253   }
254   marker->SetMarkerColor(GetQualityColor());
255   marker->SetY(*fNumberPtr);
256 }
257 //_________________________________________________________________________
258 void AliTPCCalibQAChecker::ProcessEntries()
259 {
260   //
261   // Processing function which analyses the number of affected rows of a tree draw
262   //
263   TString draw=fStrDraw;
264   if (draw.IsNull()) return;
265   
266   TString cuts=fStrCuts;
267   
268   TString opt=fStrDrawOpt;
269   opt+="goff";
270   
271   //draw and get the histogram
272   Int_t res=(*fTreePtr)->Draw(draw.Data(),cuts.Data(),opt.Data());
273   fQualityLevel=GetQuality(res);
274 }
275 //_________________________________________________________________________
276 void AliTPCCalibQAChecker::ProcessMean()
277 {
278   //
279   // Processing function which analyses the mean of the resulting histogram
280   //
281
282   TH1* h=(*fHistPtr);
283   Double_t value=h->GetMean();
284   if (fAlarmType==kNentries) value=h->GetEntries();
285   fQualityLevel=GetQuality(value);
286 }
287 //_________________________________________________________________________
288 void AliTPCCalibQAChecker::ProcessBin()
289 {
290   //
291   // Process a histogram bin by bin and check for thresholds
292   //
293
294   //bin quality counters
295   Int_t nquality[kNQualityFlags];
296   for (Int_t iquality=(Int_t)kINFO; iquality<kNQualityFlags; ++iquality) nquality[iquality]=0;
297     
298   TH1 *h=(*fHistPtr);
299   
300   Int_t nbinsX=h->GetNbinsX();
301   Int_t nbinsY=h->GetNbinsY();
302   Int_t nbinsZ=h->GetNbinsZ();
303   Int_t nbinsTotal=nbinsX*nbinsY*nbinsZ;
304
305   //loop over all bins
306   for (Int_t ibinZ=1;ibinZ<nbinsZ+1;++ibinZ){
307     for (Int_t ibinY=1;ibinY<nbinsY+1;++ibinY){
308       for (Int_t ibinX=1;ibinX<nbinsX+1;++ibinX){
309         Double_t value = (*fHistPtr)->GetBinContent(ibinX, ibinY, ibinZ);
310         QualityFlag_t quality=GetQuality(value);
311         nquality[quality]++;
312       }
313     }
314   }
315
316   //loop over Quality levels and set quality
317   for (Int_t iquality=(Int_t)kINFO; iquality<kNQualityFlags; ++iquality){
318     if (fAlarmType==kBinAny){
319       if (nquality[iquality]) fQualityLevel=(QualityFlag_t)iquality;
320     } else if (fAlarmType==kBinAll){
321       if (nquality[iquality]==nbinsTotal) fQualityLevel=(QualityFlag_t)iquality;
322     }
323   }
324 }
325 //_________________________________________________________________________
326 void AliTPCCalibQAChecker::CreateRepresentationHist()
327 {
328   //
329   // Create the representation histogram which will be shown in the draw function
330   //
331   ResetRepresentationHist();
332
333   TString draw=fStrDrawRep;
334   if (draw.IsNull()) {
335     draw=fStrDraw;
336     fStrDrawRepOpt=fStrDrawOpt;
337   } else {
338     draw.ReplaceAll("%alarm%",fStrDraw.Data());
339   }
340   if (draw.IsNull()) return;
341   
342   TString cuts=fStrCuts;
343   
344   TString opt=fStrDrawRepOpt;
345   opt+="goff";
346   
347   Int_t res=(*fTreePtr)->Draw(draw.Data(),cuts.Data(),opt.Data());
348   TH1 *hist=(*fTreePtr)->GetHistogram();
349   if (res<0 || !hist){
350     AliError(Form("Could not create representation histogram of alarm '%s'",GetName()));
351     return;
352   }
353   fHistRep=hist->Clone();
354   TH1 *h=(TH1*)fHistRep;
355   h->SetDirectory(0);
356   h->SetNameTitle(Form("h%s",GetName()),GetTitle());
357 }
358 //_________________________________________________________________________
359 void AliTPCCalibQAChecker::CreateAlarmHist()
360 {
361   //
362   // create alarm histogram from the tree
363   //
364     
365   TString draw=fStrDraw;
366   if (draw.IsNull()) return;
367   
368   TString cuts=fStrCuts;
369   
370   TString opt=fStrDrawOpt;
371   opt+="goff";
372   
373   //draw and get the histogram
374   Int_t res=(*fTreePtr)->Draw(draw.Data(),cuts.Data(),opt.Data());
375   fHist=(*fTreePtr)->GetHistogram();
376   if (res<0 || !fHist){
377     AliError(Form("Could not create alarm histogram of alarm '%s'",GetName()));
378     return;
379   }
380   fHist->SetDirectory(0);
381   fHistPtr=&fHist;
382 }
383 //_________________________________________________________________________
384 void AliTPCCalibQAChecker::ResetAlarmHist()
385 {
386   //
387   // delete the alarm histogram and reset the pointer
388   //
389   if (fHistPtr){
390     if (*fHistPtr) delete *fHistPtr;
391     fHistPtr=0x0;
392   }
393 }
394 //_________________________________________________________________________
395 void AliTPCCalibQAChecker::Draw(Option_t *option)
396 {
397   //
398   // object draw function
399   // by default the pad backgound color is set to the quality level color
400   // use 'nobc' to change this
401   //
402
403   //case of a node with subcheckers and no specific representation histogram
404   if (HasSubCheckers()&&!fHistRep) {DrawSubNodes(option); return;}
405   if (fHistRep) {DrawRepresentationHist(option); return;}
406 }
407 //_________________________________________________________________________
408 void AliTPCCalibQAChecker::Print(Option_t *option) const
409 {
410   //
411   // print the quality status. If we have sub checkers print recursively
412   //
413   TString sOpt(option);
414   cout << sOpt << GetName() << ": " << GetQualityName() << endl;
415   if (fArrSubCheckers && fArrSubCheckers->GetEntries()>0){
416     sOpt.ReplaceAll("+-","  ");
417     sOpt+="+-";
418     TIter next(fArrSubCheckers);
419     TObject *o=0x0;
420     while ( (o=next()) ) o->Print(sOpt.Data());
421   }
422 }
423 //_________________________________________________________________________
424 void AliTPCCalibQAChecker::SetAlarmThreshold(const Double_t min, const Double_t max, const QualityFlag_t quality)
425 {
426   //
427   //set the alarm thresholds for a specific quality level
428   //
429   if ((Int_t)quality<(Int_t)kINFO||(Int_t)quality>=kNQualityFlags) return;
430   fThresMin[quality]=min;
431   fThresMax[quality]=max;
432 }
433 //_________________________________________________________________________
434 void AliTPCCalibQAChecker::ResetAlarmThreshold(const QualityFlag_t quality)
435 {
436   //
437   //set the alarm thresholds for a specific quality level
438   //
439   if ((Int_t)quality<(Int_t)kINFO||(Int_t)quality>=kNQualityFlags) return;
440   fThresMin[quality]=0;
441   fThresMax[quality]=0;
442 }
443 //_________________________________________________________________________
444 void AliTPCCalibQAChecker::ResetAlarmThresholds()
445 {
446   //
447   //reset all the alarm thresholds
448   //
449   for (Int_t i=0;i<kNQualityFlags;++i){
450     fThresMin[i]=0;
451     fThresMax[i]=0;
452   }
453 }
454 //_________________________________________________________________________
455 void AliTPCCalibQAChecker::SetQualityDescription(const char* text, const QualityFlag_t quality)
456 {
457   //
458   // set an description for the quality level
459   // %min and %max will be replaced by the min and max values of the alarm, when the quality
460   // description is queried (see QualityDescription)
461   //
462
463   if (quality<kINFO||quality>kFATAL) return;
464   if (! fArrAlarmDescriptions ) fArrAlarmDescriptions=new TObjArray(kNQualityFlags);
465   TObjString *s=(TObjString*)fArrAlarmDescriptions->At(quality);
466   if (!s) fArrAlarmDescriptions->AddAt(s=new TObjString,quality);
467   s->SetString(text);
468 }
469
470 //_________________________________________________________________________
471 const AliTPCCalibQAChecker* AliTPCCalibQAChecker::GetSubChecker(const char* name, Bool_t recursive) const
472 {
473   //
474   //
475   //
476   TString sname(name);
477   if (sname==GetName()) return this;
478   if (!fArrSubCheckers || !fArrSubCheckers->GetEntries()) return 0x0;
479   const AliTPCCalibQAChecker *al=0x0;
480   if (recursive){
481     TIter next(fArrSubCheckers);
482     TObject *o=0x0;
483     while ( (o=next()) ){
484       AliTPCCalibQAChecker *sal=(AliTPCCalibQAChecker*)o;
485       al=sal->GetSubChecker(name);
486       if (al) break;
487     }
488   }else{
489     al=dynamic_cast<AliTPCCalibQAChecker*>(fArrSubCheckers->FindObject(name));
490   }
491   return al;
492 }
493 //_________________________________________________________________________
494 Int_t AliTPCCalibQAChecker::GetNumberOfSubCheckers(Bool_t recursive) const
495 {
496   //
497   // get the number of sub checkers
498   // if recursive get total number of non subchecker type sub checkers
499   //
500   Int_t nsub=0;
501   if (recursive){
502     if (!fArrSubCheckers) return 1;
503     if (!fArrSubCheckers->GetEntries()) return 0;
504     TIter next(fArrSubCheckers);
505     TObject *o=0x0;
506     while ( (o=next()) ){
507       AliTPCCalibQAChecker *al=(AliTPCCalibQAChecker*)o;
508       nsub+=al->GetNumberOfSubCheckers();
509     }
510   } else {
511     if (fArrSubCheckers) nsub=fArrSubCheckers->GetEntries();
512   }
513   return nsub;
514 }
515 //_________________________________________________________________________
516 AliTPCCalibQAChecker* AliTPCCalibQAChecker::NextSubChecker()
517 {
518   //
519   // loop over sub checkers
520   // if recursive, recursively return the pointers of non subchecker type sub checkers
521   //
522   if (!fArrSubCheckers && !fArrSubCheckers->GetEntries()) return 0;
523   if (!fIterSubCheckers) fIterSubCheckers=fArrSubCheckers->MakeIterator();
524   AliTPCCalibQAChecker *al=(AliTPCCalibQAChecker*)fIterSubCheckers->Next();
525   if (!al){
526     delete fIterSubCheckers;
527     fIterSubCheckers=0x0;
528   }
529 //   if (recursive && al->GetNumberOfSubCheckers(kFALSE)) al=al->NextSubChecker();
530   return al;
531 }
532 //_________________________________________________________________________
533 const char* AliTPCCalibQAChecker::QualityName(const AliTPCCalibQAChecker::QualityFlag_t quality)
534 {
535   //
536   // get quality name for quality
537   //
538   switch (quality){
539   case kINFO:
540     return "Info";
541     break;
542   case kWARNING:
543     return "Warning";
544     break;
545   case kERROR:
546     return "Error";
547     break;
548   case kFATAL:
549     return "Fatal";
550     break;
551   default:
552     return "";
553   }
554 }
555 //_________________________________________________________________________
556 Color_t AliTPCCalibQAChecker::QualityColor(const AliTPCCalibQAChecker::QualityFlag_t quality)
557 {
558   //
559   // get quality color for quality
560   //
561   Color_t info = kSpring-8;
562   Color_t warning = kOrange;
563   Color_t error = kRed;
564   Color_t fatal = kRed+2;
565   Color_t none = kWhite;
566   
567   switch(quality) {
568   case kINFO :
569     return info;
570     break;
571   case kWARNING :
572     return warning;
573     break;
574   case kERROR :
575     return error;
576     break;
577   case kFATAL :
578     return fatal;
579     break;
580   default:
581     return none;
582   }
583   return none;
584   
585 }
586 //_________________________________________________________________________
587 const char* AliTPCCalibQAChecker::QualityDescription(const QualityFlag_t quality) const
588 {
589   //
590   // return description for quality
591   //
592   if (!fArrAlarmDescriptions || !fArrAlarmDescriptions->At(quality)) return "";
593   TString s(fArrAlarmDescriptions->At(quality)->GetName());
594   TString min, max;
595   min+=fThresMin[quality];
596   max+=fThresMax[quality];
597   s.ReplaceAll("%min",min);
598   s.ReplaceAll("%max",max);
599   return s.Data();
600 }
601 //_________________________________________________________________________
602 Int_t AliTPCCalibQAChecker::DrawInPad(TVirtualPad *pad, Int_t sub)
603 {
604   //
605   //
606   //
607   
608   if (fArrSubCheckers){
609     if (fArrSubCheckers->GetEntries()>0){
610       TIter next(fArrSubCheckers);
611       TObject *o=0x0;
612       while ( (o=next()) ) {
613         AliTPCCalibQAChecker *al=(AliTPCCalibQAChecker*)o;
614         sub=al->DrawInPad(pad,sub);
615       }
616     }
617   } else {
618     pad->cd(sub);
619     ++sub;
620     Draw();
621   }
622   return sub;
623 }
624 //_________________________________________________________________________
625 void AliTPCCalibQAChecker::DrawSubNodes(Option_t */*option*/)
626 {
627   //
628   // divide the current pad in sub pads and draw all sub nodes
629   //
630   if (!gPad) new TCanvas;
631   TPad *mother=(TPad*)gPad;
632   mother->Clear();
633   //find number of sub pads needed
634   Int_t nPads = GetNumberOfSubCheckers();
635   Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
636   Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
637   gPad->Divide(nCols,nRows);
638   
639   DrawInPad(gPad);
640   mother->Update();
641   TPad *pad=0;
642   Int_t ipad=1;
643   while ( (pad=(TPad*)mother->GetPad(ipad)) ){
644     TFrame* frame=(TFrame*)pad->GetPrimitive("TFrame");
645     if (frame) frame->SetFillColor(kWhite);
646     else {printf("no frame\n");}
647     pad->Modified();
648     ++ipad;
649   }
650   mother->Update();
651 }
652 //_________________________________________________________________________
653 void AliTPCCalibQAChecker::DrawRepresentationHist(Option_t *option)
654 {
655   //
656   // Draw the representation histogram
657   //
658   if (!fHistRep) return;
659   if (!gPad) new TCanvas;
660   Bool_t withBackColor=kTRUE;
661   
662   TString opt=option;
663   opt.ToLower();
664   
665   if (opt.Contains("nobc")) withBackColor=kFALSE;
666   opt.ReplaceAll("nobc","");
667   
668   if (opt.IsNull()) opt=fStrDrawRepOpt;
669   opt.ToLower();
670   
671   opt.ReplaceAll("prof","");
672   
673   fHistRep->Draw(opt.Data());
674   if (withBackColor) {
675     gPad->SetFillColor(GetQualityColor());
676   }
677   
678   gPad->Modified();
679 }
680 //_________________________________________________________________________
681 void AliTPCCalibQAChecker::AddQualityLines(TH1 *hist)
682 {
683   //
684   // add lines indicating the quality level to hist
685   //
686
687   Double_t xmin=hist->GetXaxis()->GetXmin();
688   Double_t xmax=hist->GetXaxis()->GetXmax();
689   Double_t min=0;
690   Double_t max=0;
691   
692   for (Int_t i=(Int_t)kINFO; i<kNQualityFlags; ++i){
693     if (fThresMin[i]>=fThresMax[i]) continue;
694     min=fThresMin[i];
695     max=fThresMax[i];
696     TLine *lmin=new TLine(xmin,min,xmax,min);
697     lmin->SetLineWidth(2);
698     lmin->SetLineColor(QualityColor((QualityFlag_t)i));
699     TLine *lmax=new TLine(xmin,max,xmax,max);
700     lmax->SetLineWidth(2);
701     lmax->SetLineColor(QualityColor((QualityFlag_t)i));
702     hist->GetListOfFunctions()->Add(lmin);
703     hist->GetListOfFunctions()->Add(lmax);
704   }
705   hist->SetAxisRange(min,max,"Y");
706 }