ROOT logo
///////////////////////////////////////////////////////////////////////////////
//                                                                           //
//  Base class for the AliTPCCalibViewer and AliTRDCalibViewer               //
//  used for the calibration monitor                                         //
//                                                                           //
//  Authors:     Marian Ivanov (Marian.Ivanov@cern.ch)                       //
//               Jens Wiechula (Jens.Wiechula@cern.ch)                       //
//               Ionut Arsene  (iarsene@cern.ch)                             //
//                                                                           //
///////////////////////////////////////////////////////////////////////////////


#include <iostream>
#include <fstream>
#include <TString.h>
#include <TRandom.h>
#include <TLegend.h>
#include <TLine.h>
//#include <TCanvas.h>
#include <TROOT.h>
#include <TStyle.h>
#include <TH1.h> 
#include <TH1F.h>
#include <TMath.h>
#include <THashTable.h>
#include <TObjString.h>
#include <TLinearFitter.h>
#include <TFile.h>
#include <TKey.h>
#include <TGraph.h>
#include <TDirectory.h>
#include <TFriendElement.h>

#include "TTreeStream.h"
#include "AliBaseCalibViewer.h"

ClassImp(AliBaseCalibViewer)

AliBaseCalibViewer::AliBaseCalibViewer()
:TObject(),
  fTree(0),
  fFile(0),
  fListOfObjectsToBeDeleted(0),
  fTreeMustBeDeleted(0), 
  fAbbreviation(0), 
  fAppendString(0)
{
  //
  // Default constructor
  //
}

//_____________________________________________________________________________
AliBaseCalibViewer::AliBaseCalibViewer(const AliBaseCalibViewer &c)
:TObject(c),
  fTree(0),
  fFile(0),
  fListOfObjectsToBeDeleted(0),
  fTreeMustBeDeleted(0),
  fAbbreviation(0), 
  fAppendString(0)
{
  //
  // dummy AliBaseCalibViewer copy constructor
  // not yet working!!!
  //
  fTree = c.fTree;
  fTreeMustBeDeleted = c.fTreeMustBeDeleted;
  fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
  fAbbreviation = c.fAbbreviation;
  fAppendString = c.fAppendString;
}

//_____________________________________________________________________________
AliBaseCalibViewer::AliBaseCalibViewer(TTree* tree)
:TObject(),
  fTree(0),
  fFile(0),
  fListOfObjectsToBeDeleted(0),
  fTreeMustBeDeleted(0),
  fAbbreviation(0), 
  fAppendString(0)
{
  //
  // Constructor that initializes the calibration viewer
  //
  fTree = tree;
  fTreeMustBeDeleted = kFALSE;
  fListOfObjectsToBeDeleted = new TObjArray();
  fAbbreviation = "~";
  fAppendString = ".fElements";
}

//_____________________________________________________________________________
AliBaseCalibViewer::AliBaseCalibViewer(const Char_t* fileName, const Char_t* treeName)
:TObject(),
  fTree(0),
  fFile(0),
  fListOfObjectsToBeDeleted(0),
  fTreeMustBeDeleted(0),
  fAbbreviation(0), 
  fAppendString(0)

{
  //
  // Constructor to initialize the calibration viewer
  // the file 'fileName' contains the tree 'treeName'
  //
  fFile = new TFile(fileName, "read");
  fTree = (TTree*) fFile->Get(treeName);
  fTreeMustBeDeleted = kTRUE;
  fListOfObjectsToBeDeleted = new TObjArray();
  fAbbreviation = "~";
  fAppendString = ".fElements";
}

//____________________________________________________________________________
AliBaseCalibViewer & AliBaseCalibViewer::operator =(const AliBaseCalibViewer & param)
{
  //
  // assignment operator - dummy
  // not yet working!!!
  //
  if(&param == this) return *this;
  TObject::operator=(param);


  fTree = param.fTree;
  fTreeMustBeDeleted = param.fTreeMustBeDeleted;
  fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
  fAbbreviation = param.fAbbreviation;
  fAppendString = param.fAppendString;
  return (*this);
}

//_____________________________________________________________________________
AliBaseCalibViewer::~AliBaseCalibViewer()
{
  //
  // AliBaseCalibViewer destructor
  // all objects will be deleted, the file will be closed, the pictures will disappear
  //
  if (fTree && fTreeMustBeDeleted) {
    fTree->SetCacheSize(0);
    fTree->Delete();
  }
  if (fFile) {
    fFile->Close();
    fFile = 0;
  }

  for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
    delete fListOfObjectsToBeDeleted->At(i);
  }
  delete fListOfObjectsToBeDeleted;
}

//_____________________________________________________________________________
void AliBaseCalibViewer::Delete(Option_t* /*option*/) {
  //
  // Should be called from AliBaseCalibViewerGUI class only.
  // If you use Delete() do not call the destructor.
  // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
  //

  if (fTree && fTreeMustBeDeleted) {
    fTree->SetCacheSize(0);
    fTree->Delete();
  }
  delete fFile;
  delete fListOfObjectsToBeDeleted;
}

//_____________________________________________________________________________
void AliBaseCalibViewer::FormatHistoLabels(TH1 *histo) const {
  // 
  // formats title and axis labels of histo 
  // removes '.fElements'
  // 
  if (!histo) return;
  TString replaceString(fAppendString.Data());
  TString *str = new TString(histo->GetTitle());
  str->ReplaceAll(replaceString, "");
  histo->SetTitle(str->Data());
  delete str;
  if (histo->GetXaxis()) {
    str = new TString(histo->GetXaxis()->GetTitle());
    str->ReplaceAll(replaceString, "");
    histo->GetXaxis()->SetTitle(str->Data());
    delete str;
  }
  if (histo->GetYaxis()) {
    str = new TString(histo->GetYaxis()->GetTitle());
    str->ReplaceAll(replaceString, "");
    histo->GetYaxis()->SetTitle(str->Data());
    delete str;
  }
  if (histo->GetZaxis()) {
    str = new TString(histo->GetZaxis()->GetTitle());
    str->ReplaceAll(replaceString, "");
    histo->GetZaxis()->SetTitle(str->Data());
    delete str;
  }
}

//_____________________________________________________________________________
TFriendElement* AliBaseCalibViewer::AddReferenceTree(const Char_t* filename, const Char_t* treename, const Char_t* refname){
  //
  // add a reference tree to the current tree
  // by default the treename is 'tree' and the reference treename is 'R'
  //
  TFile *file = new TFile(filename);
  fListOfObjectsToBeDeleted->Add(file);
  TTree * tree = (TTree*)file->Get(treename);
  return AddFriend(tree, refname);
}

//_____________________________________________________________________________
TString* AliBaseCalibViewer::Fit(const Char_t* drawCommand, const Char_t* formula, const Char_t* cuts, 
    Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
  //
  // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
  // returns chi2, fitParam and covMatrix
  // returns TString with fitted formula
  //

  TString formulaStr(formula); 
  TString drawStr(drawCommand);
  TString cutStr(cuts);

  // abbreviations:
  drawStr.ReplaceAll(fAbbreviation, fAppendString);
  cutStr.ReplaceAll(fAbbreviation, fAppendString);
  formulaStr.ReplaceAll(fAbbreviation, fAppendString);

  formulaStr.ReplaceAll("++", fAbbreviation);
  TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data()); 
  Int_t dim = formulaTokens->GetEntriesFast();

  fitParam.ResizeTo(dim);
  covMatrix.ResizeTo(dim,dim);

  TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
  fitter->StoreData(kTRUE);   
  fitter->ClearPoints();

  Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
  if (entries == -1) {
    delete fitter;
    delete formulaTokens;
    return new TString("An ERROR has occured during fitting!");
  }
  Double_t **values = new Double_t*[dim+1] ; 

  for (Int_t i = 0; i < dim + 1; i++){
    Int_t centries = 0;
    if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
    else  centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");

    if (entries != centries) {
      delete fitter;
      delete [] values;
      delete formulaTokens;
      return new TString("An ERROR has occured during fitting!");
    }
    values[i] = new Double_t[entries];
    memcpy(values[i],  fTree->GetV1(), entries*sizeof(Double_t)); 
  }

  // add points to the fitter
  for (Int_t i = 0; i < entries; i++){
    Double_t x[1000];
    for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
    fitter->AddPoint(x, values[dim][i], 1);
  }

  fitter->Eval();
  fitter->GetParameters(fitParam);
  fitter->GetCovarianceMatrix(covMatrix);
  chi2 = fitter->GetChisquare();

  TString *preturnFormula = new TString(Form("( %e+",fitParam[0])), &returnFormula = *preturnFormula;

  for (Int_t iparam = 0; iparam < dim; iparam++) {
    returnFormula.Append(Form("%s*(%e)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
    if (iparam < dim-1) returnFormula.Append("+");
  }
  returnFormula.Append(" )");
  delete formulaTokens;
  delete fitter;
  for (Int_t i = 0; i < dim + 1; i++) delete [] values[i];
  delete[] values;
  return preturnFormula;
}

//_____________________________________________________________________________
Double_t AliBaseCalibViewer::GetLTM(Int_t n, Double_t *array, Double_t *sigma, Double_t fraction){
  //
  //  returns the LTM and sigma
  //
  Double_t *ddata = new Double_t[n];
  Double_t mean = 0, lsigma = 0;
  UInt_t nPoints = 0;
  for (UInt_t i = 0; i < (UInt_t)n; i++) {
    ddata[nPoints]= array[nPoints];
    nPoints++;
  }
  Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
  AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
  if (sigma) *sigma = lsigma;
  delete [] ddata;
  return mean;
}

//_____________________________________________________________________________
Int_t AliBaseCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
  // Returns the 'bin' for 'value'
  // The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
  // avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
  /* Begin_Latex
     GetBin(value) = #frac{nbins - 1}{binUp - binLow} #upoint (value - binLow) +1
     End_Latex
     */

  Int_t bin =  TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
  // avoid index out of bounds:   
  if (value < binLow) bin = 0;
  if (value > binUp)  bin = nbins + 1;
  return bin;

}

