]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/CDB/AliBaseCalibViewer.cxx
adding getters for the LHC period (OCDB folder) and its first/last run
[u/mrichter/AliRoot.git] / STEER / CDB / AliBaseCalibViewer.cxx
1 ///////////////////////////////////////////////////////////////////////////////
2 //                                                                           //
3 //  Base class for the AliTPCCalibViewer and AliTRDCalibViewer               //
4 //  used for the calibration monitor                                         //
5 //                                                                           //
6 //  Authors:     Marian Ivanov (Marian.Ivanov@cern.ch)                       //
7 //               Jens Wiechula (Jens.Wiechula@cern.ch)                       //
8 //               Ionut Arsene  (iarsene@cern.ch)                             //
9 //                                                                           //
10 ///////////////////////////////////////////////////////////////////////////////
11
12
13 #include <iostream>
14 #include <fstream>
15 #include <TString.h>
16 #include <TRandom.h>
17 #include <TLegend.h>
18 #include <TLine.h>
19 //#include <TCanvas.h>
20 #include <TROOT.h>
21 #include <TStyle.h>
22 #include <TH1.h> 
23 #include <TH1F.h>
24 #include <TMath.h>
25 #include <THashTable.h>
26 #include <TObjString.h>
27 #include <TLinearFitter.h>
28 #include <TFile.h>
29 #include <TKey.h>
30 #include <TGraph.h>
31 #include <TDirectory.h>
32 #include <TFriendElement.h>
33
34 #include "TTreeStream.h"
35 #include "AliBaseCalibViewer.h"
36
37 ClassImp(AliBaseCalibViewer)
38
39 AliBaseCalibViewer::AliBaseCalibViewer()
40                   :TObject(),
41                    fTree(0),
42                    fFile(0),
43                    fListOfObjectsToBeDeleted(0),
44                    fTreeMustBeDeleted(0), 
45                    fAbbreviation(0), 
46                    fAppendString(0)
47 {
48   //
49   // Default constructor
50   //
51 }
52
53 //_____________________________________________________________________________
54 AliBaseCalibViewer::AliBaseCalibViewer(const AliBaseCalibViewer &c)
55                   :TObject(c),
56                    fTree(0),
57                    fFile(0),
58                    fListOfObjectsToBeDeleted(0),
59                    fTreeMustBeDeleted(0),
60                    fAbbreviation(0), 
61                    fAppendString(0)
62 {
63   //
64   // dummy AliBaseCalibViewer copy constructor
65   // not yet working!!!
66   //
67   fTree = c.fTree;
68   fTreeMustBeDeleted = c.fTreeMustBeDeleted;
69   fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
70   fAbbreviation = c.fAbbreviation;
71   fAppendString = c.fAppendString;
72 }
73
74 //_____________________________________________________________________________
75 AliBaseCalibViewer::AliBaseCalibViewer(TTree* tree)
76                   :TObject(),
77                    fTree(0),
78                    fFile(0),
79                    fListOfObjectsToBeDeleted(0),
80                    fTreeMustBeDeleted(0),
81                    fAbbreviation(0), 
82                    fAppendString(0)
83 {
84   //
85   // Constructor that initializes the calibration viewer
86   //
87   fTree = tree;
88   fTreeMustBeDeleted = kFALSE;
89   fListOfObjectsToBeDeleted = new TObjArray();
90   fAbbreviation = "~";
91   fAppendString = ".fElements";
92 }
93
94 //_____________________________________________________________________________
95 AliBaseCalibViewer::AliBaseCalibViewer(const Char_t* fileName, const Char_t* treeName)
96                   :TObject(),
97                    fTree(0),
98                    fFile(0),
99                    fListOfObjectsToBeDeleted(0),
100                    fTreeMustBeDeleted(0),
101                    fAbbreviation(0), 
102                    fAppendString(0)
103                    
104 {
105    //
106    // Constructor to initialize the calibration viewer
107    // the file 'fileName' contains the tree 'treeName'
108    //
109    fFile = new TFile(fileName, "read");
110    fTree = (TTree*) fFile->Get(treeName);
111    fTreeMustBeDeleted = kTRUE;
112    fListOfObjectsToBeDeleted = new TObjArray();
113    fAbbreviation = "~";
114    fAppendString = ".fElements";
115 }
116
117 //____________________________________________________________________________
118 AliBaseCalibViewer & AliBaseCalibViewer::operator =(const AliBaseCalibViewer & param)
119 {
120    //
121    // assignment operator - dummy
122    // not yet working!!!
123    //
124   if(&param == this) return *this;
125   TObject::operator=(param);
126
127
128    fTree = param.fTree;
129    fTreeMustBeDeleted = param.fTreeMustBeDeleted;
130    fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
131    fAbbreviation = param.fAbbreviation;
132    fAppendString = param.fAppendString;
133    return (*this);
134 }
135
136 //_____________________________________________________________________________
137 AliBaseCalibViewer::~AliBaseCalibViewer()
138 {
139    //
140    // AliBaseCalibViewer destructor
141    // all objects will be deleted, the file will be closed, the pictures will disappear
142    //
143    if (fTree && fTreeMustBeDeleted) {
144       fTree->SetCacheSize(0);
145       fTree->Delete();
146    }
147    if (fFile) {
148       fFile->Close();
149       fFile = 0;
150    }
151
152    for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
153       delete fListOfObjectsToBeDeleted->At(i);
154    }
155    delete fListOfObjectsToBeDeleted;
156 }
157
158 //_____________________________________________________________________________
159 void AliBaseCalibViewer::Delete(Option_t* /*option*/) {
160    //
161    // Should be called from AliBaseCalibViewerGUI class only.
162    // If you use Delete() do not call the destructor.
163    // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
164    //
165    
166    if (fTree && fTreeMustBeDeleted) {
167       fTree->SetCacheSize(0);
168       fTree->Delete();
169    }
170    delete fFile;
171    delete fListOfObjectsToBeDeleted;
172 }
173
174 //_____________________________________________________________________________
175 void AliBaseCalibViewer::FormatHistoLabels(TH1 *histo) const {
176    // 
177    // formats title and axis labels of histo 
178    // removes '.fElements'
179    // 
180    if (!histo) return;
181    TString replaceString(fAppendString.Data());
182    TString *str = new TString(histo->GetTitle());
183    str->ReplaceAll(replaceString, "");
184    histo->SetTitle(str->Data());
185    delete str;
186    if (histo->GetXaxis()) {
187       str = new TString(histo->GetXaxis()->GetTitle());
188       str->ReplaceAll(replaceString, "");
189       histo->GetXaxis()->SetTitle(str->Data());
190       delete str;
191    }
192    if (histo->GetYaxis()) {
193       str = new TString(histo->GetYaxis()->GetTitle());
194       str->ReplaceAll(replaceString, "");
195       histo->GetYaxis()->SetTitle(str->Data());
196       delete str;
197    }
198    if (histo->GetZaxis()) {
199       str = new TString(histo->GetZaxis()->GetTitle());
200       str->ReplaceAll(replaceString, "");
201       histo->GetZaxis()->SetTitle(str->Data());
202       delete str;
203    }
204 }
205
206 //_____________________________________________________________________________
207 TFriendElement* AliBaseCalibViewer::AddReferenceTree(const Char_t* filename, const Char_t* treename, const Char_t* refname){
208   //
209   // add a reference tree to the current tree
210   // by default the treename is 'tree' and the reference treename is 'R'
211   //
212    TFile *file = new TFile(filename);
213    fListOfObjectsToBeDeleted->Add(file);
214    TTree * tree = (TTree*)file->Get(treename);
215    return AddFriend(tree, refname);
216 }
217
218 //_____________________________________________________________________________
219 TString* AliBaseCalibViewer::Fit(const Char_t* drawCommand, const Char_t* formula, const Char_t* cuts, 
220                                  Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
221    //
222    // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
223    // returns chi2, fitParam and covMatrix
224    // returns TString with fitted formula
225    //
226    
227    TString formulaStr(formula); 
228    TString drawStr(drawCommand);
229    TString cutStr(cuts);
230    
231    // abbreviations:
232    drawStr.ReplaceAll(fAbbreviation, fAppendString);
233    cutStr.ReplaceAll(fAbbreviation, fAppendString);
234    formulaStr.ReplaceAll(fAbbreviation, fAppendString);
235    
236    formulaStr.ReplaceAll("++", fAbbreviation);
237    TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data()); 
238    Int_t dim = formulaTokens->GetEntriesFast();
239    
240    fitParam.ResizeTo(dim);
241    covMatrix.ResizeTo(dim,dim);
242    
243    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
244    fitter->StoreData(kTRUE);   
245    fitter->ClearPoints();
246    
247    Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
248    if (entries == -1) {
249      delete fitter;
250      return new TString("An ERROR has occured during fitting!");
251    }
252    Double_t **values = new Double_t*[dim+1] ; 
253    
254    for (Int_t i = 0; i < dim + 1; i++){
255       Int_t centries = 0;
256       if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
257       else  centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
258       
259       if (entries != centries) {
260         delete fitter;
261         delete [] values;
262         return new TString("An ERROR has occured during fitting!");
263       }
264       values[i] = new Double_t[entries];
265       memcpy(values[i],  fTree->GetV1(), entries*sizeof(Double_t)); 
266    }
267    
268    // add points to the fitter
269    for (Int_t i = 0; i < entries; i++){
270       Double_t x[1000];
271       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
272       fitter->AddPoint(x, values[dim][i], 1);
273    }
274
275    fitter->Eval();
276    fitter->GetParameters(fitParam);
277    fitter->GetCovarianceMatrix(covMatrix);
278    chi2 = fitter->GetChisquare();
279    
280    TString *preturnFormula = new TString(Form("( %e+",fitParam[0])), &returnFormula = *preturnFormula;
281    
282    for (Int_t iparam = 0; iparam < dim; iparam++) {
283      returnFormula.Append(Form("%s*(%e)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
284      if (iparam < dim-1) returnFormula.Append("+");
285    }
286    returnFormula.Append(" )");
287    delete formulaTokens;
288    delete fitter;
289    for (Int_t i = 0; i < dim + 1; i++) delete [] values[i];
290    delete[] values;
291    return preturnFormula;
292 }
293
294 //_____________________________________________________________________________
295 Double_t AliBaseCalibViewer::GetLTM(Int_t n, Double_t *array, Double_t *sigma, Double_t fraction){
296    //
297    //  returns the LTM and sigma
298    //
299    Double_t *ddata = new Double_t[n];
300    Double_t mean = 0, lsigma = 0;
301    UInt_t nPoints = 0;
302    for (UInt_t i = 0; i < (UInt_t)n; i++) {
303          ddata[nPoints]= array[nPoints];
304          nPoints++;
305    }
306    Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
307    AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
308    if (sigma) *sigma = lsigma;
309    delete [] ddata;
310    return mean;
311 }
312
313 //_____________________________________________________________________________
314 Int_t AliBaseCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
315    // Returns the 'bin' for 'value'
316    // The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
317    // avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
318    /* Begin_Latex
319          GetBin(value) = #frac{nbins - 1}{binUp - binLow} #upoint (value - binLow) +1
320       End_Latex
321    */
322    
323    Int_t bin =  TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
324    // avoid index out of bounds:   
325    if (value < binLow) bin = 0;
326    if (value > binUp)  bin = nbins + 1;
327    return bin;
328    
329 }
330
331 //_____________________________________________________________________________
332 TH1F* AliBaseCalibViewer::SigmaCut(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, 
333                                    Float_t sigmaStep, Bool_t pm) {
334    //
335    // Creates a cumulative histogram Begin_Latex S(t, #mu, #sigma) End_Latex, where you can see, how much of the data are inside sigma-intervals around the mean value
336    // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'histogram'
337    // 'mean' and 'sigma' are Begin_Latex #mu End_Latex and  Begin_Latex #sigma End_Latex of the distribution in 'histogram', to be specified by the user
338    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
339    // sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
340    // pm: Decide weather Begin_Latex t > 0 End_Latex (first case) or Begin_Latex t End_Latex arbitrary (secound case)
341    // The actual work is done on the array.
342    /* Begin_Latex 
343          f(x, #mu, #sigma)     #Rightarrow       S(t, #mu, #sigma) = (#int_{#mu}^{#mu + t #sigma} f(x, #mu, #sigma) dx + #int_{#mu}^{#mu - t #sigma} f(x, #mu, #sigma) dx) / (#int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx),    for  t > 0    
344          or      
345          f(x, #mu, #sigma)     #Rightarrow       S(t, #mu, #sigma) = #int_{#mu}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
346       End_Latex  
347       Begin_Macro(source)
348       {
349          Float_t mean = 0;
350          Float_t sigma = 1.5;
351          Float_t sigmaMax = 4;
352          gROOT->SetStyle("Plain");
353          TH1F *distribution = new TH1F("Distribution1", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
354          TRandom rand(23);
355          for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
356          Float_t *ar = distribution->GetArray();
357          
358          TCanvas* macro_example_canvas = new TCanvas("cAliBaseCalibViewer", "", 350, 350);
359          macro_example_canvas->Divide(0,3);
360          TVirtualPad *pad1 = macro_example_canvas->cd(1);
361          pad1->SetGridy();
362          pad1->SetGridx();
363          distribution->Draw();
364          TVirtualPad *pad2 = macro_example_canvas->cd(2);
365          pad2->SetGridy();
366          pad2->SetGridx();
367          
368          TH1F *shist = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax);
369          shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
370          shist->Draw();  
371          TVirtualPad *pad3 = macro_example_canvas->cd(3);
372          pad3->SetGridy();
373          pad3->SetGridx();
374          TH1F *shistPM = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax, -1, kTRUE);
375          shistPM->Draw();   
376          return macro_example_canvas;
377       }  
378       End_Macro
379    */ 
380    
381    Float_t *array = histogram->GetArray();
382    Int_t    nbins = histogram->GetXaxis()->GetNbins();
383    Float_t binLow = histogram->GetXaxis()->GetXmin();
384    Float_t binUp  = histogram->GetXaxis()->GetXmax();
385    return AliBaseCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
386 }   
387    
388 //_____________________________________________________________________________
389 TH1F* AliBaseCalibViewer::SigmaCut(Int_t n, Float_t *array, Float_t mean, Float_t sigma, Int_t nbins, Float_t binLow, Float_t binUp, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm){
390    //
391    // Creates a histogram Begin_Latex S(t, #mu, #sigma) End_Latex, where you can see, how much of the data are inside sigma-intervals around the mean value
392    // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
393    // 'mean' and 'sigma' are Begin_Latex #mu End_Latex and  Begin_Latex #sigma End_Latex of the distribution in 'array', to be specified by the user
394    // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
395    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
396    // sigmaStep: the binsize of the generated histogram
397    // Here the actual work is done.
398    
399    if (TMath::Abs(sigma) < 1.e-10) return 0;
400    Float_t binWidth = (binUp-binLow)/(nbins - 1);
401    if (sigmaStep <= 0) sigmaStep = binWidth;
402    Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1  due to overflow bin in histograms
403    if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
404    Float_t kbinLow = !pm ? 0 : -sigmaMax;
405    Float_t kbinUp  = sigmaMax;
406    TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
407    hist->SetDirectory(0);
408    hist->Reset();
409    
410    // calculate normalization
411    Double_t normalization = 0;
412    for (Int_t i = 0; i <= n; i++) {
413         normalization += array[i];
414    }
415    
416    // given units: units from given histogram
417    // sigma units: in units of sigma
418    // iDelta: integrate in interval (mean +- iDelta), given units
419    // x:      ofset from mean for integration, given units
420    // hist:   needs 
421    
422    // fill histogram
423    for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
424       // integrate array
425       Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
426       Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
427       // add bin of mean value only once to the histogram
428       for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
429          valueP += (mean + x <= binUp)  ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
430          valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0; 
431       }
432
433       if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueP, normalization);
434       if (valueP / normalization > 100) return hist;
435       if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueM, normalization);
436       if (valueM / normalization > 100) return hist;
437       valueP = (valueP / normalization);
438       valueM = (valueM / normalization);
439       if (pm) {
440          Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
441          hist->SetBinContent(bin, valueP);
442          bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
443          hist->SetBinContent(bin, valueM);
444       }
445       else { // if (!pm)
446          Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
447          hist->SetBinContent(bin, valueP + valueM);
448       }
449    }
450    if (!pm) hist->SetMaximum(1.2);
451    return hist;
452 }
453
454 //_____________________________________________________________________________
455 TH1F* AliBaseCalibViewer::SigmaCut(Int_t /*n*/, Double_t */*array*/, Double_t /*mean*/, Double_t /*sigma*/, 
456                                    Int_t /*nbins*/, Double_t */*xbins*/, Double_t /*sigmaMax*/){
457    // 
458    // SigmaCut for variable binsize
459    // NOT YET IMPLEMENTED !!!
460    // 
461    printf("SigmaCut with variable binsize, Not yet implemented\n");
462
463    return 0;
464 }
465
466 //_____________________________________________________________________________
467 Int_t  AliBaseCalibViewer::DrawHisto1D(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
468                                        const Char_t *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const
469 {
470    //
471    // Easy drawing of data, in principle the same as EasyDraw1D
472    // Difference: A line for the mean / median / LTM is drawn
473    // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
474    // example: sigmas = "2; 4; 6;"  at Begin_Latex 2 #sigma End_Latex, Begin_Latex 4 #sigma End_Latex and Begin_Latex 6 #sigma End_Latex  a line is drawn.
475    // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
476    //
477   Int_t oldOptStat = gStyle->GetOptStat();
478   gStyle->SetOptStat(0000000);
479   Double_t ltmFraction = 0.8;
480   
481   TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
482   TVectorF nsigma(sigmasTokens->GetEntriesFast());
483   for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
484     TString str(((TObjString*)sigmasTokens->At(i))->GetString());
485     Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
486     nsigma[i] = sig;
487   }
488   
489   TString drawStr(drawCommand);
490   Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
491   if (dangerousToDraw) {
492     Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
493     return -1;
494   }
495   drawStr += " >> tempHist";
496   Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
497   TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
498    // FIXME is this histogram deleted automatically?
499   Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
500   
501   Double_t mean = TMath::Mean(entries, values);
502   Double_t median = TMath::Median(entries, values);
503   Double_t sigma = TMath::RMS(entries, values);
504   Double_t maxY = htemp->GetMaximum();
505   
506   TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
507    //fListOfObjectsToBeDeleted->Add(legend);
508   
509   if (plotMean) {
510       // draw Mean
511     TLine* line = new TLine(mean, 0, mean, maxY);
512       //fListOfObjectsToBeDeleted->Add(line);
513     line->SetLineColor(kRed);
514     line->SetLineWidth(2);
515     line->SetLineStyle(1);
516     line->Draw();
517     legend->AddEntry(line, Form("Mean: %f", mean), "l");
518       // draw sigma lines
519     for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
520       TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
521          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
522       linePlusSigma->SetLineColor(kRed);
523       linePlusSigma->SetLineStyle(2 + i);
524       linePlusSigma->Draw();
525       TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
526          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
527       lineMinusSigma->SetLineColor(kRed);
528       lineMinusSigma->SetLineStyle(2 + i);
529       lineMinusSigma->Draw();
530       legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
531     }
532   }
533   if (plotMedian) {
534       // draw median
535     TLine* line = new TLine(median, 0, median, maxY);
536       //fListOfObjectsToBeDeleted->Add(line);
537     line->SetLineColor(kBlue);
538     line->SetLineWidth(2);
539     line->SetLineStyle(1);
540     line->Draw();
541     legend->AddEntry(line, Form("Median: %f", median), "l");
542       // draw sigma lines
543     for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
544       TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
545          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
546       linePlusSigma->SetLineColor(kBlue);
547       linePlusSigma->SetLineStyle(2 + i);
548       linePlusSigma->Draw();
549       TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
550          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
551       lineMinusSigma->SetLineColor(kBlue);
552       lineMinusSigma->SetLineStyle(2 + i);
553       lineMinusSigma->Draw();
554       legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
555     }
556   }
557   if (plotLTM) {
558       // draw LTM
559     Double_t ltmRms = 0;
560     Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
561     TLine* line = new TLine(ltm, 0, ltm, maxY);
562       //fListOfObjectsToBeDeleted->Add(line);
563     line->SetLineColor(kGreen+2);
564     line->SetLineWidth(2);
565     line->SetLineStyle(1);
566     line->Draw();
567     legend->AddEntry(line, Form("LTM: %f", ltm), "l");
568       // draw sigma lines
569     for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
570       TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
571          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
572       linePlusSigma->SetLineColor(kGreen+2);
573       linePlusSigma->SetLineStyle(2+i);
574       linePlusSigma->Draw();
575       
576       TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
577          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
578       lineMinusSigma->SetLineColor(kGreen+2);
579       lineMinusSigma->SetLineStyle(2+i);
580       lineMinusSigma->Draw();
581       legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms)), "l");
582     }
583   }
584   if (!plotMean && !plotMedian && !plotLTM) return -1;
585   legend->Draw();
586   gStyle->SetOptStat(oldOptStat);
587   return 1;
588 }
589
590 //_____________________________________________________________________________
591 Int_t AliBaseCalibViewer::SigmaCut(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
592                                    Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, 
593                                    const Char_t *sigmas, Float_t sigmaStep) const {
594    //
595    // Creates a histogram, where you can see, how much of the data are inside sigma-intervals 
596    // around the mean/median/LTM
597    // with drawCommand, sector and cuts you specify your input data, see EasyDraw
598    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
599    // sigmaStep: the binsize of the generated histogram
600    // plotMean/plotMedian/plotLTM: specifies where to put the center
601    //
602   
603    Double_t ltmFraction = 0.8;
604    
605    TString drawStr(drawCommand);
606    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
607    if (dangerousToDraw) {
608       Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
609       return -1;
610    }
611    drawStr += " >> tempHist";
612    
613    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
614    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
615    // FIXME is this histogram deleted automatically?
616    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
617    
618    Double_t mean = TMath::Mean(entries, values);
619    Double_t median = TMath::Median(entries, values);
620    Double_t sigma = TMath::RMS(entries, values);
621    
622    TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
623    //fListOfObjectsToBeDeleted->Add(legend);
624    TH1F *cutHistoMean = 0;
625    TH1F *cutHistoMedian = 0;
626    TH1F *cutHistoLTM = 0;
627    
628    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
629    TVectorF nsigma(sigmasTokens->GetEntriesFast());
630    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
631       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
632       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
633       nsigma[i] = sig;
634    }
635   
636    if (plotMean) {
637       cutHistoMean = SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
638       if (cutHistoMean) {
639          //fListOfObjectsToBeDeleted->Add(cutHistoMean);
640          cutHistoMean->SetLineColor(kRed);
641          legend->AddEntry(cutHistoMean, "Mean", "l");
642          cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
643          cutHistoMean->Draw();
644          DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
645       } // if (cutHistoMean)
646        
647    }
648    if (plotMedian) {
649       cutHistoMedian = SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
650       if (cutHistoMedian) {
651          //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
652          cutHistoMedian->SetLineColor(kBlue);
653          legend->AddEntry(cutHistoMedian, "Median", "l");
654          cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
655          if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
656             else cutHistoMedian->Draw();
657          DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
658       }  // if (cutHistoMedian)
659    }
660    if (plotLTM) {
661       Double_t ltmRms = 0;
662       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
663       cutHistoLTM = SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
664       if (cutHistoLTM) {
665          //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
666          cutHistoLTM->SetLineColor(kGreen+2);
667          legend->AddEntry(cutHistoLTM, "LTM", "l");
668          cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
669          if ((plotMean && cutHistoMean) || (plotMedian && cutHistoMedian)) cutHistoLTM->Draw("same");
670             else cutHistoLTM->Draw();
671          DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
672       }
673    }
674    if (!plotMean && !plotMedian && !plotLTM) return -1;
675    legend->Draw();
676    return 1;
677 }
678
679 //_____________________________________________________________________________
680 Int_t AliBaseCalibViewer::Integrate(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
681                                     Float_t /*sigmaMax*/, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, 
682                                     const Char_t *sigmas, Float_t /*sigmaStep*/) const {
683    //
684    // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"   
685    // "mean" and "sigma" are Begin_Latex #mu End_Latex and  Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
686    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
687    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
688    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
689    // The actual work is done on the array.
690    /* Begin_Latex 
691          f(x, #mu, #sigma)     #Rightarrow       S(t, #mu, #sigma) = #int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
692       End_Latex  
693    */
694    
695    Double_t ltmFraction = 0.8;
696    
697    TString drawStr(drawCommand);
698    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
699    if (dangerousToDraw) {
700       Warning("Integrate", "The draw string must not contain ':' or '>>'.");
701       return -1;
702    }
703    drawStr += " >> tempHist";
704    
705    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
706    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
707    TGraph *integralGraphMean   = 0;
708    TGraph *integralGraphMedian = 0;
709    TGraph *integralGraphLTM    = 0;
710    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
711    Int_t    *index  = new Int_t[entries];
712    Float_t  *xarray = new Float_t[entries];
713    Float_t  *yarray = new Float_t[entries];
714    TMath::Sort(entries, values, index, kFALSE);
715    
716    Double_t mean = TMath::Mean(entries, values);
717    Double_t median = TMath::Median(entries, values);
718    Double_t sigma = TMath::RMS(entries, values);
719    
720    // parse sigmas string
721    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
722    TVectorF nsigma(sigmasTokens->GetEntriesFast());
723    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
724       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
725       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
726       nsigma[i] = sig;
727    }
728   
729    TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
730    //fListOfObjectsToBeDeleted->Add(legend);
731   
732    if (plotMean) {
733       for (Int_t i = 0; i < entries; i++) {
734          xarray[i] = (values[index[i]] - mean) / sigma; 
735          yarray[i] = float(i) / float(entries);
736       }
737       integralGraphMean = new TGraph(entries, xarray, yarray);
738       if (integralGraphMean) {
739          //fListOfObjectsToBeDeleted->Add(integralGraphMean);
740          integralGraphMean->SetLineColor(kRed);
741          legend->AddEntry(integralGraphMean, "Mean", "l");
742          integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
743          integralGraphMean->Draw("alu");
744          DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
745       }
746    }
747    if (plotMedian) {
748       for (Int_t i = 0; i < entries; i++) {
749          xarray[i] = (values[index[i]] - median) / sigma; 
750          yarray[i] = float(i) / float(entries);
751       }
752       integralGraphMedian = new TGraph(entries, xarray, yarray);
753       if (integralGraphMedian) {
754          //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
755          integralGraphMedian->SetLineColor(kBlue);
756          legend->AddEntry(integralGraphMedian, "Median", "l");
757          integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
758          if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
759             else integralGraphMedian->Draw("alu");
760          DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
761       }
762    }
763    if (plotLTM) {
764       Double_t ltmRms = 0;
765       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
766       for (Int_t i = 0; i < entries; i++) {
767          xarray[i] = (values[index[i]] - ltm) / ltmRms; 
768          yarray[i] = float(i) / float(entries);
769       }
770       integralGraphLTM = new TGraph(entries, xarray, yarray);
771       if (integralGraphLTM) {
772          //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
773          integralGraphLTM->SetLineColor(kGreen+2);
774          legend->AddEntry(integralGraphLTM, "LTM", "l");
775          integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
776          if ((plotMean && integralGraphMean) || (plotMedian && integralGraphMedian)) integralGraphLTM->Draw("samelu");
777             else integralGraphLTM->Draw("alu");
778          DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
779       }
780    }
781    delete [] index;
782    delete [] xarray;
783    delete [] yarray;
784    if (!plotMean && !plotMedian && !plotLTM) return -1;
785    legend->Draw();
786    return entries;
787 }
788
789 //_____________________________________________________________________________
790 TH1F* AliBaseCalibViewer::Integrate(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
791    //
792    // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"   
793    // "mean" and "sigma" are Begin_Latex #mu End_Latex and  Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
794    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
795    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
796    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
797    // The actual work is done on the array.
798    /* Begin_Latex 
799          f(x, #mu, #sigma)     #Rightarrow       S(t, #mu, #sigma) = #int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
800       End_Latex  
801       Begin_Macro(source)
802       {
803          Float_t mean = 0;
804          Float_t sigma = 1.5;
805          Float_t sigmaMax = 4;
806          gROOT->SetStyle("Plain");
807          TH1F *distribution = new TH1F("Distribution2", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
808          TRandom rand(23);
809          for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
810          Float_t *ar = distribution->GetArray();
811          
812          TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_Integrate", "", 350, 350);
813          macro_example_canvas->Divide(0,2);
814          TVirtualPad *pad1 = macro_example_canvas->cd(1);
815          pad1->SetGridy();
816          pad1->SetGridx();
817          distribution->Draw();
818          TVirtualPad *pad2 = macro_example_canvas->cd(2);
819          pad2->SetGridy();
820          pad2->SetGridx();
821          TH1F *shist = AliTPCCalibViewer::Integrate(distribution, mean, sigma, sigmaMax);
822          shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
823          shist->Draw();  
824          
825          return macro_example_canvas_Integrate;
826       }  
827       End_Macro
828    */ 
829
830    
831    Float_t *array = histogram->GetArray();
832    Int_t    nbins = histogram->GetXaxis()->GetNbins();
833    Float_t binLow = histogram->GetXaxis()->GetXmin();
834    Float_t binUp  = histogram->GetXaxis()->GetXmax();
835    return Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
836 }
837
838 //_____________________________________________________________________________
839 TH1F* AliBaseCalibViewer::Integrate(Int_t n, Float_t *array, Int_t nbins, Float_t binLow, Float_t binUp, 
840                                     Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
841    // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"   
842    // "mean" and "sigma" are Begin_Latex #mu End_Latex and  Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
843    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
844    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
845    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
846    // Here the actual work is done.
847       
848    Bool_t givenUnits = kTRUE;
849    if (TMath::Abs(sigma) < 1.e-10 && TMath::Abs(sigmaMax) < 1.e-10) givenUnits = kFALSE;
850    if (givenUnits) {
851       sigma = 1;
852       sigmaMax = (binUp - binLow) / 2.;
853    }
854    
855    Float_t binWidth = (binUp-binLow)/(nbins - 1);
856    if (sigmaStep <= 0) sigmaStep = binWidth;
857    Int_t kbins =  (Int_t)(sigmaMax * sigma / sigmaStep) + 1;  // + 1  due to overflow bin in histograms
858    Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
859    Float_t kbinUp  = givenUnits ? binUp  : sigmaMax;
860    TH1F *hist = 0; 
861    if (givenUnits)  hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp); 
862    if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
863    hist->SetDirectory(0);
864    hist->Reset();
865    
866    // calculate normalization
867  //  printf("calculating normalization, integrating from bin 1 to %i \n", n);
868    Double_t normalization = 0;
869    for (Int_t i = 1; i <= n; i++) {
870         normalization += array[i];
871    }
872  //  printf("normalization: %f \n", normalization);
873    
874    // given units: units from given histogram
875    // sigma units: in units of sigma
876    // iDelta: integrate in interval (mean +- iDelta), given units
877    // x:      ofset from mean for integration, given units
878    // hist:   needs 
879    
880    // fill histogram
881    for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
882       // integrate array
883       Double_t value = 0;
884       for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
885          value += (x <= binUp && x >= binLow)  ? array[GetBin(x, nbins, binLow, binUp)] : 0;
886       }
887       if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", value, normalization);
888       if (value / normalization > 100) return hist;
889       Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
890     //  printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
891     //  printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
892       value = (value / normalization);
893       hist->SetBinContent(bin, value);
894    }
895    return hist;
896 }
897
898 //_____________________________________________________________________________
899 void AliBaseCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
900    //
901    // Private function for SigmaCut(...) and Integrate(...)
902    // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
903    //
904   
905    // start to draw the lines, loop over requested sigmas
906   for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
907     if (!pm) {
908       Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
909       TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
910          //fListOfObjectsToBeDeleted->Add(lineUp);
911       lineUp->SetLineColor(color);
912       lineUp->SetLineStyle(2 + i);
913       lineUp->Draw();
914       TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
915          //fListOfObjectsToBeDeleted->Add(lineLeft);
916       lineLeft->SetLineColor(color);
917       lineLeft->SetLineStyle(2 + i);
918       lineLeft->Draw();
919       legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
920     }
921     else { // if (pm)
922       Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
923       TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
924          //fListOfObjectsToBeDeleted->Add(lineUp1);
925       lineUp1->SetLineColor(color);
926       lineUp1->SetLineStyle(2 + i);
927       lineUp1->Draw();
928       TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
929          //fListOfObjectsToBeDeleted->Add(lineLeft1);
930       lineLeft1->SetLineColor(color);
931       lineLeft1->SetLineStyle(2 + i);
932       lineLeft1->Draw();
933       legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
934       bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
935       TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
936          //fListOfObjectsToBeDeleted->Add(lineUp2);
937       lineUp2->SetLineColor(color);
938       lineUp2->SetLineStyle(2 + i);
939       lineUp2->Draw();
940       TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
941          //fListOfObjectsToBeDeleted->Add(lineLeft2);
942       lineLeft2->SetLineColor(color);
943       lineLeft2->SetLineStyle(2 + i);
944       lineLeft2->Draw();
945       legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
946     }
947   }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
948 }
949
950 //_____________________________________________________________________________
951 void AliBaseCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
952    //
953    // Private function for SigmaCut(...) and Integrate(...)
954    // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
955    //
956   
957    // start to draw the lines, loop over requested sigmas
958   for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
959     if (!pm) {
960       TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
961          //fListOfObjectsToBeDeleted->Add(lineUp);
962       lineUp->SetLineColor(color);
963       lineUp->SetLineStyle(2 + i);
964       lineUp->Draw();
965       TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
966          //fListOfObjectsToBeDeleted->Add(lineLeft);
967       lineLeft->SetLineColor(color);
968       lineLeft->SetLineStyle(2 + i);
969       lineLeft->Draw();
970       legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
971     }
972     else { // if (pm)
973       TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
974          //fListOfObjectsToBeDeleted->Add(lineUp1);
975       lineUp1->SetLineColor(color);
976       lineUp1->SetLineStyle(2 + i);
977       lineUp1->Draw();
978       TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
979          //fListOfObjectsToBeDeleted->Add(lineLeft1);
980       lineLeft1->SetLineColor(color);
981       lineLeft1->SetLineStyle(2 + i);
982       lineLeft1->Draw();
983       legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
984       TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
985          //fListOfObjectsToBeDeleted->Add(lineUp2);
986       lineUp2->SetLineColor(color);
987       lineUp2->SetLineStyle(2 + i);
988       lineUp2->Draw();
989       TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
990          //fListOfObjectsToBeDeleted->Add(lineLeft2);
991       lineLeft2->SetLineColor(color);
992       lineLeft2->SetLineStyle(2 + i);
993       lineLeft2->Draw();
994       legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i])), "l");
995     }
996   }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
997 }