ROOT logo
/* $Id$ */

// ------------------------------------------------------
//
// Class to handle corrections.
//
// ------------------------------------------------------
//

#include <TFile.h>
#include <TCanvas.h>
#include <TH3F.h>
#include <TH2F.h>
#include <TMath.h>

#include <AliLog.h>
#include "AliCorrectionMatrix2D.h"
#include "AliCorrectionMatrix3D.h"

#include "AliCorrection.h"

//____________________________________________________________________
ClassImp(AliCorrection)

//____________________________________________________________________
AliCorrection::AliCorrection() : TNamed(),
  fEventCorr(0),
  fTrackCorr(0)
{
  // default constructor
}

//____________________________________________________________________
AliCorrection::AliCorrection(const Char_t* name, const Char_t* title, AliPWG0Helper::AnalysisMode analysisMode) : TNamed(name, title),
  fEventCorr(0),
  fTrackCorr(0)
{
  // constructor initializing tnamed

  Float_t* binLimitsPt = 0;
  Int_t nBinsPt = 0;

  // different binnings, better solution could be anticipated...
  if ((analysisMode & AliPWG0Helper::kTPC || analysisMode & AliPWG0Helper::kTPCITS) && (analysisMode & AliPWG0Helper::kFieldOn))
  {
    static Float_t binLimitsPtTmp[] = {0.0, 0.05, 0.1, 0.125, 0.15, 0.175, 0.2, 0.225, 0.25, 0.275, 0.3, 0.325, 0.35, 0.375, 0.4, 0.425, 0.45, 0.475, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.5, 2.0, 5.0, 10.0, 100.0};
    //static Float_t binLimitsPtTmp[] = {0.0, 0.05, 0.1, 0.125, 0.15, 0.175, 0.2, 0.225, 0.25, 0.275, 0.3, 0.325, 0.35, 0.375, 0.4, 0.425, 0.45, 0.475, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 2.0, 5.0, 10.0, 100.0};
    binLimitsPt = (Float_t*) binLimitsPtTmp;
    nBinsPt = 28;
  }
  else if (analysisMode & AliPWG0Helper::kSPD || !(analysisMode & AliPWG0Helper::kFieldOn))
  {
    static Float_t binLimitsPtTmp[] = {-0.5, 100.0};
    binLimitsPt = (Float_t*) binLimitsPtTmp;
    nBinsPt = 1;
  }

  if (!binLimitsPt)
  {
    Printf("FATAL: Invalid binning");
    return;
  }

  // third axis for track histogram
  Int_t nBinsN2 = 1;
  Float_t binLimitsN2[]   = {-0.5, 1000};
  const char* title3 = "Ntracks";
  
  //Int_t nBinsN2 = 21;
  //Float_t binLimitsN2[]   = {-0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, 12.5, 14.5, 16.5, 18.5, 20.5, 25.5, 30.5, 40.5, 50.5, 100.5, 1000.5};
  
  // phi
  //Int_t nBinsN2 = 36;
  //Float_t binLimitsN2[]   = {0.000, 0.175, 0.349, 0.524, 0.698, 0.873, 1.047, 1.222, 1.396, 1.571, 1.745, 1.920, 2.094, 2.269, 2.443, 2.618, 2.793, 2.967, 3.142, 3.316, 3.491, 3.665, 3.840, 4.014, 4.189, 4.363, 4.538, 4.712, 4.887, 5.061, 5.236, 5.411, 5.585, 5.760, 5.934, 6.109, TMath::Pi() * 2};
  //const char* title3 = "#phi";

  // mult axis for event histogram
  Int_t nBinsN = 22;
  Float_t binLimitsN[]   = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, 12.5, 14.5, 16.5, 18.5, 20.5, 25.5, 30.5, 40.5, 50.5, 100.5, 300.5};
  //Int_t nBinsN = 8;
  //Float_t binLimitsN[]   = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 300.5};
  //Int_t nBinsN = 52;
  //Float_t binLimitsN[]   = {-0.5, 0.5, 1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, 11.5, 12.5, 13.5, 14.5, 15.5, 16.5, 17.5, 18.5, 19.5, 20.5, 21.5, 22.5, 23.5, 24.5, 25.5, 26.5, 27.5, 28.5, 29.5, 30.5, 31.5, 32.5, 33.5, 34.5, 35.5, 36.5, 37.5, 38.5, 39.5, 40.5, 41.5, 42.5, 43.5, 44.5, 45.5, 46.5, 47.5, 48.5, 49.5, 50.5, 100.5, 300.5};

  //Float_t binLimitsVtx[] = {-20,-15,-10,-6,-3,0,3,6,10,15,20};
  //Float_t binLimitsVtx[] = {-20,-19,-18,-17,-16,-15,-14,-13,-12,-11,-10,-9,-8,-7,-6,-5,-4,-3,-2,-1,0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20};
  //Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-8,-6,-4,-2,0,2,4,6,8,10,15,20,25,30};
  Float_t binLimitsVtx[] = {-30,-25,-20,-15,-10,-7,-5.5,-4,-3,-2,-1,0,1,2,3,4,5.5,7,10,15,20,25,30};
  /*Float_t binLimitsEta[] = {-3.0,-2.8,-2.6,-2.4,-2.2,
                            -2.0,-1.8,-1.6,-1.4,-1.2,
                            -1.0,-0.8,-0.6,-0.4,-0.2, 0.0,
                             0.2, 0.4, 0.6, 0.8, 1.0,
                             1.2, 1.4, 1.6, 1.8, 2.0,
                             2.2, 2.4, 2.6, 2.8, 3.0};*/
  //Float_t binLimitsEta[] = { -3.0, -2.9, -2.8, -2.7, -2.6, -2.5, -2.4, -2.3, -2.2, -2.1, -2.0, -1.9, -1.8, -1.7, -1.6, -1.5, -1.4, -1.3, -1.2, -1.1, -1.0, -0.9, -0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, -0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 2.7, 2.8, 2.9, 3.0 };
  /*Float_t binLimitsEta[] = {-3,-2.75,-2.5,-2.25,-2.0,
                            -1.75,-1.5,-1.25,
                            -1,-0.75,-0.5,-0.25,
                            0,0.25,0.5,0.75,
                            1.0,1.25,1.5,1.75,
                            2.0,2.25,2.5,2.75,3.0}; */
  //Float_t binLimitsEta[] = {-2.8,-2.4,-2,-1.6,-1.2,-0.8,-0.4,0,0.4,0.8,1.2,1.6,2.0,2.4,2.8}; 
  //Float_t binLimitsEta[] = {-2.8,-2.4,-2,-1.6,-1.2,-0.8,-0.5,0,0.5,0.8,1.2,1.6,2.0,2.4,2.8}; 

  //Float_t binLimitsEta[] = {-3,-2.6,-2.2,-1.8,-1.4,-1,-0.6,-0.2,0.2,0.6,1,1.4,1.8,2.2,2.6,3.0};
  //1. standard
  //  Float_t binLimitsEta[] = {-2,-1.9,-1.8,-1.7,-1.6,-1.5,-1.4,-1.3,-1.2,-1.1,-1.0,-0.9,-0.8,-0.7,-0.6,-0.5,-0.4,-0.3,-0.2,-0.1,0.0,0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0};
  // 2. MB Working Group definition
  Float_t binLimitsEta[] = {-2, -1.8, -1.6, -1.4, -1.2, -1., -0.8, -0.6, -0.4, -0.2, 0.0, 0.2, 0.4, 0.6, 0.8, 1.0, 1.2, 1.4, 1.6, 1.8, 2.0};