//_____________________________________________________________________________
TH1F* AliBaseCalibViewer::SigmaCut(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, 
    Float_t sigmaStep, Bool_t pm) {
  //
  // 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
  // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'histogram'
  // '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
  // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
  // sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
  // pm: Decide weather Begin_Latex t > 0 End_Latex (first case) or Begin_Latex t End_Latex arbitrary (secound case)
  // The actual work is done on the array.
  /* Begin_Latex 
     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    
     or      
     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
     End_Latex  
     Begin_Macro(source)
     {
     Float_t mean = 0;
     Float_t sigma = 1.5;
     Float_t sigmaMax = 4;
     gROOT->SetStyle("Plain");
     TH1F *distribution = new TH1F("Distribution1", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
     TRandom rand(23);
     for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
     Float_t *ar = distribution->GetArray();

     TCanvas* macro_example_canvas = new TCanvas("cAliBaseCalibViewer", "", 350, 350);
     macro_example_canvas->Divide(0,3);
     TVirtualPad *pad1 = macro_example_canvas->cd(1);
     pad1->SetGridy();
     pad1->SetGridx();
     distribution->Draw();
     TVirtualPad *pad2 = macro_example_canvas->cd(2);
     pad2->SetGridy();
     pad2->SetGridx();

     TH1F *shist = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax);
     shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
     shist->Draw();  
     TVirtualPad *pad3 = macro_example_canvas->cd(3);
     pad3->SetGridy();
     pad3->SetGridx();
     TH1F *shistPM = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax, -1, kTRUE);
     shistPM->Draw();   
     return macro_example_canvas;
     }  
     End_Macro
     */ 

  Float_t *array = histogram->GetArray();
  Int_t    nbins = histogram->GetXaxis()->GetNbins();
  Float_t binLow = histogram->GetXaxis()->GetXmin();
  Float_t binUp  = histogram->GetXaxis()->GetXmax();
  return AliBaseCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
}   

//_____________________________________________________________________________
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){
  //
  // 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
  // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
  // '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
  // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
  // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
  // sigmaStep: the binsize of the generated histogram
  // Here the actual work is done.

  if (TMath::Abs(sigma) < 1.e-10) return 0;
  Float_t binWidth = (binUp-binLow)/(nbins - 1);
  if (sigmaStep <= 0) sigmaStep = binWidth;
  Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1  due to overflow bin in histograms
  if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
  Float_t kbinLow = !pm ? 0 : -sigmaMax;
  Float_t kbinUp  = sigmaMax;
  TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
  hist->SetDirectory(0);
  hist->Reset();

  // calculate normalization
  Double_t normalization = 0;
  for (Int_t i = 0; i <= n; i++) {
    normalization += array[i];
  }

  // given units: units from given histogram
  // sigma units: in units of sigma
  // iDelta: integrate in interval (mean +- iDelta), given units
  // x:      ofset from mean for integration, given units
  // hist:   needs 

  // fill histogram
  for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
    // integrate array
    Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
    Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
    // add bin of mean value only once to the histogram
    for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
      valueP += (mean + x <= binUp)  ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
      valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0; 
    }

    if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueP, normalization);
    if (valueP / normalization > 100) return hist;
    if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueM, normalization);
    if (valueM / normalization > 100) return hist;
    valueP = (valueP / normalization);
    valueM = (valueM / normalization);
    if (pm) {
      Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
      hist->SetBinContent(bin, valueP);
      bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
      hist->SetBinContent(bin, valueM);
    }
    else { // if (!pm)
      Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
      hist->SetBinContent(bin, valueP + valueM);
    }
  }
  if (!pm) hist->SetMaximum(1.2);
  return hist;
}

//_____________________________________________________________________________
TH1F* AliBaseCalibViewer::SigmaCut(Int_t /*n*/, Double_t */*array*/, Double_t /*mean*/, Double_t /*sigma*/, 
    Int_t /*nbins*/, Double_t */*xbins*/, Double_t /*sigmaMax*/){
  // 
  // SigmaCut for variable binsize
  // NOT YET IMPLEMENTED !!!
  // 
  printf("SigmaCut with variable binsize, Not yet implemented\n");

  return 0;
}

//_____________________________________________________________________________
Int_t  AliBaseCalibViewer::DrawHisto1D(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
    const Char_t *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const
{
  //
  // Easy drawing of data, in principle the same as EasyDraw1D
  // Difference: A line for the mean / median / LTM is drawn
  // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
  // 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.
  // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
  //
  Int_t oldOptStat = gStyle->GetOptStat();
  gStyle->SetOptStat(0000000);
  Double_t ltmFraction = 0.8;

  TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
  TVectorF nsigma(sigmasTokens->GetEntriesFast());
  for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
    TString str(((TObjString*)sigmasTokens->At(i))->GetString());
    Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
    nsigma[i] = sig;
  }
  delete sigmasTokens;
  //
  TString drawStr(drawCommand);
  Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
  if (dangerousToDraw) {
    Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
    return -1;
  }
  drawStr += " >> tempHist";
  Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
  TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
  // FIXME is this histogram deleted automatically?
  Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers

  Double_t mean = TMath::Mean(entries, values);
  Double_t median = TMath::Median(entries, values);
  Double_t sigma = TMath::RMS(entries, values);
  Double_t maxY = htemp->GetMaximum();

  TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
  //fListOfObjectsToBeDeleted->Add(legend);

  if (plotMean) {
    // draw Mean
    TLine* line = new TLine(mean, 0, mean, maxY);
    //fListOfObjectsToBeDeleted->Add(line);
    line->SetLineColor(kRed);
    line->SetLineWidth(2);
    line->SetLineStyle(1);
    line->Draw();
    legend->AddEntry(line, Form("Mean: %f", mean), "l");
    // draw sigma lines
    for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
      TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
      //fListOfObjectsToBeDeleted->Add(linePlusSigma);
      linePlusSigma->SetLineColor(kRed);
      linePlusSigma->SetLineStyle(2 + i);
      linePlusSigma->Draw();
      TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
      //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
      lineMinusSigma->SetLineColor(kRed);
      lineMinusSigma->SetLineStyle(2 + i);
      lineMinusSigma->Draw();
      legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
    }
  }
  if (plotMedian) {
    // draw median
    TLine* line = new TLine(median, 0, median, maxY);
    //fListOfObjectsToBeDeleted->Add(line);
    line->SetLineColor(kBlue);
    line->SetLineWidth(2);
    line->SetLineStyle(1);
    line->Draw();
    legend->AddEntry(line, Form("Median: %f", median), "l");
    // draw sigma lines
    for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
      TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
      //fListOfObjectsToBeDeleted->Add(linePlusSigma);
      linePlusSigma->SetLineColor(kBlue);
      linePlusSigma->SetLineStyle(2 + i);
      linePlusSigma->Draw();
      TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
      //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
      lineMinusSigma->SetLineColor(kBlue);
      lineMinusSigma->SetLineStyle(2 + i);
      lineMinusSigma->Draw();
      legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
    }
  }
  if (plotLTM) {
    // draw LTM
    Double_t ltmRms = 0;
    Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
    TLine* line = new TLine(ltm, 0, ltm, maxY);
    //fListOfObjectsToBeDeleted->Add(line);
    line->SetLineColor(kGreen+2);
    line->SetLineWidth(2);
    line->SetLineStyle(1);
    line->Draw();
    legend->AddEntry(line, Form("LTM: %f", ltm), "l");
    // draw sigma lines
    for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
      TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
      //fListOfObjectsToBeDeleted->Add(linePlusSigma);
      linePlusSigma->SetLineColor(kGreen+2);
      linePlusSigma->SetLineStyle(2+i);
      linePlusSigma->Draw();

      TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
      //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
      lineMinusSigma->SetLineColor(kGreen+2);
      lineMinusSigma->SetLineStyle(2+i);
      lineMinusSigma->Draw();
      legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms)), "l");
    }
  }
  if (!plotMean && !plotMedian && !plotLTM) return -1;
  legend->Draw();
  gStyle->SetOptStat(oldOptStat);
  return 1;
}

