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