//  Float_t binLimitsEta[] = { -7.0, -6.9, -6.8, -6.7, -6.6, -6.5, -6.4, -6.3, -6.2, -6.1, -6.0, -5.9, -5.8, -5.7, -5.6, -5.5, -5.4, -5.3, -5.2, -5.1, -5.0, -4.9, -4.8, -4.7, -4.6, -4.5, -4.4, -4.3, -4.2, -4.1, -4.0, -3.9, -3.8, -3.7, -3.6, -3.5, -3.4, -3.3, -3.2, -3.1, -3.0, -2.9, -2.8, -2.7, -2.6, -2.5, -2.4, -2.3, -2.2, -2.1, -2.0, -1.9, -1.8, -1.7, -1.6, -1.5, -1.4, -1.3, -1.2, -1.1, -1.0, -0.9, -0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, -0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 2.7, 2.8, 2.9, 3.0, 3.1, 3.2, 3.3, 3.4, 3.5, 3.6, 3.7, 3.8, 3.9, 4.0, 4.1, 4.2, 4.3, 4.4, 4.5, 4.6, 4.7, 4.8, 4.9, 5.0, 5.1, 5.2, 5.3, 5.4, 5.5, 5.6, 5.7, 5.8, 5.9, 6.0, 6.1, 6.2, 6.3, 6.4, 6.5, 6.6, 6.7, 6.8, 6.9, 7.0 };

  fEventCorr = new AliCorrectionMatrix2D("EventCorrection", Form("%s EventCorrection", fTitle.Data()), 22, binLimitsVtx, nBinsN, binLimitsN);

  if (analysisMode & AliPWG0Helper::kSPD)
  {
    // 1. standard
    //TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 40, binLimitsEta, nBinsN2, binLimitsN2);
    // 2. MB Working Group definition
    TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 20, binLimitsEta, nBinsN2, binLimitsN2);
    fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
    fTrackCorr->SetAxisTitles("vtx-z (cm)", "#eta", title3);
    delete dummyBinning;
  }
  else
  {
    // 1. standard
    //TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 40, binLimitsEta , nBinsPt, binLimitsPt);
    // 2. MB Working Group definition
    TH3F* dummyBinning = new TH3F("dummyBinning","dummyBinning",22, binLimitsVtx, 20, binLimitsEta , nBinsPt, binLimitsPt);
    fTrackCorr = new AliCorrectionMatrix3D("TrackCorrection", Form("%s TrackCorrection", fTitle.Data()), dummyBinning);
    fTrackCorr->SetAxisTitles("vtx-z (cm)", "#eta", "p_{T} (GeV/c)");
    delete dummyBinning;
  }


  fEventCorr->SetAxisTitles("vtx-z (cm)", "Ntracks");
}