//_____________________________________________________________________________
Int_t AliBaseCalibViewer::SigmaCut(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
    Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, 
    const Char_t *sigmas, Float_t sigmaStep) const {
  //
  // Creates a histogram, where you can see, how much of the data are inside sigma-intervals 
  // around the mean/median/LTM
  // with drawCommand, sector and cuts you specify your input data, see EasyDraw
  // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
  // sigmaStep: the binsize of the generated histogram
  // plotMean/plotMedian/plotLTM: specifies where to put the center
  //

  Double_t ltmFraction = 0.8;

  TString drawStr(drawCommand);
  Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
  if (dangerousToDraw) {
    Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
    return -1;
  }
  drawStr += " >> tempHist";

  Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
  TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
  // FIXME is this histogram deleted automatically?
  Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers

  Double_t mean = TMath::Mean(entries, values);
  Double_t median = TMath::Median(entries, values);
  Double_t sigma = TMath::RMS(entries, values);

  TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
  //fListOfObjectsToBeDeleted->Add(legend);
  TH1F *cutHistoMean = 0;
  TH1F *cutHistoMedian = 0;
  TH1F *cutHistoLTM = 0;

  TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
  TVectorF nsigma(sigmasTokens->GetEntriesFast());
  for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
    TString str(((TObjString*)sigmasTokens->At(i))->GetString());
    Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
    nsigma[i] = sig;
  }
  delete sigmasTokens;
  //
  if (plotMean) {
    cutHistoMean = SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
    if (cutHistoMean) {
      //fListOfObjectsToBeDeleted->Add(cutHistoMean);
      cutHistoMean->SetLineColor(kRed);
      legend->AddEntry(cutHistoMean, "Mean", "l");
      cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      cutHistoMean->Draw();
      DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
    } // if (cutHistoMean)

  }
  if (plotMedian) {
    cutHistoMedian = SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
    if (cutHistoMedian) {
      //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
      cutHistoMedian->SetLineColor(kBlue);
      legend->AddEntry(cutHistoMedian, "Median", "l");
      cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
      else cutHistoMedian->Draw();
      DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
    }  // if (cutHistoMedian)
  }
  if (plotLTM) {
    Double_t ltmRms = 0;
    Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
    cutHistoLTM = SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
    if (cutHistoLTM) {
      //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
      cutHistoLTM->SetLineColor(kGreen+2);
      legend->AddEntry(cutHistoLTM, "LTM", "l");
      cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      if ((plotMean && cutHistoMean) || (plotMedian && cutHistoMedian)) cutHistoLTM->Draw("same");
      else cutHistoLTM->Draw();
      DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
    }
  }
  if (!plotMean && !plotMedian && !plotLTM) return -1;
  legend->Draw();
  return 1;
}

//_____________________________________________________________________________
Int_t AliBaseCalibViewer::Integrate(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts, 
    Float_t /*sigmaMax*/, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, 
    const Char_t *sigmas, Float_t /*sigmaStep*/) const {
  //
  // 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"   
  // "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
  // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
  // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
  // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
  // The actual work is done on the array.
  /* Begin_Latex 
     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
     End_Latex  
     */

  Double_t ltmFraction = 0.8;

  TString drawStr(drawCommand);
  Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
  if (dangerousToDraw) {
    Warning("Integrate", "The draw string must not contain ':' or '>>'.");
    return -1;
  }
  drawStr += " >> tempHist";

  Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
  TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
  TGraph *integralGraphMean   = 0;
  TGraph *integralGraphMedian = 0;
  TGraph *integralGraphLTM    = 0;
  Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
  Int_t    *index  = new Int_t[entries];
  Float_t  *xarray = new Float_t[entries];
  Float_t  *yarray = new Float_t[entries];
  TMath::Sort(entries, values, index, kFALSE);

  Double_t mean = TMath::Mean(entries, values);
  Double_t median = TMath::Median(entries, values);
  Double_t sigma = TMath::RMS(entries, values);

  // parse sigmas string
  TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
  TVectorF nsigma(sigmasTokens->GetEntriesFast());
  for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
    TString str(((TObjString*)sigmasTokens->At(i))->GetString());
    Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
    nsigma[i] = sig;
  }
  delete sigmasTokens;
  TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
  //fListOfObjectsToBeDeleted->Add(legend);

  if (plotMean) {
    for (Int_t i = 0; i < entries; i++) {
      xarray[i] = (values[index[i]] - mean) / sigma; 
      yarray[i] = float(i) / float(entries);
    }
    integralGraphMean = new TGraph(entries, xarray, yarray);
    if (integralGraphMean) {
      //fListOfObjectsToBeDeleted->Add(integralGraphMean);
      integralGraphMean->SetLineColor(kRed);
      legend->AddEntry(integralGraphMean, "Mean", "l");
      integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      integralGraphMean->Draw("alu");
      DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
    }
  }
  if (plotMedian) {
    for (Int_t i = 0; i < entries; i++) {
      xarray[i] = (values[index[i]] - median) / sigma; 
      yarray[i] = float(i) / float(entries);
    }
    integralGraphMedian = new TGraph(entries, xarray, yarray);
    if (integralGraphMedian) {
      //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
      integralGraphMedian->SetLineColor(kBlue);
      legend->AddEntry(integralGraphMedian, "Median", "l");
      integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
      else integralGraphMedian->Draw("alu");
      DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
    }
  }
  if (plotLTM) {
    Double_t ltmRms = 0;
    Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
    for (Int_t i = 0; i < entries; i++) {
      xarray[i] = (values[index[i]] - ltm) / ltmRms; 
      yarray[i] = float(i) / float(entries);
    }
    integralGraphLTM = new TGraph(entries, xarray, yarray);
    if (integralGraphLTM) {
      //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
      integralGraphLTM->SetLineColor(kGreen+2);
      legend->AddEntry(integralGraphLTM, "LTM", "l");
      integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
      if ((plotMean && integralGraphMean) || (plotMedian && integralGraphMedian)) integralGraphLTM->Draw("samelu");
      else integralGraphLTM->Draw("alu");
      DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
    }
  }
  delete [] index;
  delete [] xarray;
  delete [] yarray;
  if (!plotMean && !plotMedian && !plotLTM) return -1;
  legend->Draw();
  return entries;
}

//_____________________________________________________________________________
TH1F* AliBaseCalibViewer::Integrate(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
  //
  // 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"   
  // "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
  // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
  // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
  // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
  // The actual work is done on the array.
  /* Begin_Latex 
     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
     End_Latex  
     Begin_Macro(source)
     {
     Float_t mean = 0;
     Float_t sigma = 1.5;
     Float_t sigmaMax = 4;
     gROOT->SetStyle("Plain");
     TH1F *distribution = new TH1F("Distribution2", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
     TRandom rand(23);
     for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
     Float_t *ar = distribution->GetArray();

     TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_Integrate", "", 350, 350);
     macro_example_canvas->Divide(0,2);
     TVirtualPad *pad1 = macro_example_canvas->cd(1);
     pad1->SetGridy();
     pad1->SetGridx();
     distribution->Draw();
     TVirtualPad *pad2 = macro_example_canvas->cd(2);
     pad2->SetGridy();
     pad2->SetGridx();
     TH1F *shist = AliTPCCalibViewer::Integrate(distribution, mean, sigma, sigmaMax);
     shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
     shist->Draw();  

     return macro_example_canvas_Integrate;
     }  
     End_Macro
     */ 


  Float_t *array = histogram->GetArray();
  Int_t    nbins = histogram->GetXaxis()->GetNbins();
  Float_t binLow = histogram->GetXaxis()->GetXmin();
  Float_t binUp  = histogram->GetXaxis()->GetXmax();
  return Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
}

//_____________________________________________________________________________
TH1F* AliBaseCalibViewer::Integrate(Int_t n, Float_t *array, Int_t nbins, Float_t binLow, Float_t binUp, 
    Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
  // 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"   
  // "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
  // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
  // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
  // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
  // Here the actual work is done.

  Bool_t givenUnits = kTRUE;
  if (TMath::Abs(sigma) < 1.e-10 && TMath::Abs(sigmaMax) < 1.e-10) givenUnits = kFALSE;
  if (givenUnits) {
    sigma = 1;
    sigmaMax = (binUp - binLow) / 2.;
  }

  Float_t binWidth = (binUp-binLow)/(nbins - 1);
  if (sigmaStep <= 0) sigmaStep = binWidth;
  Int_t kbins =  (Int_t)(sigmaMax * sigma / sigmaStep) + 1;  // + 1  due to overflow bin in histograms
  Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
  Float_t kbinUp  = givenUnits ? binUp  : sigmaMax;
  TH1F *hist = 0; 
  if (givenUnits)  hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp); 
  if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
  hist->SetDirectory(0);
  hist->Reset();

  // calculate normalization
  //  printf("calculating normalization, integrating from bin 1 to %i \n", n);
  Double_t normalization = 0;
  for (Int_t i = 1; i <= n; i++) {
    normalization += array[i];
  }
  //  printf("normalization: %f \n", normalization);

  // given units: units from given histogram
  // sigma units: in units of sigma
  // iDelta: integrate in interval (mean +- iDelta), given units
  // x:      ofset from mean for integration, given units
  // hist:   needs 

  // fill histogram
  for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
    // integrate array
    Double_t value = 0;
    for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
      value += (x <= binUp && x >= binLow)  ? array[GetBin(x, nbins, binLow, binUp)] : 0;
    }
    if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", value, normalization);
    if (value / normalization > 100) return hist;
    Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
    //  printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
    //  printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
    value = (value / normalization);
    hist->SetBinContent(bin, value);
  }
  return hist;
}

//_____________________________________________________________________________
void AliBaseCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
  //
  // Private function for SigmaCut(...) and Integrate(...)
  // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
  //

  // start to draw the lines, loop over requested sigmas
  for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
    if (!pm) {
      Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
      TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineUp);
      lineUp->SetLineColor(color);
      lineUp->SetLineStyle(2 + i);
      lineUp->Draw();
      TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineLeft);
      lineLeft->SetLineColor(color);
      lineLeft->SetLineStyle(2 + i);
      lineLeft->Draw();
      legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
    }
    else { // if (pm)
      Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
      TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineUp1);
      lineUp1->SetLineColor(color);
      lineUp1->SetLineStyle(2 + i);
      lineUp1->Draw();
      TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineLeft1);
      lineLeft1->SetLineColor(color);
      lineLeft1->SetLineStyle(2 + i);
      lineLeft1->Draw();
      legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
      bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
      TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineUp2);
      lineUp2->SetLineColor(color);
      lineUp2->SetLineStyle(2 + i);
      lineUp2->Draw();
      TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
      //fListOfObjectsToBeDeleted->Add(lineLeft2);
      lineLeft2->SetLineColor(color);
      lineLeft2->SetLineStyle(2 + i);
      lineLeft2->Draw();
      legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
    }
  }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
}

//_____________________________________________________________________________
void AliBaseCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
  //
  // Private function for SigmaCut(...) and Integrate(...)
  // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
  //

  // start to draw the lines, loop over requested sigmas
  for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
    if (!pm) {
      TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineUp);
      lineUp->SetLineColor(color);
      lineUp->SetLineStyle(2 + i);
      lineUp->Draw();
      TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineLeft);
      lineLeft->SetLineColor(color);
      lineLeft->SetLineStyle(2 + i);
      lineLeft->Draw();
      legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
    }
    else { // if (pm)
      TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineUp1);
      lineUp1->SetLineColor(color);
      lineUp1->SetLineStyle(2 + i);
      lineUp1->Draw();
      TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineLeft1);
      lineLeft1->SetLineColor(color);
      lineLeft1->SetLineStyle(2 + i);
      lineLeft1->Draw();
      legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
      TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineUp2);
      lineUp2->SetLineColor(color);
      lineUp2->SetLineStyle(2 + i);
      lineUp2->Draw();
      TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
      //fListOfObjectsToBeDeleted->Add(lineLeft2);
      lineLeft2->SetLineColor(color);
      lineLeft2->SetLineStyle(2 + i);
      lineLeft2->Draw();
      legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i])), "l");
    }
  }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
}
 AliBaseCalibViewer.cxx:1
 AliBaseCalibViewer.cxx:2
 AliBaseCalibViewer.cxx:3
 AliBaseCalibViewer.cxx:4
 AliBaseCalibViewer.cxx:5
 AliBaseCalibViewer.cxx:6
 AliBaseCalibViewer.cxx:7
 AliBaseCalibViewer.cxx:8
 AliBaseCalibViewer.cxx:9
 AliBaseCalibViewer.cxx:10
 AliBaseCalibViewer.cxx:11
 AliBaseCalibViewer.cxx:12
 AliBaseCalibViewer.cxx:13
 AliBaseCalibViewer.cxx:14
 AliBaseCalibViewer.cxx:15
 AliBaseCalibViewer.cxx:16
 AliBaseCalibViewer.cxx:17
 AliBaseCalibViewer.cxx:18
 AliBaseCalibViewer.cxx:19
 AliBaseCalibViewer.cxx:20
 AliBaseCalibViewer.cxx:21
 AliBaseCalibViewer.cxx:22
 AliBaseCalibViewer.cxx:23
 AliBaseCalibViewer.cxx:24
 AliBaseCalibViewer.cxx:25
 AliBaseCalibViewer.cxx:26
 AliBaseCalibViewer.cxx:27
 AliBaseCalibViewer.cxx:28
 AliBaseCalibViewer.cxx:29
 AliBaseCalibViewer.cxx:30
 AliBaseCalibViewer.cxx:31
 AliBaseCalibViewer.cxx:32
 AliBaseCalibViewer.cxx:33
 AliBaseCalibViewer.cxx:34
 AliBaseCalibViewer.cxx:35
 AliBaseCalibViewer.cxx:36
 AliBaseCalibViewer.cxx:37
 AliBaseCalibViewer.cxx:38
 AliBaseCalibViewer.cxx:39
 AliBaseCalibViewer.cxx:40
 AliBaseCalibViewer.cxx:41
 AliBaseCalibViewer.cxx:42
 AliBaseCalibViewer.cxx:43
 AliBaseCalibViewer.cxx:44
 AliBaseCalibViewer.cxx:45
 AliBaseCalibViewer.cxx:46
 AliBaseCalibViewer.cxx:47
 AliBaseCalibViewer.cxx:48
 AliBaseCalibViewer.cxx:49
 AliBaseCalibViewer.cxx:50
 AliBaseCalibViewer.cxx:51
 AliBaseCalibViewer.cxx:52
 AliBaseCalibViewer.cxx:53
 AliBaseCalibViewer.cxx:54
 AliBaseCalibViewer.cxx:55
 AliBaseCalibViewer.cxx:56
 AliBaseCalibViewer.cxx:57
 AliBaseCalibViewer.cxx:58
 AliBaseCalibViewer.cxx:59
 AliBaseCalibViewer.cxx:60
 AliBaseCalibViewer.cxx:61
 AliBaseCalibViewer.cxx:62
 AliBaseCalibViewer.cxx:63
 AliBaseCalibViewer.cxx:64
 AliBaseCalibViewer.cxx:65
 AliBaseCalibViewer.cxx:66
 AliBaseCalibViewer.cxx:67
 AliBaseCalibViewer.cxx:68
 AliBaseCalibViewer.cxx:69
 AliBaseCalibViewer.cxx:70
 AliBaseCalibViewer.cxx:71
 AliBaseCalibViewer.cxx:72
 AliBaseCalibViewer.cxx:73
 AliBaseCalibViewer.cxx:74
 AliBaseCalibViewer.cxx:75
 AliBaseCalibViewer.cxx:76
 AliBaseCalibViewer.cxx:77
 AliBaseCalibViewer.cxx:78
 AliBaseCalibViewer.cxx:79
 AliBaseCalibViewer.cxx:80
 AliBaseCalibViewer.cxx:81
 AliBaseCalibViewer.cxx:82
 AliBaseCalibViewer.cxx:83
 AliBaseCalibViewer.cxx:84
 AliBaseCalibViewer.cxx:85
 AliBaseCalibViewer.cxx:86
 AliBaseCalibViewer.cxx:87
 AliBaseCalibViewer.cxx:88
 AliBaseCalibViewer.cxx:89
 AliBaseCalibViewer.cxx:90
 AliBaseCalibViewer.cxx:91
 AliBaseCalibViewer.cxx:92
 AliBaseCalibViewer.cxx:93
 AliBaseCalibViewer.cxx:94
 AliBaseCalibViewer.cxx:95
 AliBaseCalibViewer.cxx:96
 AliBaseCalibViewer.cxx:97
 AliBaseCalibViewer.cxx:98
 AliBaseCalibViewer.cxx:99
 AliBaseCalibViewer.cxx:100
 AliBaseCalibViewer.cxx:101
 AliBaseCalibViewer.cxx:102
 AliBaseCalibViewer.cxx:103
 AliBaseCalibViewer.cxx:104
 AliBaseCalibViewer.cxx:105
 AliBaseCalibViewer.cxx:106
 AliBaseCalibViewer.cxx:107
 AliBaseCalibViewer.cxx:108
 AliBaseCalibViewer.cxx:109
 AliBaseCalibViewer.cxx:110
 AliBaseCalibViewer.cxx:111
 AliBaseCalibViewer.cxx:112
 AliBaseCalibViewer.cxx:113
 AliBaseCalibViewer.cxx:114
 AliBaseCalibViewer.cxx:115
 AliBaseCalibViewer.cxx:116
 AliBaseCalibViewer.cxx:117
 AliBaseCalibViewer.cxx:118
 AliBaseCalibViewer.cxx:119
 AliBaseCalibViewer.cxx:120
 AliBaseCalibViewer.cxx:121
 AliBaseCalibViewer.cxx:122
 AliBaseCalibViewer.cxx:123
 AliBaseCalibViewer.cxx:124
 AliBaseCalibViewer.cxx:125
 AliBaseCalibViewer.cxx:126
 AliBaseCalibViewer.cxx:127
 AliBaseCalibViewer.cxx:128
 AliBaseCalibViewer.cxx:129
 AliBaseCalibViewer.cxx:130
 AliBaseCalibViewer.cxx:131
 AliBaseCalibViewer.cxx:132
 AliBaseCalibViewer.cxx:133
 AliBaseCalibViewer.cxx:134
 AliBaseCalibViewer.cxx:135
 AliBaseCalibViewer.cxx:136
 AliBaseCalibViewer.cxx:137
 AliBaseCalibViewer.cxx:138
 AliBaseCalibViewer.cxx:139
 AliBaseCalibViewer.cxx:140
 AliBaseCalibViewer.cxx:141
 AliBaseCalibViewer.cxx:142
 AliBaseCalibViewer.cxx:143
 AliBaseCalibViewer.cxx:144
 AliBaseCalibViewer.cxx:145
 AliBaseCalibViewer.cxx:146
 AliBaseCalibViewer.cxx:147
 AliBaseCalibViewer.cxx:148
 AliBaseCalibViewer.cxx:149
 AliBaseCalibViewer.cxx:150
 AliBaseCalibViewer.cxx:151
 AliBaseCalibViewer.cxx:152
 AliBaseCalibViewer.cxx:153
 AliBaseCalibViewer.cxx:154
 AliBaseCalibViewer.cxx:155
 AliBaseCalibViewer.cxx:156
 AliBaseCalibViewer.cxx:157
 AliBaseCalibViewer.cxx:158
 AliBaseCalibViewer.cxx:159
 AliBaseCalibViewer.cxx:160
 AliBaseCalibViewer.cxx:161
 AliBaseCalibViewer.cxx:162
 AliBaseCalibViewer.cxx:163
 AliBaseCalibViewer.cxx:164
 AliBaseCalibViewer.cxx:165
 AliBaseCalibViewer.cxx:166
 AliBaseCalibViewer.cxx:167
 AliBaseCalibViewer.cxx:168
 AliBaseCalibViewer.cxx:169
 AliBaseCalibViewer.cxx:170
 AliBaseCalibViewer.cxx:171
 AliBaseCalibViewer.cxx:172
 AliBaseCalibViewer.cxx:173
 AliBaseCalibViewer.cxx:174
 AliBaseCalibViewer.cxx:175
 AliBaseCalibViewer.cxx:176
 AliBaseCalibViewer.cxx:177
 AliBaseCalibViewer.cxx:178
 AliBaseCalibViewer.cxx:179
 AliBaseCalibViewer.cxx:180
 AliBaseCalibViewer.cxx:181
 AliBaseCalibViewer.cxx:182
 AliBaseCalibViewer.cxx:183
 AliBaseCalibViewer.cxx:184
 AliBaseCalibViewer.cxx:185
 AliBaseCalibViewer.cxx:186
 AliBaseCalibViewer.cxx:187
 AliBaseCalibViewer.cxx:188
 AliBaseCalibViewer.cxx:189
 AliBaseCalibViewer.cxx:190
 AliBaseCalibViewer.cxx:191
 AliBaseCalibViewer.cxx:192
 AliBaseCalibViewer.cxx:193
 AliBaseCalibViewer.cxx:194
 AliBaseCalibViewer.cxx:195
 AliBaseCalibViewer.cxx:196
 AliBaseCalibViewer.cxx:197
 AliBaseCalibViewer.cxx:198
 AliBaseCalibViewer.cxx:199
 AliBaseCalibViewer.cxx:200
 AliBaseCalibViewer.cxx:201
 AliBaseCalibViewer.cxx:202
 AliBaseCalibViewer.cxx:203
 AliBaseCalibViewer.cxx:204
 AliBaseCalibViewer.cxx:205
 AliBaseCalibViewer.cxx:206
 AliBaseCalibViewer.cxx:207
 AliBaseCalibViewer.cxx:208
 AliBaseCalibViewer.cxx:209
 AliBaseCalibViewer.cxx:210
 AliBaseCalibViewer.cxx:211
 AliBaseCalibViewer.cxx:212
 AliBaseCalibViewer.cxx:213
 AliBaseCalibViewer.cxx:214
 AliBaseCalibViewer.cxx:215
 AliBaseCalibViewer.cxx:216
 AliBaseCalibViewer.cxx:217
 AliBaseCalibViewer.cxx:218
 AliBaseCalibViewer.cxx:219
 AliBaseCalibViewer.cxx:220
 AliBaseCalibViewer.cxx:221
 AliBaseCalibViewer.cxx:222
 AliBaseCalibViewer.cxx:223
 AliBaseCalibViewer.cxx:224
 AliBaseCalibViewer.cxx:225
 AliBaseCalibViewer.cxx:226
 AliBaseCalibViewer.cxx:227
 AliBaseCalibViewer.cxx:228
 AliBaseCalibViewer.cxx:229
 AliBaseCalibViewer.cxx:230
 AliBaseCalibViewer.cxx:231
 AliBaseCalibViewer.cxx:232
 AliBaseCalibViewer.cxx:233
 AliBaseCalibViewer.cxx:234
 AliBaseCalibViewer.cxx:235
 AliBaseCalibViewer.cxx:236
 AliBaseCalibViewer.cxx:237
 AliBaseCalibViewer.cxx:238
 AliBaseCalibViewer.cxx:239
 AliBaseCalibViewer.cxx:240
 AliBaseCalibViewer.cxx:241
 AliBaseCalibViewer.cxx:242
 AliBaseCalibViewer.cxx:243
 AliBaseCalibViewer.cxx:244
 AliBaseCalibViewer.cxx:245
 AliBaseCalibViewer.cxx:246
 AliBaseCalibViewer.cxx:247
 AliBaseCalibViewer.cxx:248
 AliBaseCalibViewer.cxx:249
 AliBaseCalibViewer.cxx:250
 AliBaseCalibViewer.cxx:251
 AliBaseCalibViewer.cxx:252
 AliBaseCalibViewer.cxx:253
 AliBaseCalibViewer.cxx:254
 AliBaseCalibViewer.cxx:255
 AliBaseCalibViewer.cxx:256
 AliBaseCalibViewer.cxx:257
 AliBaseCalibViewer.cxx:258
 AliBaseCalibViewer.cxx:259
 AliBaseCalibViewer.cxx:260
 AliBaseCalibViewer.cxx:261
 AliBaseCalibViewer.cxx:262
 AliBaseCalibViewer.cxx:263
 AliBaseCalibViewer.cxx:264
 AliBaseCalibViewer.cxx:265
 AliBaseCalibViewer.cxx:266
 AliBaseCalibViewer.cxx:267
 AliBaseCalibViewer.cxx:268
 AliBaseCalibViewer.cxx:269
 AliBaseCalibViewer.cxx:270
 AliBaseCalibViewer.cxx:271
 AliBaseCalibViewer.cxx:272
 AliBaseCalibViewer.cxx:273
 AliBaseCalibViewer.cxx:274
 AliBaseCalibViewer.cxx:275
 AliBaseCalibViewer.cxx:276
 AliBaseCalibViewer.cxx:277
 AliBaseCalibViewer.cxx:278
 AliBaseCalibViewer.cxx:279
 AliBaseCalibViewer.cxx:280
 AliBaseCalibViewer.cxx:281
 AliBaseCalibViewer.cxx:282
 AliBaseCalibViewer.cxx:283
 AliBaseCalibViewer.cxx:284
 AliBaseCalibViewer.cxx:285
 AliBaseCalibViewer.cxx:286
 AliBaseCalibViewer.cxx:287
 AliBaseCalibViewer.cxx:288
 AliBaseCalibViewer.cxx:289
 AliBaseCalibViewer.cxx:290
 AliBaseCalibViewer.cxx:291
 AliBaseCalibViewer.cxx:292
 AliBaseCalibViewer.cxx:293
 AliBaseCalibViewer.cxx:294
 AliBaseCalibViewer.cxx:295
 AliBaseCalibViewer.cxx:296
 AliBaseCalibViewer.cxx:297
 AliBaseCalibViewer.cxx:298
 AliBaseCalibViewer.cxx:299
 AliBaseCalibViewer.cxx:300
 AliBaseCalibViewer.cxx:301
 AliBaseCalibViewer.cxx:302
 AliBaseCalibViewer.cxx:303
 AliBaseCalibViewer.cxx:304
 AliBaseCalibViewer.cxx:305
 AliBaseCalibViewer.cxx:306
 AliBaseCalibViewer.cxx:307
 AliBaseCalibViewer.cxx:308
 AliBaseCalibViewer.cxx:309
 AliBaseCalibViewer.cxx:310
 AliBaseCalibViewer.cxx:311
 AliBaseCalibViewer.cxx:312
 AliBaseCalibViewer.cxx:313
 AliBaseCalibViewer.cxx:314
 AliBaseCalibViewer.cxx:315
 AliBaseCalibViewer.cxx:316
 AliBaseCalibViewer.cxx:317
 AliBaseCalibViewer.cxx:318
 AliBaseCalibViewer.cxx:319
 AliBaseCalibViewer.cxx:320
 AliBaseCalibViewer.cxx:321
 AliBaseCalibViewer.cxx:322
 AliBaseCalibViewer.cxx:323
 AliBaseCalibViewer.cxx:324
 AliBaseCalibViewer.cxx:325
 AliBaseCalibViewer.cxx:326
 AliBaseCalibViewer.cxx:327
 AliBaseCalibViewer.cxx:328
 AliBaseCalibViewer.cxx:329
 AliBaseCalibViewer.cxx:330
 AliBaseCalibViewer.cxx:331
 AliBaseCalibViewer.cxx:332
 AliBaseCalibViewer.cxx:333
 AliBaseCalibViewer.cxx:334
 AliBaseCalibViewer.cxx:335
 AliBaseCalibViewer.cxx:336
 AliBaseCalibViewer.cxx:337
 AliBaseCalibViewer.cxx:338
 AliBaseCalibViewer.cxx:339
 AliBaseCalibViewer.cxx:340
 AliBaseCalibViewer.cxx:341
 AliBaseCalibViewer.cxx:342
 AliBaseCalibViewer.cxx:343
 AliBaseCalibViewer.cxx:344
 AliBaseCalibViewer.cxx:345
 AliBaseCalibViewer.cxx:346
 AliBaseCalibViewer.cxx:347
 AliBaseCalibViewer.cxx:348
 AliBaseCalibViewer.cxx:349
 AliBaseCalibViewer.cxx:350
 AliBaseCalibViewer.cxx:351
 AliBaseCalibViewer.cxx:352
 AliBaseCalibViewer.cxx:353
 AliBaseCalibViewer.cxx:354
 AliBaseCalibViewer.cxx:355
 AliBaseCalibViewer.cxx:356
 AliBaseCalibViewer.cxx:357
 AliBaseCalibViewer.cxx:358
 AliBaseCalibViewer.cxx:359
 AliBaseCalibViewer.cxx:360
 AliBaseCalibViewer.cxx:361
 AliBaseCalibViewer.cxx:362
 AliBaseCalibViewer.cxx:363
 AliBaseCalibViewer.cxx:364
 AliBaseCalibViewer.cxx:365
 AliBaseCalibViewer.cxx:366
 AliBaseCalibViewer.cxx:367
 AliBaseCalibViewer.cxx:368
 AliBaseCalibViewer.cxx:369
 AliBaseCalibViewer.cxx:370
 AliBaseCalibViewer.cxx:371
 AliBaseCalibViewer.cxx:372
 AliBaseCalibViewer.cxx:373
 AliBaseCalibViewer.cxx:374
 AliBaseCalibViewer.cxx:375
 AliBaseCalibViewer.cxx:376
 AliBaseCalibViewer.cxx:377
 AliBaseCalibViewer.cxx:378
 AliBaseCalibViewer.cxx:379
 AliBaseCalibViewer.cxx:380
 AliBaseCalibViewer.cxx:381
 AliBaseCalibViewer.cxx:382
 AliBaseCalibViewer.cxx:383
 AliBaseCalibViewer.cxx:384
 AliBaseCalibViewer.cxx:385
 AliBaseCalibViewer.cxx:386
 AliBaseCalibViewer.cxx:387
 AliBaseCalibViewer.cxx:388
 AliBaseCalibViewer.cxx:389
 AliBaseCalibViewer.cxx:390
 AliBaseCalibViewer.cxx:391
 AliBaseCalibViewer.cxx:392
 AliBaseCalibViewer.cxx:393
 AliBaseCalibViewer.cxx:394
 AliBaseCalibViewer.cxx:395
 AliBaseCalibViewer.cxx:396
 AliBaseCalibViewer.cxx:397
 AliBaseCalibViewer.cxx:398
 AliBaseCalibViewer.cxx:399
 AliBaseCalibViewer.cxx:400
 AliBaseCalibViewer.cxx:401
 AliBaseCalibViewer.cxx:402
 AliBaseCalibViewer.cxx:403
 AliBaseCalibViewer.cxx:404
 AliBaseCalibViewer.cxx:405
 AliBaseCalibViewer.cxx:406
 AliBaseCalibViewer.cxx:407
 AliBaseCalibViewer.cxx:408
 AliBaseCalibViewer.cxx:409
 AliBaseCalibViewer.cxx:410
 AliBaseCalibViewer.cxx:411
 AliBaseCalibViewer.cxx:412
 AliBaseCalibViewer.cxx:413
 AliBaseCalibViewer.cxx:414
 AliBaseCalibViewer.cxx:415
 AliBaseCalibViewer.cxx:416
 AliBaseCalibViewer.cxx:417
 AliBaseCalibViewer.cxx:418
 AliBaseCalibViewer.cxx:419
 AliBaseCalibViewer.cxx:420
 AliBaseCalibViewer.cxx:421
 AliBaseCalibViewer.cxx:422
 AliBaseCalibViewer.cxx:423
 AliBaseCalibViewer.cxx:424
 AliBaseCalibViewer.cxx:425
 AliBaseCalibViewer.cxx:426
 AliBaseCalibViewer.cxx:427
 AliBaseCalibViewer.cxx:428
 AliBaseCalibViewer.cxx:429
 AliBaseCalibViewer.cxx:430
 AliBaseCalibViewer.cxx:431
 AliBaseCalibViewer.cxx:432
 AliBaseCalibViewer.cxx:433
 AliBaseCalibViewer.cxx:434
 AliBaseCalibViewer.cxx:435
 AliBaseCalibViewer.cxx:436
 AliBaseCalibViewer.cxx:437
 AliBaseCalibViewer.cxx:438
 AliBaseCalibViewer.cxx:439
 AliBaseCalibViewer.cxx:440
 AliBaseCalibViewer.cxx:441
 AliBaseCalibViewer.cxx:442
 AliBaseCalibViewer.cxx:443
 AliBaseCalibViewer.cxx:444
 AliBaseCalibViewer.cxx:445
 AliBaseCalibViewer.cxx:446
 AliBaseCalibViewer.cxx:447
 AliBaseCalibViewer.cxx:448
 AliBaseCalibViewer.cxx:449
 AliBaseCalibViewer.cxx:450
 AliBaseCalibViewer.cxx:451
 AliBaseCalibViewer.cxx:452
 AliBaseCalibViewer.cxx:453
 AliBaseCalibViewer.cxx:454
 AliBaseCalibViewer.cxx:455
 AliBaseCalibViewer.cxx:456
 AliBaseCalibViewer.cxx:457
 AliBaseCalibViewer.cxx:458
 AliBaseCalibViewer.cxx:459
 AliBaseCalibViewer.cxx:460
 AliBaseCalibViewer.cxx:461
 AliBaseCalibViewer.cxx:462
 AliBaseCalibViewer.cxx:463
 AliBaseCalibViewer.cxx:464
 AliBaseCalibViewer.cxx:465
 AliBaseCalibViewer.cxx:466
 AliBaseCalibViewer.cxx:467
 AliBaseCalibViewer.cxx:468
 AliBaseCalibViewer.cxx:469
 AliBaseCalibViewer.cxx:470
 AliBaseCalibViewer.cxx:471
 AliBaseCalibViewer.cxx:472
 AliBaseCalibViewer.cxx:473
 AliBaseCalibViewer.cxx:474
 AliBaseCalibViewer.cxx:475
 AliBaseCalibViewer.cxx:476
 AliBaseCalibViewer.cxx:477
 AliBaseCalibViewer.cxx:478
 AliBaseCalibViewer.cxx:479
 AliBaseCalibViewer.cxx:480
 AliBaseCalibViewer.cxx:481
 AliBaseCalibViewer.cxx:482
 AliBaseCalibViewer.cxx:483
 AliBaseCalibViewer.cxx:484
 AliBaseCalibViewer.cxx:485
 AliBaseCalibViewer.cxx:486
 AliBaseCalibViewer.cxx:487
 AliBaseCalibViewer.cxx:488
 AliBaseCalibViewer.cxx:489
 AliBaseCalibViewer.cxx:490
 AliBaseCalibViewer.cxx:491
 AliBaseCalibViewer.cxx:492
 AliBaseCalibViewer.cxx:493
 AliBaseCalibViewer.cxx:494
 AliBaseCalibViewer.cxx:495
 AliBaseCalibViewer.cxx:496
 AliBaseCalibViewer.cxx:497
 AliBaseCalibViewer.cxx:498
 AliBaseCalibViewer.cxx:499
 AliBaseCalibViewer.cxx:500
 AliBaseCalibViewer.cxx:501
 AliBaseCalibViewer.cxx:502
 AliBaseCalibViewer.cxx:503
 AliBaseCalibViewer.cxx:504
 AliBaseCalibViewer.cxx:505
 AliBaseCalibViewer.cxx:506
 AliBaseCalibViewer.cxx:507
 AliBaseCalibViewer.cxx:508
 AliBaseCalibViewer.cxx:509
 AliBaseCalibViewer.cxx:510
 AliBaseCalibViewer.cxx:511
 AliBaseCalibViewer.cxx:512
 AliBaseCalibViewer.cxx:513
 AliBaseCalibViewer.cxx:514
 AliBaseCalibViewer.cxx:515
 AliBaseCalibViewer.cxx:516
 AliBaseCalibViewer.cxx:517
 AliBaseCalibViewer.cxx:518
 AliBaseCalibViewer.cxx:519
 AliBaseCalibViewer.cxx:520
 AliBaseCalibViewer.cxx:521
 AliBaseCalibViewer.cxx:522
 AliBaseCalibViewer.cxx:523
 AliBaseCalibViewer.cxx:524
 AliBaseCalibViewer.cxx:525
 AliBaseCalibViewer.cxx:526
 AliBaseCalibViewer.cxx:527
 AliBaseCalibViewer.cxx:528
 AliBaseCalibViewer.cxx:529
 AliBaseCalibViewer.cxx:530
 AliBaseCalibViewer.cxx:531
 AliBaseCalibViewer.cxx:532
 AliBaseCalibViewer.cxx:533
 AliBaseCalibViewer.cxx:534
 AliBaseCalibViewer.cxx:535
 AliBaseCalibViewer.cxx:536
 AliBaseCalibViewer.cxx:537
 AliBaseCalibViewer.cxx:538
 AliBaseCalibViewer.cxx:539
 AliBaseCalibViewer.cxx:540
 AliBaseCalibViewer.cxx:541
 AliBaseCalibViewer.cxx:542
 AliBaseCalibViewer.cxx:543
 AliBaseCalibViewer.cxx:544
 AliBaseCalibViewer.cxx:545
 AliBaseCalibViewer.cxx:546
 AliBaseCalibViewer.cxx:547
 AliBaseCalibViewer.cxx:548
 AliBaseCalibViewer.cxx:549
 AliBaseCalibViewer.cxx:550
 AliBaseCalibViewer.cxx:551
 AliBaseCalibViewer.cxx:552
 AliBaseCalibViewer.cxx:553
 AliBaseCalibViewer.cxx:554
 AliBaseCalibViewer.cxx:555
 AliBaseCalibViewer.cxx:556
 AliBaseCalibViewer.cxx:557
 AliBaseCalibViewer.cxx:558
 AliBaseCalibViewer.cxx:559
 AliBaseCalibViewer.cxx:560
 AliBaseCalibViewer.cxx:561
 AliBaseCalibViewer.cxx:562
 AliBaseCalibViewer.cxx:563
 AliBaseCalibViewer.cxx:564
 AliBaseCalibViewer.cxx:565
 AliBaseCalibViewer.cxx:566
 AliBaseCalibViewer.cxx:567
 AliBaseCalibViewer.cxx:568
 AliBaseCalibViewer.cxx:569
 AliBaseCalibViewer.cxx:570
 AliBaseCalibViewer.cxx:571
 AliBaseCalibViewer.cxx:572
 AliBaseCalibViewer.cxx:573
 AliBaseCalibViewer.cxx:574
 AliBaseCalibViewer.cxx:575
 AliBaseCalibViewer.cxx:576
 AliBaseCalibViewer.cxx:577
 AliBaseCalibViewer.cxx:578
 AliBaseCalibViewer.cxx:579
 AliBaseCalibViewer.cxx:580
 AliBaseCalibViewer.cxx:581
 AliBaseCalibViewer.cxx:582
 AliBaseCalibViewer.cxx:583
 AliBaseCalibViewer.cxx:584
 AliBaseCalibViewer.cxx:585
 AliBaseCalibViewer.cxx:586
 AliBaseCalibViewer.cxx:587
 AliBaseCalibViewer.cxx:588
 AliBaseCalibViewer.cxx:589
 AliBaseCalibViewer.cxx:590
 AliBaseCalibViewer.cxx:591
 AliBaseCalibViewer.cxx:592
 AliBaseCalibViewer.cxx:593
 AliBaseCalibViewer.cxx:594
 AliBaseCalibViewer.cxx:595
 AliBaseCalibViewer.cxx:596
 AliBaseCalibViewer.cxx:597
 AliBaseCalibViewer.cxx:598
 AliBaseCalibViewer.cxx:599
 AliBaseCalibViewer.cxx:600
 AliBaseCalibViewer.cxx:601
 AliBaseCalibViewer.cxx:602
 AliBaseCalibViewer.cxx:603
 AliBaseCalibViewer.cxx:604
 AliBaseCalibViewer.cxx:605
 AliBaseCalibViewer.cxx:606
 AliBaseCalibViewer.cxx:607
 AliBaseCalibViewer.cxx:608
 AliBaseCalibViewer.cxx:609
 AliBaseCalibViewer.cxx:610
 AliBaseCalibViewer.cxx:611
 AliBaseCalibViewer.cxx:612
 AliBaseCalibViewer.cxx:613
 AliBaseCalibViewer.cxx:614
 AliBaseCalibViewer.cxx:615
 AliBaseCalibViewer.cxx:616
 AliBaseCalibViewer.cxx:617
 AliBaseCalibViewer.cxx:618
 AliBaseCalibViewer.cxx:619
 AliBaseCalibViewer.cxx:620
 AliBaseCalibViewer.cxx:621
 AliBaseCalibViewer.cxx:622
 AliBaseCalibViewer.cxx:623
 AliBaseCalibViewer.cxx:624
 AliBaseCalibViewer.cxx:625
 AliBaseCalibViewer.cxx:626
 AliBaseCalibViewer.cxx:627
 AliBaseCalibViewer.cxx:628
 AliBaseCalibViewer.cxx:629
 AliBaseCalibViewer.cxx:630
 AliBaseCalibViewer.cxx:631
 AliBaseCalibViewer.cxx:632
 AliBaseCalibViewer.cxx:633
 AliBaseCalibViewer.cxx:634
 AliBaseCalibViewer.cxx:635
 AliBaseCalibViewer.cxx:636
 AliBaseCalibViewer.cxx:637
 AliBaseCalibViewer.cxx:638
 AliBaseCalibViewer.cxx:639
 AliBaseCalibViewer.cxx:640
 AliBaseCalibViewer.cxx:641
 AliBaseCalibViewer.cxx:642
 AliBaseCalibViewer.cxx:643
 AliBaseCalibViewer.cxx:644
 AliBaseCalibViewer.cxx:645
 AliBaseCalibViewer.cxx:646
 AliBaseCalibViewer.cxx:647
 AliBaseCalibViewer.cxx:648
 AliBaseCalibViewer.cxx:649
 AliBaseCalibViewer.cxx:650
 AliBaseCalibViewer.cxx:651
 AliBaseCalibViewer.cxx:652
 AliBaseCalibViewer.cxx:653
 AliBaseCalibViewer.cxx:654
 AliBaseCalibViewer.cxx:655
 AliBaseCalibViewer.cxx:656
 AliBaseCalibViewer.cxx:657
 AliBaseCalibViewer.cxx:658
 AliBaseCalibViewer.cxx:659
 AliBaseCalibViewer.cxx:660
 AliBaseCalibViewer.cxx:661
 AliBaseCalibViewer.cxx:662
 AliBaseCalibViewer.cxx:663
 AliBaseCalibViewer.cxx:664
 AliBaseCalibViewer.cxx:665
 AliBaseCalibViewer.cxx:666
 AliBaseCalibViewer.cxx:667
 AliBaseCalibViewer.cxx:668
 AliBaseCalibViewer.cxx:669
 AliBaseCalibViewer.cxx:670
 AliBaseCalibViewer.cxx:671
 AliBaseCalibViewer.cxx:672
 AliBaseCalibViewer.cxx:673
 AliBaseCalibViewer.cxx:674
 AliBaseCalibViewer.cxx:675
 AliBaseCalibViewer.cxx:676
 AliBaseCalibViewer.cxx:677
 AliBaseCalibViewer.cxx:678
 AliBaseCalibViewer.cxx:679
 AliBaseCalibViewer.cxx:680
 AliBaseCalibViewer.cxx:681
 AliBaseCalibViewer.cxx:682
 AliBaseCalibViewer.cxx:683
 AliBaseCalibViewer.cxx:684
 AliBaseCalibViewer.cxx:685
 AliBaseCalibViewer.cxx:686
 AliBaseCalibViewer.cxx:687
 AliBaseCalibViewer.cxx:688
 AliBaseCalibViewer.cxx:689
 AliBaseCalibViewer.cxx:690
 AliBaseCalibViewer.cxx:691
 AliBaseCalibViewer.cxx:692
 AliBaseCalibViewer.cxx:693
 AliBaseCalibViewer.cxx:694
 AliBaseCalibViewer.cxx:695
 AliBaseCalibViewer.cxx:696
 AliBaseCalibViewer.cxx:697
 AliBaseCalibViewer.cxx:698
 AliBaseCalibViewer.cxx:699
 AliBaseCalibViewer.cxx:700
 AliBaseCalibViewer.cxx:701
 AliBaseCalibViewer.cxx:702
 AliBaseCalibViewer.cxx:703
 AliBaseCalibViewer.cxx:704
 AliBaseCalibViewer.cxx:705
 AliBaseCalibViewer.cxx:706
 AliBaseCalibViewer.cxx:707
 AliBaseCalibViewer.cxx:708
 AliBaseCalibViewer.cxx:709
 AliBaseCalibViewer.cxx:710
 AliBaseCalibViewer.cxx:711
 AliBaseCalibViewer.cxx:712
 AliBaseCalibViewer.cxx:713
 AliBaseCalibViewer.cxx:714
 AliBaseCalibViewer.cxx:715
 AliBaseCalibViewer.cxx:716
 AliBaseCalibViewer.cxx:717
 AliBaseCalibViewer.cxx:718
 AliBaseCalibViewer.cxx:719
 AliBaseCalibViewer.cxx:720
 AliBaseCalibViewer.cxx:721
 AliBaseCalibViewer.cxx:722
 AliBaseCalibViewer.cxx:723
 AliBaseCalibViewer.cxx:724
 AliBaseCalibViewer.cxx:725
 AliBaseCalibViewer.cxx:726
 AliBaseCalibViewer.cxx:727
 AliBaseCalibViewer.cxx:728
 AliBaseCalibViewer.cxx:729
 AliBaseCalibViewer.cxx:730
 AliBaseCalibViewer.cxx:731
 AliBaseCalibViewer.cxx:732
 AliBaseCalibViewer.cxx:733
 AliBaseCalibViewer.cxx:734
 AliBaseCalibViewer.cxx:735
 AliBaseCalibViewer.cxx:736
 AliBaseCalibViewer.cxx:737
 AliBaseCalibViewer.cxx:738
 AliBaseCalibViewer.cxx:739
 AliBaseCalibViewer.cxx:740
 AliBaseCalibViewer.cxx:741
 AliBaseCalibViewer.cxx:742
 AliBaseCalibViewer.cxx:743
 AliBaseCalibViewer.cxx:744
 AliBaseCalibViewer.cxx:745
 AliBaseCalibViewer.cxx:746
 AliBaseCalibViewer.cxx:747
 AliBaseCalibViewer.cxx:748
 AliBaseCalibViewer.cxx:749
 AliBaseCalibViewer.cxx:750
 AliBaseCalibViewer.cxx:751
 AliBaseCalibViewer.cxx:752
 AliBaseCalibViewer.cxx:753
 AliBaseCalibViewer.cxx:754
 AliBaseCalibViewer.cxx:755
 AliBaseCalibViewer.cxx:756
 AliBaseCalibViewer.cxx:757
 AliBaseCalibViewer.cxx:758
 AliBaseCalibViewer.cxx:759
 AliBaseCalibViewer.cxx:760
 AliBaseCalibViewer.cxx:761
 AliBaseCalibViewer.cxx:762
 AliBaseCalibViewer.cxx:763
 AliBaseCalibViewer.cxx:764
 AliBaseCalibViewer.cxx:765
 AliBaseCalibViewer.cxx:766
 AliBaseCalibViewer.cxx:767
 AliBaseCalibViewer.cxx:768
 AliBaseCalibViewer.cxx:769
 AliBaseCalibViewer.cxx:770
 AliBaseCalibViewer.cxx:771
 AliBaseCalibViewer.cxx:772
 AliBaseCalibViewer.cxx:773
 AliBaseCalibViewer.cxx:774
 AliBaseCalibViewer.cxx:775
 AliBaseCalibViewer.cxx:776
 AliBaseCalibViewer.cxx:777
 AliBaseCalibViewer.cxx:778
 AliBaseCalibViewer.cxx:779
 AliBaseCalibViewer.cxx:780
 AliBaseCalibViewer.cxx:781
 AliBaseCalibViewer.cxx:782
 AliBaseCalibViewer.cxx:783
 AliBaseCalibViewer.cxx:784
 AliBaseCalibViewer.cxx:785
 AliBaseCalibViewer.cxx:786
 AliBaseCalibViewer.cxx:787
 AliBaseCalibViewer.cxx:788
 AliBaseCalibViewer.cxx:789
 AliBaseCalibViewer.cxx:790
 AliBaseCalibViewer.cxx:791
 AliBaseCalibViewer.cxx:792
 AliBaseCalibViewer.cxx:793
 AliBaseCalibViewer.cxx:794
 AliBaseCalibViewer.cxx:795
 AliBaseCalibViewer.cxx:796
 AliBaseCalibViewer.cxx:797
 AliBaseCalibViewer.cxx:798
 AliBaseCalibViewer.cxx:799
 AliBaseCalibViewer.cxx:800
 AliBaseCalibViewer.cxx:801
 AliBaseCalibViewer.cxx:802
 AliBaseCalibViewer.cxx:803
 AliBaseCalibViewer.cxx:804
 AliBaseCalibViewer.cxx:805
 AliBaseCalibViewer.cxx:806
 AliBaseCalibViewer.cxx:807
 AliBaseCalibViewer.cxx:808
 AliBaseCalibViewer.cxx:809
 AliBaseCalibViewer.cxx:810
 AliBaseCalibViewer.cxx:811
 AliBaseCalibViewer.cxx:812
 AliBaseCalibViewer.cxx:813
 AliBaseCalibViewer.cxx:814
 AliBaseCalibViewer.cxx:815
 AliBaseCalibViewer.cxx:816
 AliBaseCalibViewer.cxx:817
 AliBaseCalibViewer.cxx:818
 AliBaseCalibViewer.cxx:819
 AliBaseCalibViewer.cxx:820
 AliBaseCalibViewer.cxx:821
 AliBaseCalibViewer.cxx:822
 AliBaseCalibViewer.cxx:823
 AliBaseCalibViewer.cxx:824
 AliBaseCalibViewer.cxx:825
 AliBaseCalibViewer.cxx:826
 AliBaseCalibViewer.cxx:827
 AliBaseCalibViewer.cxx:828
 AliBaseCalibViewer.cxx:829
 AliBaseCalibViewer.cxx:830
 AliBaseCalibViewer.cxx:831
 AliBaseCalibViewer.cxx:832
 AliBaseCalibViewer.cxx:833
 AliBaseCalibViewer.cxx:834
 AliBaseCalibViewer.cxx:835
 AliBaseCalibViewer.cxx:836
 AliBaseCalibViewer.cxx:837
 AliBaseCalibViewer.cxx:838
 AliBaseCalibViewer.cxx:839
 AliBaseCalibViewer.cxx:840
 AliBaseCalibViewer.cxx:841
 AliBaseCalibViewer.cxx:842
 AliBaseCalibViewer.cxx:843
 AliBaseCalibViewer.cxx:844
 AliBaseCalibViewer.cxx:845
 AliBaseCalibViewer.cxx:846
 AliBaseCalibViewer.cxx:847
 AliBaseCalibViewer.cxx:848
 AliBaseCalibViewer.cxx:849
 AliBaseCalibViewer.cxx:850
 AliBaseCalibViewer.cxx:851
 AliBaseCalibViewer.cxx:852
 AliBaseCalibViewer.cxx:853
 AliBaseCalibViewer.cxx:854
 AliBaseCalibViewer.cxx:855
 AliBaseCalibViewer.cxx:856
 AliBaseCalibViewer.cxx:857
 AliBaseCalibViewer.cxx:858
 AliBaseCalibViewer.cxx:859
 AliBaseCalibViewer.cxx:860
 AliBaseCalibViewer.cxx:861
 AliBaseCalibViewer.cxx:862
 AliBaseCalibViewer.cxx:863
 AliBaseCalibViewer.cxx:864
 AliBaseCalibViewer.cxx:865
 AliBaseCalibViewer.cxx:866
 AliBaseCalibViewer.cxx:867
 AliBaseCalibViewer.cxx:868
 AliBaseCalibViewer.cxx:869
 AliBaseCalibViewer.cxx:870
 AliBaseCalibViewer.cxx:871
 AliBaseCalibViewer.cxx:872
 AliBaseCalibViewer.cxx:873
 AliBaseCalibViewer.cxx:874
 AliBaseCalibViewer.cxx:875
 AliBaseCalibViewer.cxx:876
 AliBaseCalibViewer.cxx:877
 AliBaseCalibViewer.cxx:878
 AliBaseCalibViewer.cxx:879
 AliBaseCalibViewer.cxx:880
 AliBaseCalibViewer.cxx:881
 AliBaseCalibViewer.cxx:882
 AliBaseCalibViewer.cxx:883
 AliBaseCalibViewer.cxx:884
 AliBaseCalibViewer.cxx:885
 AliBaseCalibViewer.cxx:886
 AliBaseCalibViewer.cxx:887
 AliBaseCalibViewer.cxx:888
 AliBaseCalibViewer.cxx:889
 AliBaseCalibViewer.cxx:890
 AliBaseCalibViewer.cxx:891
 AliBaseCalibViewer.cxx:892
 AliBaseCalibViewer.cxx:893
 AliBaseCalibViewer.cxx:894
 AliBaseCalibViewer.cxx:895
 AliBaseCalibViewer.cxx:896
 AliBaseCalibViewer.cxx:897
 AliBaseCalibViewer.cxx:898
 AliBaseCalibViewer.cxx:899
 AliBaseCalibViewer.cxx:900
 AliBaseCalibViewer.cxx:901
 AliBaseCalibViewer.cxx:902
 AliBaseCalibViewer.cxx:903
 AliBaseCalibViewer.cxx:904
 AliBaseCalibViewer.cxx:905
 AliBaseCalibViewer.cxx:906
 AliBaseCalibViewer.cxx:907
 AliBaseCalibViewer.cxx:908
 AliBaseCalibViewer.cxx:909
 AliBaseCalibViewer.cxx:910
 AliBaseCalibViewer.cxx:911
 AliBaseCalibViewer.cxx:912
 AliBaseCalibViewer.cxx:913
 AliBaseCalibViewer.cxx:914
 AliBaseCalibViewer.cxx:915
 AliBaseCalibViewer.cxx:916
 AliBaseCalibViewer.cxx:917
 AliBaseCalibViewer.cxx:918
 AliBaseCalibViewer.cxx:919
 AliBaseCalibViewer.cxx:920
 AliBaseCalibViewer.cxx:921
 AliBaseCalibViewer.cxx:922
 AliBaseCalibViewer.cxx:923
 AliBaseCalibViewer.cxx:924
 AliBaseCalibViewer.cxx:925
 AliBaseCalibViewer.cxx:926
 AliBaseCalibViewer.cxx:927
 AliBaseCalibViewer.cxx:928
 AliBaseCalibViewer.cxx:929
 AliBaseCalibViewer.cxx:930
 AliBaseCalibViewer.cxx:931
 AliBaseCalibViewer.cxx:932
 AliBaseCalibViewer.cxx:933
 AliBaseCalibViewer.cxx:934
 AliBaseCalibViewer.cxx:935
 AliBaseCalibViewer.cxx:936
 AliBaseCalibViewer.cxx:937
 AliBaseCalibViewer.cxx:938
 AliBaseCalibViewer.cxx:939
 AliBaseCalibViewer.cxx:940
 AliBaseCalibViewer.cxx:941
 AliBaseCalibViewer.cxx:942
 AliBaseCalibViewer.cxx:943
 AliBaseCalibViewer.cxx:944
 AliBaseCalibViewer.cxx:945
 AliBaseCalibViewer.cxx:946
 AliBaseCalibViewer.cxx:947
 AliBaseCalibViewer.cxx:948
 AliBaseCalibViewer.cxx:949
 AliBaseCalibViewer.cxx:950
 AliBaseCalibViewer.cxx:951
 AliBaseCalibViewer.cxx:952
 AliBaseCalibViewer.cxx:953
 AliBaseCalibViewer.cxx:954
 AliBaseCalibViewer.cxx:955
 AliBaseCalibViewer.cxx:956
 AliBaseCalibViewer.cxx:957
 AliBaseCalibViewer.cxx:958
 AliBaseCalibViewer.cxx:959
 AliBaseCalibViewer.cxx:960
 AliBaseCalibViewer.cxx:961
 AliBaseCalibViewer.cxx:962
 AliBaseCalibViewer.cxx:963
 AliBaseCalibViewer.cxx:964
 AliBaseCalibViewer.cxx:965
 AliBaseCalibViewer.cxx:966
 AliBaseCalibViewer.cxx:967
 AliBaseCalibViewer.cxx:968
 AliBaseCalibViewer.cxx:969
 AliBaseCalibViewer.cxx:970
 AliBaseCalibViewer.cxx:971
 AliBaseCalibViewer.cxx:972
 AliBaseCalibViewer.cxx:973
 AliBaseCalibViewer.cxx:974
 AliBaseCalibViewer.cxx:975
 AliBaseCalibViewer.cxx:976
 AliBaseCalibViewer.cxx:977
 AliBaseCalibViewer.cxx:978
 AliBaseCalibViewer.cxx:979
 AliBaseCalibViewer.cxx:980
 AliBaseCalibViewer.cxx:981
 AliBaseCalibViewer.cxx:982
 AliBaseCalibViewer.cxx:983
 AliBaseCalibViewer.cxx:984
 AliBaseCalibViewer.cxx:985
 AliBaseCalibViewer.cxx:986
 AliBaseCalibViewer.cxx:987
 AliBaseCalibViewer.cxx:988
 AliBaseCalibViewer.cxx:989
 AliBaseCalibViewer.cxx:990
 AliBaseCalibViewer.cxx:991
 AliBaseCalibViewer.cxx:992
 AliBaseCalibViewer.cxx:993
 AliBaseCalibViewer.cxx:994
 AliBaseCalibViewer.cxx:995
 AliBaseCalibViewer.cxx:996
 AliBaseCalibViewer.cxx:997
 AliBaseCalibViewer.cxx:998
 AliBaseCalibViewer.cxx:999
 AliBaseCalibViewer.cxx:1000
 AliBaseCalibViewer.cxx:1001