//____________________________________________________________________
AliCorrection::AliCorrection(const AliCorrection& c) : TNamed(c),
  fEventCorr(0),
  fTrackCorr(0)
{
  // copy constructor
  ((AliCorrection &)c).Copy(*this);
}

//____________________________________________________________________
AliCorrection::~AliCorrection()
{
  //
  // destructor
  //

  if (fEventCorr)
  {
    delete fEventCorr;
    fEventCorr = 0;
  }

  if (fTrackCorr)
  {
    delete fTrackCorr;
    fTrackCorr = 0;
  }
}

//____________________________________________________________________
AliCorrection &AliCorrection::operator=(const AliCorrection &c)
{
  // assigment operator

  if (this != &c)
    ((AliCorrection &) c).Copy(*this);

  return *this;
}

//____________________________________________________________________
void AliCorrection::Copy(TObject& c) const
{
  // copy function

  AliCorrection& target = (AliCorrection &) c;

  if (fEventCorr)
    target.fEventCorr = dynamic_cast<AliCorrectionMatrix2D*> (fEventCorr->Clone());

  if (fTrackCorr)
    target.fTrackCorr = dynamic_cast<AliCorrectionMatrix3D*> (fTrackCorr->Clone());
}

//____________________________________________________________________
Long64_t AliCorrection::Merge(TCollection* list)
{
  // Merge a list of AliCorrection objects with this (needed for
  // PROOF). 
  // Returns the number of merged objects (including this).

  if (!list)
    return 0;
  
  if (list->IsEmpty())
    return 1;

  TIterator* iter = list->MakeIterator();
  TObject* obj;

  // collections of measured and generated histograms
  TList* collectionEvent = new TList;
  TList* collectionTrack = new TList;
  
  Int_t count = 0;
  while ((obj = iter->Next())) {
    
    AliCorrection* entry = dynamic_cast<AliCorrection*> (obj);
    if (entry == 0) 
      continue;

    collectionEvent->Add(entry->fEventCorr);
    collectionTrack->Add(entry->fTrackCorr);

    count++;
  }
  fEventCorr->Merge(collectionEvent);
  fTrackCorr->Merge(collectionTrack);

  delete collectionEvent;
  delete collectionTrack;

  return count+1;
}

//____________________________________________________________________
void AliCorrection::Divide()
{
  //
  // divide the histograms to get the correction
  //
  
  if (!fEventCorr || !fTrackCorr)
    return;
    
  fEventCorr->Divide();
  fTrackCorr->Divide();

  Int_t emptyBins = fTrackCorr->CheckEmptyBins(-9.99, 9.99, -0.79, 0.79, 0.2, 4.9);
  printf("INFO: In the central region the track correction of %s has %d empty bins\n", GetName(), emptyBins);
}

//____________________________________________________________________
void AliCorrection::Add(AliCorrection* aCorrectionToAdd, Float_t c)
{
  //
  // add to measured and generated the measured and generated of aCorrectionToAdd
  // with the weight c

  fEventCorr->Add(aCorrectionToAdd->GetEventCorrection(),c);
  fTrackCorr->Add(aCorrectionToAdd->GetTrackCorrection(),c);
}


//____________________________________________________________________
Bool_t AliCorrection::LoadHistograms(const Char_t* dir)
{
  //
  // loads the histograms from a file
  // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
  //

  if (!fEventCorr || !fTrackCorr)
    return kFALSE;

  if (!dir)
    dir = GetName();

  if (!gDirectory->cd(dir))
    return kFALSE;

  Bool_t success = fEventCorr->LoadHistograms();
  success &= fTrackCorr->LoadHistograms();

  gDirectory->cd("..");

  return success;
}

//____________________________________________________________________
void AliCorrection::SaveHistograms()
{
  //
  // saves the histograms in a directory with the name of this object (GetName)
  //
  
  gDirectory->mkdir(GetName());
  gDirectory->cd(GetName());

  if (fEventCorr)
    fEventCorr->SaveHistograms();

  if (fTrackCorr)
    fTrackCorr->SaveHistograms();
    
  gDirectory->cd("..");
}

//____________________________________________________________________
void AliCorrection::ReduceInformation()
{
  // this function deletes the measured and generated histograms to reduce the amount of data
  // in memory

  if (!fEventCorr || !fTrackCorr)
    return;

  fEventCorr->ReduceInformation();
  fTrackCorr->ReduceInformation();
}

//____________________________________________________________________
void AliCorrection::Reset(Option_t* option)
{
  // resets the histograms

  if (fEventCorr)
    fEventCorr->Reset(option);

  if (fTrackCorr)
    fTrackCorr->Reset(option);
}

//____________________________________________________________________
void AliCorrection::DrawHistograms(const Char_t* name)
{
  // draws the corrections

  if (!name)
    name = GetName();

  if (fEventCorr)
    fEventCorr->DrawHistograms(Form("%s event", name));

  if (fTrackCorr)
    fTrackCorr->DrawHistograms(Form("%s track", name));
}

void AliCorrection::DrawOverview(const char* canvasName)
{
  // draw projection of the corrections
  //   to the 3 axis of fTrackCorr and 2 axis of fEventCorr

  TString canvasNameTmp(GetName());
  if (canvasName)
    canvasNameTmp = canvasName;

  TCanvas* canvas = new TCanvas(canvasNameTmp, canvasNameTmp, 1200, 1200);
  canvas->Divide(3, 3);

  if (fTrackCorr) {
    canvas->cd(1);
    fTrackCorr->Get2DCorrectionHistogram("xy", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
    
    canvas->cd(2);
    fTrackCorr->Get2DCorrectionHistogram("yz", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
    
    canvas->cd(3);
    fTrackCorr->Get2DCorrectionHistogram("zx", 0, 0)->DrawCopy("COLZ")->GetZaxis()->SetRangeUser(0, 10);
    
    canvas->cd(4);
    fTrackCorr->Get1DCorrectionHistogram("x", 0.3, 5, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);

    canvas->cd(5);
    fTrackCorr->Get1DCorrectionHistogram("y", 0.3, 5, 0, 0)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);

    canvas->cd(6);
    fTrackCorr->Get1DCorrectionHistogram("z", 0, -1, -1, 1)->DrawCopy()->GetYaxis()->SetRangeUser(0, 10);
  }

  if (fEventCorr)
  {
    canvas->cd(7);
    fEventCorr->GetCorrectionHistogram()->DrawCopy("COLZ")->GetYaxis()->SetRangeUser(0, 30);
    
    canvas->cd(8);
    fEventCorr->Get1DCorrectionHistogram("x")->DrawCopy();

    canvas->cd(9);
    fEventCorr->Get1DCorrectionHistogram("y")->DrawCopy()->GetXaxis()->SetRangeUser(0, 30);
  }
}

//____________________________________________________________________
void AliCorrection::SetCorrectionToUnity()
{
  // set the corrections to unity

  if (fEventCorr)
    fEventCorr->SetCorrectionToUnity();

  if (fTrackCorr)
    fTrackCorr->SetCorrectionToUnity();
}

//____________________________________________________________________
void AliCorrection::Multiply()
{
  // call Multiply

  if (fEventCorr)
  {
    fEventCorr->Multiply();
    // now we manually copy the overflow bin of the y axis (multiplicity) over. This is important to get the event count correct
    TH2* hist = fEventCorr->GetMeasuredHistogram();
    for (Int_t x = 1; x <= hist->GetNbinsX(); ++x)
      fEventCorr->GetGeneratedHistogram()->SetBinContent(x, hist->GetNbinsY() + 1, hist->GetBinContent(x, hist->GetNbinsY() + 1));
  }

  if (fTrackCorr)
    fTrackCorr->Multiply();
}

//____________________________________________________________________
void AliCorrection::Scale(Double_t factor)
{
  // scales the two contained corrections

  fEventCorr->Scale(factor);
  fTrackCorr->Scale(factor);
}

//____________________________________________________________________
void AliCorrection::PrintStats(Float_t zRange, Float_t etaRange, Float_t ptCut)
{
  // prints statistics and effective correction factors

  Printf("AliCorrection::PrintInfo: Values in |eta| < %.2f, |vtx-z| < %.2f and pt > %.2f:", etaRange, zRange, ptCut);

  // prevent to be on bin border
  zRange -= 0.1;
  etaRange -= 0.1;

  TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
  TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();

  TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
  TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();

  Float_t tracksM = measured->Integral(measured->GetXaxis()->FindBin(-zRange), measured->GetXaxis()->FindBin(zRange), measured->GetYaxis()->FindBin(-etaRange), measured->GetYaxis()->FindBin(etaRange), measured->GetZaxis()->FindBin(ptCut), measured->GetZaxis()->GetNbins());
  Float_t tracksG = generated->Integral(generated->GetXaxis()->FindBin(-zRange), generated->GetXaxis()->FindBin(zRange), generated->GetYaxis()->FindBin(-etaRange), generated->GetYaxis()->FindBin(etaRange), generated->GetZaxis()->FindBin(ptCut), generated->GetZaxis()->GetNbins());

  Float_t eventsM = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-zRange), measuredEvents->GetXaxis()->FindBin(zRange), 1, measuredEvents->GetNbinsY());
  Float_t eventsG = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-zRange), generatedEvents->GetXaxis()->FindBin(zRange), 1, generatedEvents->GetNbinsY());

  Float_t eventsM1 = measuredEvents->Integral(measuredEvents->GetXaxis()->FindBin(-zRange), measuredEvents->GetXaxis()->FindBin(zRange), 2, measuredEvents->GetNbinsY());
  Float_t eventsG1 = generatedEvents->Integral(generatedEvents->GetXaxis()->FindBin(-zRange), generatedEvents->GetXaxis()->FindBin(zRange), 2, generatedEvents->GetNbinsY());
  
  Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f, events (N>0) measured: %.1f, events (N>0) generated %.1f", tracksM, tracksG, eventsM, eventsG, eventsM1, eventsG1);

  if (tracksM > 0)
    Printf("Effective average correction factor for TRACKS: %.3f", tracksG / tracksM);
  if (eventsM > 0)
    Printf("Effective average correction factor for EVENTS: %.3f", eventsG / eventsM);

  if (eventsM > 0 && eventsG > 0)
  {
    // normalize to number of events;
    tracksM /= eventsM;
    tracksG /= eventsG;

    Printf("%.2f tracks/event measured, %.2f tracks/event after correction --> effective average correction factor is %.3f (pt cut %.2f GeV/c)", tracksM, tracksG, (tracksM > 0) ? (tracksG / tracksM) : -1, ptCut);
  }
}

//____________________________________________________________________
void AliCorrection::PrintInfo(Float_t ptCut)
{
  // prints some stats

  TH3* measured = GetTrackCorrection()->GetMeasuredHistogram();
  TH3* generated = GetTrackCorrection()->GetGeneratedHistogram();

  TH2* measuredEvents = GetEventCorrection()->GetMeasuredHistogram();
  TH2* generatedEvents = GetEventCorrection()->GetGeneratedHistogram();

  Printf("AliCorrection::PrintInfo: Whole phasespace:");

  Printf("tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->Integral(), generated->Integral(), measuredEvents->Integral(), generatedEvents->Integral());

  Printf("Example centered bin: tracks measured: %.1f tracks generated: %.1f, events measured: %.1f, events generated %.1f", measured->GetBinContent(measured->GetNbinsX() / 2, measured->GetNbinsY() / 2, measured->GetNbinsZ() / 2), generated->GetBinContent(measured->GetNbinsX() / 2, measured->GetNbinsY() / 2, measured->GetNbinsZ() / 2), measuredEvents->GetBinContent(measuredEvents->GetNbinsX() / 2, measuredEvents->GetNbinsY() / 2), generatedEvents->GetBinContent(measuredEvents->GetNbinsX() / 2, measuredEvents->GetNbinsY() / 2));

  PrintStats(10, 0.6, ptCut);
  PrintStats(10, 0.8, ptCut);
  PrintStats(10, 1.5, ptCut);
}

//____________________________________________________________________
void AliCorrection::ResetErrorsOnCorrections()
{
  // resets the errors in the correction matrix
  
  fEventCorr->ResetErrorsOnCorrections();
  fTrackCorr->ResetErrorsOnCorrections();  
}
 AliCorrection.cxx:1
 AliCorrection.cxx:2
 AliCorrection.cxx:3
 AliCorrection.cxx:4
 AliCorrection.cxx:5
 AliCorrection.cxx:6
 AliCorrection.cxx:7
 AliCorrection.cxx:8
 AliCorrection.cxx:9
 AliCorrection.cxx:10
 AliCorrection.cxx:11
 AliCorrection.cxx:12
 AliCorrection.cxx:13
 AliCorrection.cxx:14
 AliCorrection.cxx:15
 AliCorrection.cxx:16
 AliCorrection.cxx:17
 AliCorrection.cxx:18
 AliCorrection.cxx:19
 AliCorrection.cxx:20
 AliCorrection.cxx:21
 AliCorrection.cxx:22
 AliCorrection.cxx:23
 AliCorrection.cxx:24
 AliCorrection.cxx:25
 AliCorrection.cxx:26
 AliCorrection.cxx:27
 AliCorrection.cxx:28
 AliCorrection.cxx:29
 AliCorrection.cxx:30
 AliCorrection.cxx:31
 AliCorrection.cxx:32
 AliCorrection.cxx:33
 AliCorrection.cxx:34
 AliCorrection.cxx:35
 AliCorrection.cxx:36
 AliCorrection.cxx:37
 AliCorrection.cxx:38
 AliCorrection.cxx:39
 AliCorrection.cxx:40
 AliCorrection.cxx:41
 AliCorrection.cxx:42
 AliCorrection.cxx:43
 AliCorrection.cxx:44
 AliCorrection.cxx:45
 AliCorrection.cxx:46
 AliCorrection.cxx:47
 AliCorrection.cxx:48
 AliCorrection.cxx:49
 AliCorrection.cxx:50
 AliCorrection.cxx:51
 AliCorrection.cxx:52
 AliCorrection.cxx:53
 AliCorrection.cxx:54
 AliCorrection.cxx:55
 AliCorrection.cxx:56
 AliCorrection.cxx:57
 AliCorrection.cxx:58
 AliCorrection.cxx:59
 AliCorrection.cxx:60
 AliCorrection.cxx:61
 AliCorrection.cxx:62
 AliCorrection.cxx:63
 AliCorrection.cxx:64
 AliCorrection.cxx:65
 AliCorrection.cxx:66
 AliCorrection.cxx:67
 AliCorrection.cxx:68
 AliCorrection.cxx:69
 AliCorrection.cxx:70
 AliCorrection.cxx:71
 AliCorrection.cxx:72
 AliCorrection.cxx:73
 AliCorrection.cxx:74
 AliCorrection.cxx:75
 AliCorrection.cxx:76
 AliCorrection.cxx:77
 AliCorrection.cxx:78
 AliCorrection.cxx:79
 AliCorrection.cxx:80
 AliCorrection.cxx:81
 AliCorrection.cxx:82
 AliCorrection.cxx:83
 AliCorrection.cxx:84
 AliCorrection.cxx:85
 AliCorrection.cxx:86
 AliCorrection.cxx:87
 AliCorrection.cxx:88
 AliCorrection.cxx:89
 AliCorrection.cxx:90
 AliCorrection.cxx:91
 AliCorrection.cxx:92
 AliCorrection.cxx:93
 AliCorrection.cxx:94
 AliCorrection.cxx:95
 AliCorrection.cxx:96
 AliCorrection.cxx:97
 AliCorrection.cxx:98
 AliCorrection.cxx:99
 AliCorrection.cxx:100
 AliCorrection.cxx:101
 AliCorrection.cxx:102
 AliCorrection.cxx:103
 AliCorrection.cxx:104
 AliCorrection.cxx:105
 AliCorrection.cxx:106
 AliCorrection.cxx:107
 AliCorrection.cxx:108
 AliCorrection.cxx:109
 AliCorrection.cxx:110
 AliCorrection.cxx:111
 AliCorrection.cxx:112
 AliCorrection.cxx:113
 AliCorrection.cxx:114
 AliCorrection.cxx:115
 AliCorrection.cxx:116
 AliCorrection.cxx:117
 AliCorrection.cxx:118
 AliCorrection.cxx:119
 AliCorrection.cxx:120
 AliCorrection.cxx:121
 AliCorrection.cxx:122
 AliCorrection.cxx:123
 AliCorrection.cxx:124
 AliCorrection.cxx:125
 AliCorrection.cxx:126
 AliCorrection.cxx:127
 AliCorrection.cxx:128
 AliCorrection.cxx:129
 AliCorrection.cxx:130
 AliCorrection.cxx:131
 AliCorrection.cxx:132
 AliCorrection.cxx:133
 AliCorrection.cxx:134
 AliCorrection.cxx:135
 AliCorrection.cxx:136
 AliCorrection.cxx:137
 AliCorrection.cxx:138
 AliCorrection.cxx:139
 AliCorrection.cxx:140
 AliCorrection.cxx:141
 AliCorrection.cxx:142
 AliCorrection.cxx:143
 AliCorrection.cxx:144
 AliCorrection.cxx:145
 AliCorrection.cxx:146
 AliCorrection.cxx:147
 AliCorrection.cxx:148
 AliCorrection.cxx:149
 AliCorrection.cxx:150
 AliCorrection.cxx:151
 AliCorrection.cxx:152
 AliCorrection.cxx:153
 AliCorrection.cxx:154
 AliCorrection.cxx:155
 AliCorrection.cxx:156
 AliCorrection.cxx:157
 AliCorrection.cxx:158
 AliCorrection.cxx:159
 AliCorrection.cxx:160
 AliCorrection.cxx:161
 AliCorrection.cxx:162
 AliCorrection.cxx:163
 AliCorrection.cxx:164
 AliCorrection.cxx:165
 AliCorrection.cxx:166
 AliCorrection.cxx:167
 AliCorrection.cxx:168
 AliCorrection.cxx:169
 AliCorrection.cxx:170
 AliCorrection.cxx:171
 AliCorrection.cxx:172
 AliCorrection.cxx:173
 AliCorrection.cxx:174
 AliCorrection.cxx:175
 AliCorrection.cxx:176
 AliCorrection.cxx:177
 AliCorrection.cxx:178
 AliCorrection.cxx:179
 AliCorrection.cxx:180
 AliCorrection.cxx:181
 AliCorrection.cxx:182
 AliCorrection.cxx:183
 AliCorrection.cxx:184
 AliCorrection.cxx:185
 AliCorrection.cxx:186
 AliCorrection.cxx:187
 AliCorrection.cxx:188
 AliCorrection.cxx:189
 AliCorrection.cxx:190
 AliCorrection.cxx:191
 AliCorrection.cxx:192
 AliCorrection.cxx:193
 AliCorrection.cxx:194
 AliCorrection.cxx:195
 AliCorrection.cxx:196
 AliCorrection.cxx:197
 AliCorrection.cxx:198
 AliCorrection.cxx:199
 AliCorrection.cxx:200
 AliCorrection.cxx:201
 AliCorrection.cxx:202
 AliCorrection.cxx:203
 AliCorrection.cxx:204
 AliCorrection.cxx:205
 AliCorrection.cxx:206
 AliCorrection.cxx:207
 AliCorrection.cxx:208
 AliCorrection.cxx:209
 AliCorrection.cxx:210
 AliCorrection.cxx:211
 AliCorrection.cxx:212
 AliCorrection.cxx:213
 AliCorrection.cxx:214
 AliCorrection.cxx:215
 AliCorrection.cxx:216
 AliCorrection.cxx:217
 AliCorrection.cxx:218
 AliCorrection.cxx:219
 AliCorrection.cxx:220
 AliCorrection.cxx:221
 AliCorrection.cxx:222
 AliCorrection.cxx:223
 AliCorrection.cxx:224
 AliCorrection.cxx:225
 AliCorrection.cxx:226
 AliCorrection.cxx:227
 AliCorrection.cxx:228
 AliCorrection.cxx:229
 AliCorrection.cxx:230
 AliCorrection.cxx:231
 AliCorrection.cxx:232
 AliCorrection.cxx:233
 AliCorrection.cxx:234
 AliCorrection.cxx:235
 AliCorrection.cxx:236
 AliCorrection.cxx:237
 AliCorrection.cxx:238
 AliCorrection.cxx:239
 AliCorrection.cxx:240
 AliCorrection.cxx:241
 AliCorrection.cxx:242
 AliCorrection.cxx:243
 AliCorrection.cxx:244
 AliCorrection.cxx:245
 AliCorrection.cxx:246
 AliCorrection.cxx:247
 AliCorrection.cxx:248
 AliCorrection.cxx:249
 AliCorrection.cxx:250
 AliCorrection.cxx:251
 AliCorrection.cxx:252
 AliCorrection.cxx:253
 AliCorrection.cxx:254
 AliCorrection.cxx:255
 AliCorrection.cxx:256
 AliCorrection.cxx:257
 AliCorrection.cxx:258
 AliCorrection.cxx:259
 AliCorrection.cxx:260
 AliCorrection.cxx:261
 AliCorrection.cxx:262
 AliCorrection.cxx:263
 AliCorrection.cxx:264
 AliCorrection.cxx:265
 AliCorrection.cxx:266
 AliCorrection.cxx:267
 AliCorrection.cxx:268
 AliCorrection.cxx:269
 AliCorrection.cxx:270
 AliCorrection.cxx:271
 AliCorrection.cxx:272
 AliCorrection.cxx:273
 AliCorrection.cxx:274
 AliCorrection.cxx:275
 AliCorrection.cxx:276
 AliCorrection.cxx:277
 AliCorrection.cxx:278
 AliCorrection.cxx:279
 AliCorrection.cxx:280
 AliCorrection.cxx:281
 AliCorrection.cxx:282
 AliCorrection.cxx:283
 AliCorrection.cxx:284
 AliCorrection.cxx:285
 AliCorrection.cxx:286
 AliCorrection.cxx:287
 AliCorrection.cxx:288
 AliCorrection.cxx:289
 AliCorrection.cxx:290
 AliCorrection.cxx:291
 AliCorrection.cxx:292
 AliCorrection.cxx:293
 AliCorrection.cxx:294
 AliCorrection.cxx:295
 AliCorrection.cxx:296
 AliCorrection.cxx:297
 AliCorrection.cxx:298
 AliCorrection.cxx:299
 AliCorrection.cxx:300
 AliCorrection.cxx:301
 AliCorrection.cxx:302
 AliCorrection.cxx:303
 AliCorrection.cxx:304
 AliCorrection.cxx:305
 AliCorrection.cxx:306
 AliCorrection.cxx:307
 AliCorrection.cxx:308
 AliCorrection.cxx:309
 AliCorrection.cxx:310
 AliCorrection.cxx:311
 AliCorrection.cxx:312
 AliCorrection.cxx:313
 AliCorrection.cxx:314
 AliCorrection.cxx:315
 AliCorrection.cxx:316
 AliCorrection.cxx:317
 AliCorrection.cxx:318
 AliCorrection.cxx:319
 AliCorrection.cxx:320
 AliCorrection.cxx:321
 AliCorrection.cxx:322
 AliCorrection.cxx:323
 AliCorrection.cxx:324
 AliCorrection.cxx:325
 AliCorrection.cxx:326
 AliCorrection.cxx:327
 AliCorrection.cxx:328
 AliCorrection.cxx:329
 AliCorrection.cxx:330
 AliCorrection.cxx:331
 AliCorrection.cxx:332
 AliCorrection.cxx:333
 AliCorrection.cxx:334
 AliCorrection.cxx:335
 AliCorrection.cxx:336
 AliCorrection.cxx:337
 AliCorrection.cxx:338
 AliCorrection.cxx:339
 AliCorrection.cxx:340
 AliCorrection.cxx:341
 AliCorrection.cxx:342
 AliCorrection.cxx:343
 AliCorrection.cxx:344
 AliCorrection.cxx:345
 AliCorrection.cxx:346
 AliCorrection.cxx:347
 AliCorrection.cxx:348
 AliCorrection.cxx:349
 AliCorrection.cxx:350
 AliCorrection.cxx:351
 AliCorrection.cxx:352
 AliCorrection.cxx:353
 AliCorrection.cxx:354
 AliCorrection.cxx:355
 AliCorrection.cxx:356
 AliCorrection.cxx:357
 AliCorrection.cxx:358
 AliCorrection.cxx:359
 AliCorrection.cxx:360
 AliCorrection.cxx:361
 AliCorrection.cxx:362
 AliCorrection.cxx:363
 AliCorrection.cxx:364
 AliCorrection.cxx:365
 AliCorrection.cxx:366
 AliCorrection.cxx:367
 AliCorrection.cxx:368
 AliCorrection.cxx:369
 AliCorrection.cxx:370
 AliCorrection.cxx:371
 AliCorrection.cxx:372
 AliCorrection.cxx:373
 AliCorrection.cxx:374
 AliCorrection.cxx:375
 AliCorrection.cxx:376
 AliCorrection.cxx:377
 AliCorrection.cxx:378
 AliCorrection.cxx:379
 AliCorrection.cxx:380
 AliCorrection.cxx:381
 AliCorrection.cxx:382
 AliCorrection.cxx:383
 AliCorrection.cxx:384
 AliCorrection.cxx:385
 AliCorrection.cxx:386
 AliCorrection.cxx:387
 AliCorrection.cxx:388
 AliCorrection.cxx:389
 AliCorrection.cxx:390
 AliCorrection.cxx:391
 AliCorrection.cxx:392
 AliCorrection.cxx:393
 AliCorrection.cxx:394
 AliCorrection.cxx:395
 AliCorrection.cxx:396
 AliCorrection.cxx:397
 AliCorrection.cxx:398
 AliCorrection.cxx:399
 AliCorrection.cxx:400
 AliCorrection.cxx:401
 AliCorrection.cxx:402
 AliCorrection.cxx:403
 AliCorrection.cxx:404
 AliCorrection.cxx:405
 AliCorrection.cxx:406
 AliCorrection.cxx:407
 AliCorrection.cxx:408
 AliCorrection.cxx:409
 AliCorrection.cxx:410
 AliCorrection.cxx:411
 AliCorrection.cxx:412
 AliCorrection.cxx:413
 AliCorrection.cxx:414
 AliCorrection.cxx:415
 AliCorrection.cxx:416
 AliCorrection.cxx:417
 AliCorrection.cxx:418
 AliCorrection.cxx:419
 AliCorrection.cxx:420
 AliCorrection.cxx:421
 AliCorrection.cxx:422
 AliCorrection.cxx:423
 AliCorrection.cxx:424
 AliCorrection.cxx:425
 AliCorrection.cxx:426
 AliCorrection.cxx:427
 AliCorrection.cxx:428
 AliCorrection.cxx:429
 AliCorrection.cxx:430
 AliCorrection.cxx:431
 AliCorrection.cxx:432
 AliCorrection.cxx:433
 AliCorrection.cxx:434
 AliCorrection.cxx:435
 AliCorrection.cxx:436
 AliCorrection.cxx:437
 AliCorrection.cxx:438
 AliCorrection.cxx:439
 AliCorrection.cxx:440
 AliCorrection.cxx:441
 AliCorrection.cxx:442
 AliCorrection.cxx:443
 AliCorrection.cxx:444
 AliCorrection.cxx:445
 AliCorrection.cxx:446
 AliCorrection.cxx:447
 AliCorrection.cxx:448
 AliCorrection.cxx:449
 AliCorrection.cxx:450
 AliCorrection.cxx:451
 AliCorrection.cxx:452
 AliCorrection.cxx:453
 AliCorrection.cxx:454
 AliCorrection.cxx:455
 AliCorrection.cxx:456
 AliCorrection.cxx:457
 AliCorrection.cxx:458
 AliCorrection.cxx:459
 AliCorrection.cxx:460
 AliCorrection.cxx:461
 AliCorrection.cxx:462
 AliCorrection.cxx:463
 AliCorrection.cxx:464
 AliCorrection.cxx:465
 AliCorrection.cxx:466
 AliCorrection.cxx:467
 AliCorrection.cxx:468
 AliCorrection.cxx:469
 AliCorrection.cxx:470
 AliCorrection.cxx:471
 AliCorrection.cxx:472
 AliCorrection.cxx:473
 AliCorrection.cxx:474
 AliCorrection.cxx:475
 AliCorrection.cxx:476
 AliCorrection.cxx:477
 AliCorrection.cxx:478
 AliCorrection.cxx:479
 AliCorrection.cxx:480
 AliCorrection.cxx:481
 AliCorrection.cxx:482
 AliCorrection.cxx:483
 AliCorrection.cxx:484
 AliCorrection.cxx:485
 AliCorrection.cxx:486
 AliCorrection.cxx:487
 AliCorrection.cxx:488
 AliCorrection.cxx:489
 AliCorrection.cxx:490
 AliCorrection.cxx:491
 AliCorrection.cxx:492
 AliCorrection.cxx:493
 AliCorrection.cxx:494
 AliCorrection.cxx:495
 AliCorrection.cxx:496
 AliCorrection.cxx:497
 AliCorrection.cxx:498
 AliCorrection.cxx:499
 AliCorrection.cxx:500
 AliCorrection.cxx:501
 AliCorrection.cxx:502
 AliCorrection.cxx:503