ROOT logo
// $Id$
// 
// Analysis task to estimate an event's local energy density
//
// This task is part of the emcal jet framework and should be run in the emcaljet train
// The following extensions to an accepted AliVEvent are expected:
//      - (anti-kt) jets -> necessary if one wants to exclude leading jet contribution to the event plane
//      - background estimate of rho -> this task estimates modulation, not rho itself
//      - pico tracks -> a uniform track selection is necessary to estimate the contribution of v_n harmonics
//      aod's and esd's are handled transparently
// The task will estimates a phi-dependent background density rho 
// which is added to the event as a AliLocalRhoParamter object
//
// Author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
//         (rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl)

// root includes
#include <TStyle.h>
#include <TRandom3.h>
#include <TChain.h>
#include <TMath.h>
#include <TF1.h>
#include <TF2.h>
#include <TH1F.h>
#include <TH2F.h>
#include <TProfile.h>
// aliroot includes
#include <AliAnalysisTask.h>
#include <AliAnalysisManager.h>
#include <AliCentrality.h>
#include <AliVVertex.h>
#include <AliESDEvent.h>
#include <AliAODEvent.h>
#include <AliAODTrack.h>
// emcal jet framework includes
#include <AliPicoTrack.h>
#include <AliEmcalJet.h>
#include <AliRhoParameter.h>
#include <AliLocalRhoParameter.h>
#include <AliAnalysisTaskLocalRho.h>

class AliAnalysisTaskLocalRho;
using namespace std;

ClassImp(AliAnalysisTaskLocalRho)

//_____________________________________________________________________________
AliAnalysisTaskLocalRho::AliAnalysisTaskLocalRho() : 
  AliAnalysisTaskEmcalJet("AliAnalysisTaskLocalRho", kTRUE), 
  fDebug(0), fInitialized(0), fAttachToEvent(kTRUE), fFillHistograms(kFALSE), fNoEventWeightsForQC(kTRUE), 
  fUseScaledRho(0), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), 
  fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fInCentralitySelection(-1), 
  fFitModulationType(kNoFit), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kFALSE), fDetectorType(kTPC), 
  fFitModulationOptions("WLQI"), fRunModeType(kGrid), fFitModulation(0), fMinPvalue(0.01), fMaxPvalue(1), 
  fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), 
  fSoftTrackMaxPt(5.), fHistPvalueCDF(0), fHistRhoStatusCent(0), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), 
  fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fOutputList(0), 
  fOutputListGood(0), fOutputListBad(0), fHistSwap(0), fHistAnalysisSummary(0), fProfV2(0), fProfV2Cumulant(0), 
  fProfV3(0), fProfV3Cumulant(0) 
{
  // Default constructor

  for(Int_t i(0); i < 10; i++) {
    fHistPsi2[i] = 0; 
    fHistPsi3[i] = 0;
  }
}

//_____________________________________________________________________________
AliAnalysisTaskLocalRho::AliAnalysisTaskLocalRho(const char* name, runModeType type) : 
  AliAnalysisTaskEmcalJet(name, kTRUE),
  fDebug(0), fInitialized(0), fAttachToEvent(kTRUE), fFillHistograms(kFALSE), fNoEventWeightsForQC(kTRUE), 
  fUseScaledRho(0), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), 
  fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fInCentralitySelection(-1), 
  fFitModulationType(kNoFit), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kFALSE), fDetectorType(kTPC), 
  fFitModulationOptions("WLQI"), fRunModeType(type), fFitModulation(0), fMinPvalue(0.01), fMaxPvalue(1), 
  fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), 
  fSoftTrackMaxPt(5.), fHistPvalueCDF(0), fHistRhoStatusCent(0), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), 
  fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fOutputList(0), 
  fOutputListGood(0), fOutputListBad(0), fHistSwap(0), fHistAnalysisSummary(0), fProfV2(0), fProfV2Cumulant(0), 
  fProfV3(0), fProfV3Cumulant(0) 
{
  // Constructor
  for(Int_t i(0); i < 10; i++) {
    fHistPsi2[i] = 0; 
    fHistPsi3[i] = 0;
  }

  DefineInput(0, TChain::Class());
  DefineOutput(1, TList::Class());
  switch (fRunModeType) {
  case kLocal : {
    gStyle->SetOptFit(1);
    DefineOutput(2, TList::Class());
    DefineOutput(3, TList::Class());
  } break;
  default: fDebug = -1;   // suppress debug info explicitely when not running locally
  }
}

//_____________________________________________________________________________
AliAnalysisTaskLocalRho::~AliAnalysisTaskLocalRho()
{
  // destructor
  if(fOutputList)             delete fOutputList;
  if(fOutputListGood)         delete fOutputListGood;
  if(fOutputListBad)          delete fOutputListBad;
  if(fFitModulation)          delete fFitModulation;
  if(fHistSwap)               delete fHistSwap;
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::ExecOnce()
{
  // Init the analysis
  if(fLocalRhoName=="") fLocalRhoName = Form("LocalRhoFrom_%s", GetName());
  fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0); 
  // add the local rho to the event if necessary
  if(fAttachToEvent) {
    if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
      InputEvent()->AddObject(fLocalRho);
    } else {
      AliFatal(Form("%s: Container with same name %s already present. Aborting", GetName(), fLocalRho->GetName()));
    }
  }
  AliAnalysisTaskEmcalJet::ExecOnce();        // init the base clas
  if(fUseScaledRho) {
    // unscaled rho has been retrieved by the parent class, now we retrieve rho scaled
    fRho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(Form("%s_Scaled", fRho->GetName())));
    if(!fRho) {
      AliFatal(Form("%s: Couldn't find container for scaled rho. Aborting !", GetName()));
    }
  }
  if(!GetJetContainer()) AliFatal(Form("%s: Couldn't get jet container. Aborting !", GetName()));
}

//_____________________________________________________________________________
Bool_t AliAnalysisTaskLocalRho::InitializeAnalysis() 
{
  // Initialize the anaysis

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(fLocalJetMinEta > -10 && fLocalJetMaxEta > -10) SetJetEtaLimits(fLocalJetMinEta, fLocalJetMaxEta);
  if(fLocalJetMinPhi > -10 && fLocalJetMaxPhi > -10) SetJetPhiLimits(fLocalJetMinPhi, fLocalJetMaxPhi);
  switch (fFitModulationType)  {
  case kNoFit : { SetModulationFit(new TF1("fit_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
  case kV2 : {
    SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
    fFitModulation->SetParameter(0, 0.);        // normalization
    fFitModulation->SetParameter(3, 0.2);       // v2
    fFitModulation->FixParameter(1, 1.);        // constant
    fFitModulation->FixParameter(2, 2.);        // constant
  } break;
  case kV3: {
    SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
    fFitModulation->SetParameter(0, 0.);        // normalization
    fFitModulation->SetParameter(3, 0.2);       // v3
    fFitModulation->FixParameter(1, 1.);        // constant
    fFitModulation->FixParameter(2, 3.);        // constant
  } break;
  default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
    SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
    fFitModulation->SetParameter(0, 0.);       // normalization
    fFitModulation->SetParameter(3, 0.2);      // v2
    fFitModulation->FixParameter(1, 1.);       // constant
    fFitModulation->FixParameter(2, 2.);       // constant
    fFitModulation->FixParameter(5, 3.);       // constant
    fFitModulation->SetParameter(7, 0.2);      // v3
  } break;
  }
  switch (fRunModeType) {
  case kGrid : { fFitModulationOptions += "N0"; } break;
  default : break;
  }
  FillAnalysisSummaryHistogram();
  return kTRUE;
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::UserCreateOutputObjects()
{
  // create output objects
  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
  if(!fCentralityClasses) {   // classes must be defined at this point
    Int_t c[] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100};
    fCentralityClasses = new TArrayI(sizeof(c)/sizeof(c[0]), c);
  }
  fOutputList = new TList();
  fOutputList->SetOwner(kTRUE);
  // the analysis summary histo which stores all the analysis flags is always written to file
  fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 51, -0.5, 51.5);
  if(!fFillHistograms) {
    PostData(1, fOutputList);
    return;
  }
  for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {     
    fHistPsi2[i] = BookTH1F("fHistPsi2", "#Psi_{2}", 100, -.5*TMath::Pi(), .5*TMath::Pi(), i);
    fHistPsi3[i] = BookTH1F("fHistPsi3", "#Psi_{3}", 100, -1.*TMath::Pi()/3., TMath::Pi()/3., i);
  }
  // cdf of chisquare distribution
  fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 500, 0, 1);
  fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [0=ok, 1=failed]", 101, -1, 100, 2, -.5, 1.5);
  // vn profiles
  Float_t temp[fCentralityClasses->GetSize()];
  for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
  fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
  fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
  fOutputList->Add(fProfV2);
  fOutputList->Add(fProfV3);
  switch (fFitModulationType) {
  case kQC2 : {
    fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
    fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
    fOutputList->Add(fProfV2Cumulant);
    fOutputList->Add(fProfV3Cumulant);
  } break;
  case kQC4 : {
    fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
    fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
    fOutputList->Add(fProfV2Cumulant);
    fOutputList->Add(fProfV3Cumulant);
  } break;
  default : break;
  }
  if(fUsePtWeight) fHistSwap->Sumw2();
  if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
  if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
  if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
  if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
  // increase readability of output list
  fOutputList->Sort();
  PostData(1, fOutputList);
  switch (fRunModeType) {
  case kLocal : {
    fOutputListGood = new TList();
    fOutputListGood->SetOwner(kTRUE);
    fOutputListBad = new TList();
    fOutputListBad->SetOwner(kTRUE);
    PostData(2, fOutputListGood);
    PostData(3, fOutputListBad);
  } break;
  default: break;
  }
}

//_____________________________________________________________________________
TH1F* AliAnalysisTaskLocalRho::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
{
  // Book a TH1F and connect it to the output container

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(!fOutputList) return 0x0;
  TString title(name);
  if(c!=-1) { // format centrality dependent histograms accordingly
    name = Form("%s_%i", name, c);
    title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
  }
  title += Form(";%s;[counts]", x);
  TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
  histogram->Sumw2();
  if(append) fOutputList->Add(histogram);
  return histogram;   
}

//_____________________________________________________________________________
TH2F* AliAnalysisTaskLocalRho::BookTH2F(const char* name, const char* x, const char*y, Int_t binsx, Double_t minx, Double_t maxx, 
					Int_t binsy, Double_t miny, Double_t maxy, Int_t c, Bool_t append)
{
  // Book a TH2F and connect it to the output container

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(!fOutputList) return 0x0;
  TString title(name);
  if(c!=-1) { // format centrality dependent histograms accordingly
    name = Form("%s_%i", name, c);
    title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
  }
  title += Form(";%s;%s", x, y);
  TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
  histogram->Sumw2();
  if(append) fOutputList->Add(histogram);
  return histogram;   
}

//_____________________________________________________________________________
Bool_t AliAnalysisTaskLocalRho::Run()
{
  // Execute once for each event

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(!(InputEvent()||fTracks||fJets||fRho)) return kFALSE;
  if(!fInitialized) fInitialized = InitializeAnalysis();
  // get the centrality bin (necessary for some control histograms
  fInCentralitySelection = -1;
  Double_t cent(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
  for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
    if(cent >= fCentralityClasses->At(i) && cent <= fCentralityClasses->At(1+i)) {
      fInCentralitySelection = i;
      break; }
  }
  if(fInCentralitySelection < 0) return kFALSE;
  // set the rho value 
  fLocalRho->SetVal(fRho->GetVal());
  // set the correct event plane accordign to the requested reference detector
  Double_t psi2(-1), psi3(-1);
  switch (fDetectorType) {    // determine the detector type for the rho fit
  case kTPC :     { 
    // [0] psi2         [1] psi3
    Double_t tpc[2];
    CalculateEventPlaneTPC(tpc);
    psi2 = tpc[0];         psi3 = tpc[1]; 
  } break;
  case kVZEROA :  { 
    // [0][0] psi2a     [1,0]   psi2c
    // [0][1] psi3a     [1,1]   psi3c
    Double_t vzero[2][2];
    CalculateEventPlaneVZERO(vzero);
    psi2 = vzero[0][0];    psi3 = vzero[0][1]; 
  }   break;  
  case kVZEROC :  { 
    // [0][0] psi2a     [1,0]   psi2c
    // [0][1] psi3a     [1,1]   psi3c
    Double_t vzero[2][2];
    CalculateEventPlaneVZERO(vzero);
    psi2 = vzero[1][0];    psi3 = vzero[1][1]; 
  }   break;
  case kVZEROComb : { 
    /* for the combined vzero event plane
     * [0] psi2         [1] psi3
     * not fully implmemented yet, use with caution ! */
    Double_t vzeroComb[2];
    CalculateEventPlaneCombinedVZERO(vzeroComb);
    psi2 = vzeroComb[0]; psi3 = vzeroComb[1];
  } break;
  default : break;
  }
  if(fFillHistograms) FillEventPlaneHistograms(psi2, psi3);
  switch (fFitModulationType) { // do the fits
  case kNoFit : { fFitModulation->FixParameter(0, fLocalRho->GetVal()); } break;
  case kV2 : {    // only v2
    if(CorrectRho(psi2, psi3)) {
      if(fFillHistograms) fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
      if(fUserSuppliedR2) {
	Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
	if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
      }
    }
  } break;
  case kV3 : {    // only v3
    if(CorrectRho(psi2, psi3)) {
      if(fUserSuppliedR3) {
	Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
	if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
      }
      if(fFillHistograms) fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
    }
  } break;
  case kQC2 : {   // qc2 analysis - NOTE: not a wise idea to use this !
    if(CorrectRho(psi2, psi3)) {
      if(fUserSuppliedR2 && fUserSuppliedR3) {
	// note for the qc method, resolution is REVERSED to go back to v2obs
	Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
	Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
	if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
	if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
      }
      if (fUsePtWeight) { // use weighted weights
	Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
	if(fFillHistograms) {
	  fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
	  fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11); 
	}
      } else {
	Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
	if(fFillHistograms) {
	  fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
	  fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
	}
      }
    }
  } break;
  case kQC4 : {   // NOTE: see comment at kQC2
    if(CorrectRho(psi2, psi3)) {
      if(fUserSuppliedR2 && fUserSuppliedR3) {
	// note for the qc method, resolution is REVERSED to go back to v2obs   
	Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
	Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
	if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
	if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
      }
      if (fUsePtWeight) { // use weighted weights
	if(fFillHistograms) {
	  fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
	  fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/); 
	}
      } else {
	if(fFillHistograms) {
	  fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
	  fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
	}
      }
    }
  } break;
  default : {
    if(CorrectRho(psi2, psi3)) {
      if(fUserSuppliedR2 && fUserSuppliedR3) {
	Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
	Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
	if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
	if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
      }
      if(fFillHistograms) {
	fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
	fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
      }
    }
  } break;
  }
  // if all went well, add local rho
  fLocalRho->SetLocalRho(fFitModulation);
  PostData(1, fOutputList);
  return kTRUE;
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const 
{
  // Get the vzero event plane
  if(fUseV0EventPlaneFromHeader) {
    // use the vzero event plane from the event header
    // note: to use the calibrated vzero event plane, run 
    // $ALICE_ROOT/ANALYSIS/macros/AddTaskVZEROEPSelection.C
    // prior to this task (make sure the calibration is available for the dataset
    // you want to use)
    Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
    vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
    vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
    vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
    vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
    return;
  }
  // grab the vzero event plane without recentering
  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0);    // for psi2
  Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0);    // for psi3
  for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
    Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
    //        (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
    if(iVZERO<32) {
      qxa2 += weight*TMath::Cos(2.*phi);
      qya2 += weight*TMath::Sin(2.*phi);
      qxa3 += weight*TMath::Cos(3.*phi);
      qya3 += weight*TMath::Sin(3.*phi);
    }
    else {
      qxc2 += weight*TMath::Cos(2.*phi);
      qyc2 += weight*TMath::Sin(2.*phi);
      qxc3 += weight*TMath::Cos(3.*phi);
      qyc3 += weight*TMath::Sin(3.*phi);
    }
  }
  vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
  vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
  vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
  vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::CalculateEventPlaneTPC(Double_t* tpc)
{
  // Grab the TPC event plane. if parameter fExcludeLeadingJetsFromFit is larger than 0, 
  // strip in eta of width fExcludeLeadingJetsFromFit * GetJetContainer()->GetJetRadius() around the leading jet (before
  // subtraction of rho) will be exluded from the event plane estimate

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  fNAcceptedTracks = 0;                // reset the track counter
  Double_t qx2(0), qy2(0);     // for psi2
  Double_t qx3(0), qy3(0);     // for psi3
  if(fTracks) {
    Float_t excludeInEta = -999;
    if(fExcludeLeadingJetsFromFit > 0 ) {    // remove the leading jet from ep estimate
      AliEmcalJet* leadingJet(GetJetContainer()->GetLeadingJet());
      if(leadingJet) leadingJet->Eta();
    }
    Int_t iTracks(fTracks->GetEntriesFast());
    for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
      AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
      if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
      if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
      fNAcceptedTracks++;
      qx2+= TMath::Cos(2.*track->Phi());
      qy2+= TMath::Sin(2.*track->Phi());
      qx3+= TMath::Cos(3.*track->Phi());
      qy3+= TMath::Sin(3.*track->Phi());
    }
  }
  tpc[0] = .5*TMath::ATan2(qy2, qx2);
  tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
} 

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
{
  // Grab the combined vzero event plane

  //    if(fUseV0EventPlaneFromHeader) {    // use the vzero from the header
  Double_t a(0), b(0), c(0), d(0);
  comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
  comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
  // FIXME the rest of this function isn't impelmented yet (as of 01-07-2013)
  // this means a default the combined vzero event plane from the header is used
  // to get this value 'by hand', vzeroa and vzeroc event planes have to be combined
  // according to their resolution - this will be added ...
  //
  //    } else {
  //        Double_t qx2a(0), qy2a(0), qx2c(0), qy2c(0), qx3a(0), qy3a(0), qx3c(0), qy3c(0);
  //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, qx2a, qy2a);
  //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, qx2c, qy2c);
  //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, qx3a, qy3a);
  //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, qx3c, qy3c);
  //        Double_t chi2A(-1), chi2C(-1), chi3A(-1), chi3C(-1);     // get chi from the resolution
  //        Double_t qx2(chi2A*chi2A*qx2a+chi2C*chi2C*qx2c);
  //        Double_t qy2(chi2A*chi2A*qy2a+chi2C*chi2C*qy2c);
  //        Double_t qx3(chi3A*chi3A*qx3a+chi3C*chi3C*qx3c);
  //        Double_t qy3(chi3A*chi3A*qy3a+chi3C*chi3C*qy3c);
  //        comb[0] = .5*TMath::ATan2(qy2, qx2);
  //        comb[1] = (1./3.)*TMath::ATan2(qy3, qx3);
  //    }
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::CalculateQC2(Int_t harm) 
{
  // Get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
  if(fUsePtWeight) {  // for the weighted 2-nd order q-cumulant
    QCnQnk(harm, 1, reQ, imQ);      // get the weighted 2-nd order q-vectors
    modQ = reQ*reQ+imQ*imQ;         // get abs Q-squared
    M11 = QCnM11();                 // equals S2,1 - S1,2
    return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
  } // else return the non-weighted 2-nd order q-cumulant
  QCnQnk(harm, 0, reQ, imQ);          // get the non-weighted 2-nd order q-vectors
  modQ = reQ*reQ+imQ*imQ;             // get abs Q-squared
  M = QCnM();
  return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::CalculateQC4(Int_t harm) 
{
  // Get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
  Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0);  // terms of the calculation
  if(fUsePtWeight) {  // for the weighted 4-th order q-cumulant
    QCnQnk(harm, 1, reQn1, imQn1);
    QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
    QCnQnk(harm, 3, reQn3, imQn3);
    // fill in the terms ...
    a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
    b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
    c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
    d = 8.*(reQn3*reQn1+imQn3*imQn1);
    e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
    f = -6.*QCnS(1,4);
    g = 2.*QCnS(2,2);
    M1111 = QCnM1111();
    return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
  }   // else return the unweighted case
  Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
  QCnQnk(harm, 0, reQn, imQn);
  QCnQnk(harm*2, 0, reQ2n, imQ2n);
  // fill in the terms ...
  M = QCnM();
  if(M < 4) return -999;
  a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
  b = reQ2n*reQ2n + imQ2n*imQ2n;
  c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
  e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
  f = 2.*M*(M-3);
  return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) 
{
  // Get the weighted n-th order q-vector, pass real and imaginary part as reference

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(!fTracks) return;
  fNAcceptedTracksQCn = 0;
  Int_t iTracks(fTracks->GetEntriesFast());
  for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
    AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
    if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
    fNAcceptedTracksQCn++;
    // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
    reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
    imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
  }
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::QCnS(Int_t i, Int_t j) 
{
  // Get the weighted ij-th order autocorrelation correction

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(!fTracks || i <= 0 || j <= 0) return -999;
  Int_t iTracks(fTracks->GetEntriesFast());
  Double_t Sij(0);
  for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
    AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
    if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
    Sij+=TMath::Power(track->Pt(), j);
  }
  return TMath::Power(Sij, i);
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::QCnM() 
{
  // Get multiplicity for unweighted q-cumulants. function QCnQnk should be called first

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  return (Double_t) fNAcceptedTracksQCn;
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::QCnM11() 
{
  // Get multiplicity weights for the weighted two particle cumulant

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  return (QCnS(2,1) - QCnS(1,2));
}

//_____________________________________________________________________________
Double_t AliAnalysisTaskLocalRho::QCnM1111() 
{
  // Get multiplicity weights for the weighted four particle cumulant

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  return (QCnS(4,1)-6*QCnS(1,2)*QCnS(2,1)+8*QCnS(1,3)*QCnS(1,1)+3*QCnS(2,2)-6*QCnS(1,4));
}

//_____________________________________________________________________________
Bool_t AliAnalysisTaskLocalRho::QCnRecovery(Double_t psi2, Double_t psi3) 
{
  // Decides how to deal with the situation where c2 or c3 is negative 
  // Returns kTRUE depending on whether or not a modulated rho is used for the jet background

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
    fFitModulation->SetParameter(7, 0);
    fFitModulation->SetParameter(3, 0);
    fFitModulation->SetParameter(0, fLocalRho->GetVal());
    return kTRUE;   // v2 and v3 have physical null values
  }
  switch (fQCRecovery) {
  case kFixedRho : {      // roll back to the original rho
    fFitModulation->SetParameter(7, 0);
    fFitModulation->SetParameter(3, 0);
    fFitModulation->SetParameter(0, fLocalRho->GetVal());
    return kFALSE;       // rho is forced to be fixed
  }
  case kNegativeVn : {
    Double_t c2(fFitModulation->GetParameter(3));
    Double_t c3(fFitModulation->GetParameter(7));
    if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
    if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
    fFitModulation->SetParameter(3, c2);
    fFitModulation->SetParameter(7, c3);
    return kTRUE;        // is this a physical quantity ?
  }
  case kTryFit : {
    fitModulationType tempType(fFitModulationType);  // store temporarily
    fFitModulationType = kCombined;
    fFitModulation->SetParameter(7, 0);
    fFitModulation->SetParameter(3, 0);
    Bool_t pass(CorrectRho(psi2, psi3));         // do the fit and all quality checks
    fFitModulationType = tempType;               // roll back for next event
    return pass;
  }
  default : return kFALSE;
  }
  return kFALSE;
}

//_____________________________________________________________________________
Bool_t AliAnalysisTaskLocalRho::CorrectRho(Double_t psi2, Double_t psi3) 
{
  // Get rho' -> rho(phi)
  // three routines are available, 1 and 2 can be used with or without pt weights
  //  [1] get vn from q-cumulants
  //      in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
  //      are expected. a check is performed to see if rho has no negative local minimum
  //      for full description, see Phys. Rev. C 83, 044913
  //      since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
  //      in this case one can either roll back to the 'original' fixed rho, do a fit for vn or take use
  //      vn = - sqrt(|cn|) note that because of this, use of q-cumulants is not safe ! 
  //  [2] fitting a fourier expansion to the de/dphi distribution
  //      the fit can be done with either v2, v3 or a combination.
  //      in all cases, a cut can be made on the p-value of the chi-squared value of the fit
  //      and a check can be performed to see if rho has no negative local minimum
  //  [3] get v2 and v3 from user supplied histograms
  //      in this way, a fixed value of v2 and v3 is subtracted w.r.t. whichever event plane is requested

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  switch (fFitModulationType) {       // for approaches where no fitting is required
  case kQC2 : {
    fFitModulation->FixParameter(4, psi2); 
    fFitModulation->FixParameter(6, psi3);
    fFitModulation->FixParameter(3, CalculateQC2(2));   // set here with cn, vn = sqrt(cn)
    fFitModulation->FixParameter(7, CalculateQC2(3));
    // first fill the histos of the raw cumulant distribution
    if (fUsePtWeight) { // use weighted weights
      Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
      if(fFillHistograms) {
	fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
	fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
      }
    } else {
      Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
      if(fFillHistograms) {
	fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
	fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
      }
    }
    // then see if one of the cn value is larger than zero and vn is readily available
    if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
      fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
      fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
    } else if (!QCnRecovery(psi2, psi3)) return kFALSE;  // try to recover the cumulant, this will set v2 and v3
    if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {  // general check 
      fFitModulation->SetParameter(7, 0);
      fFitModulation->SetParameter(3, 0);
      fFitModulation->SetParameter(0, fLocalRho->GetVal());
      return kFALSE;
    }
    return kTRUE;
  } break;
  case kQC4 : {
    fFitModulation->FixParameter(4, psi2); 
    fFitModulation->FixParameter(6, psi3);
    fFitModulation->FixParameter(3, CalculateQC4(2));   // set here with cn, vn = sqrt(cn)
    fFitModulation->FixParameter(7, CalculateQC4(3));
    // first fill the histos of the raw cumulant distribution
    if (fUsePtWeight) { // use weighted weights
      if(fFillHistograms) {
	fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
	fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
      }
    } else {
      if(fFillHistograms) {
	fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
	fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
      }
    }
    // then see if one of the cn value is larger than zero and vn is readily available
    if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
      fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
      fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
    } else if (!QCnRecovery(psi2, psi3)) return kFALSE;  // try to recover the cumulant, this will set v2 and v3
    if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {  // general check 
      fFitModulation->SetParameter(7, 0);
      fFitModulation->SetParameter(3, 0);
      fFitModulation->SetParameter(0, fLocalRho->GetVal());
      return kFALSE;
    }
  } break;
  case kIntegratedFlow : {
    // use v2 and v3 values from an earlier iteration over the data
    fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
    fFitModulation->FixParameter(4, psi2);
    fFitModulation->FixParameter(6, psi3);
    fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
    if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { 
      fFitModulation->SetParameter(7, 0);
      fFitModulation->SetParameter(3, 0);
      fFitModulation->SetParameter(0, fLocalRho->GetVal());
      return kFALSE;
    }
    return kTRUE;
  }
  default : break;
  }
  TString detector("");
  switch (fDetectorType) {
  case kTPC : detector+="TPC";
    break;
  case kVZEROA : detector+="VZEROA";
    break;
  case kVZEROC : detector+="VZEROC";
    break;
  case kVZEROComb : detector+="VZEROComb";
    break; 
  default: break;
  }
  Int_t iTracks(fTracks->GetEntriesFast());
  Double_t excludeInEta = -999;
  Double_t excludeInPhi = -999;
  Double_t excludeInPt  = -999;
  if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE;   // no use fitting an empty event ...
  if(fExcludeLeadingJetsFromFit > 0 ) {
    AliEmcalJet* leadingJet = GetJetContainer()->GetLeadingJet();
    if(PassesCuts(leadingJet)) {
      excludeInEta = leadingJet->Eta();
      excludeInPhi = leadingJet->Phi();
      excludeInPt  = leadingJet->Pt();
    }
  }
  fHistSwap->Reset();                 // clear the histogram
  TH1F _tempSwap;     // on stack for quick access
  TH1F _tempSwapN;    // on stack for quick access, bookkeeping histogram
  if(fRebinSwapHistoOnTheFly) {
    if(fNAcceptedTracks < 49) fNAcceptedTracks = 49;       // avoid aliasing effects
    _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), 0, TMath::TwoPi());
    if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), 0, TMath::TwoPi());
    if(fUsePtWeight) _tempSwap.Sumw2();
  }
  else _tempSwap = *fHistSwap;         // now _tempSwap holds the desired histo
  // non poissonian error when using pt weights
  Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
  for(Int_t i(0); i < iTracks; i++) {
    AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
    if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
    if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
    if(fUsePtWeight) {
      _tempSwap.Fill(track->Phi(), track->Pt());
      if(fUsePtWeightErrorPropagation) {
        totalpts += track->Pt();
        totalptsquares += track->Pt()*track->Pt();
        totalns += 1;
        _tempSwapN.Fill(track->Phi());
      }
    }
  else _tempSwap.Fill(track->Phi());
  }
  if(fUsePtWeight && fUsePtWeightErrorPropagation) {
    // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
    // the assumption here is that the bin error will be dominated by the uncertainty in the mean pt in a bin and in the uncertainty
    // of the number of tracks in a bin, the first of which will be estimated from the sample standard deviation of all tracks in the 
    // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
    if(totalns < 1) return kFALSE; // not one track passes the cuts
    for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
      if(_tempSwapN.GetBinContent(l+1) == 0) {
        _tempSwap.SetBinContent(l+1,0);
        _tempSwap.SetBinError(l+1,0);
      }
      else {
        Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
        Double_t variance = vartimesnsq/(totalns*(totalns-1.));
        Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
        Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
        Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
        Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
        Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
        if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
        else {
          _tempSwap.SetBinContent(l+1,0);
          _tempSwap.SetBinError(l+1,0);
        }
      }
    }
  }

  fFitModulation->SetParameter(0, fLocalRho->GetVal());
  switch (fFitModulationType) {
  case kNoFit : { fFitModulation->FixParameter(0, fLocalRho->GetVal() ); 
  } break;
  case kV2 : { 
    fFitModulation->FixParameter(4, psi2); 
  } break;
  case kV3 : { 
    fFitModulation->FixParameter(4, psi3); 
  } break;
  case kCombined : {
    fFitModulation->FixParameter(4, psi2); 
    fFitModulation->FixParameter(6, psi3);
  } break;
  case kFourierSeries : {
    // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
    // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
    Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
    for(Int_t i(0); i < iTracks; i++) {
      AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
      if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
      sumPt += track->Pt();
      cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2)); 
      sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
      cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3)); 
      sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
    }
    fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
    fFitModulation->SetParameter(4, psi2);
    fFitModulation->SetParameter(6, psi3);
    fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
  } break;
  default : break;
  }
  _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", 0, TMath::TwoPi());
  // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
  Double_t CDF(1.-ChiSquareCDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
  if(fFillHistograms) fHistPvalueCDF->Fill(CDF);
  if(CDF > fMinPvalue && CDF < fMaxPvalue && ( fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality
    if(fFillHistograms) fHistRhoStatusCent->Fill(fCent, 0.);
    // for LOCAL didactic purposes, save the  best and the worst fits
    // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID 
    // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
    switch (fRunModeType) {
    case kLocal : {
      if(gRandom->Uniform(0, 100) > fPercentageOfFits) break;
      static Int_t didacticCounterBest(0);
      TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
      TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
      didacticProfile->GetListOfFunctions()->Add(didactifFit);
      fOutputListGood->Add(didacticProfile);
      didacticCounterBest++;
      TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
      for(Int_t i(0); i < iTracks; i++) {
	AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
	if(PassesCuts(track)) {
	  if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
	  else didacticSurface->Fill(track->Phi(), track->Eta());
	}
      }
      if(fExcludeLeadingJetsFromFit) {       // visualize the excluded region
	TF2 *f2 = new TF2(Form("%s_LJ", didacticSurface->GetName()),"[0]*TMath::Gaus(x,[1],[2])*TMath::Gaus(y,[3],[4])", 0, TMath::TwoPi(), -1, 1);
	f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
	didacticSurface->GetListOfFunctions()->Add(f2);
      }
      fOutputListGood->Add(didacticSurface);
    } break;
    default : break;
    }
  } else {    // if the fit is of poor quality revert to the original rho estimate
    if(fFillHistograms) fHistRhoStatusCent->Fill(fCent, 1.);
    switch (fRunModeType) { // again see if we want to save the fit
    case kLocal : {
      static Int_t didacticCounterWorst(0);
      if(gRandom->Uniform(0, 100) > fPercentageOfFits) break;
      TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
      TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
      didacticProfile->GetListOfFunctions()->Add(didactifFit);
      fOutputListBad->Add(didacticProfile);
      didacticCounterWorst++;
    } break;
    default : break;
    }
    switch (fFitModulationType) {
    case kNoFit : break;        // nothing to do
    case kCombined : fFitModulation->SetParameter(7, 0);        // no break
    case kFourierSeries : fFitModulation->SetParameter(7, 0);   // no break
    default : { // needs to be done if there was a poor fit
      fFitModulation->SetParameter(3, 0);
      fFitModulation->SetParameter(0, fLocalRho->GetVal());
    } break;
    }
    return kFALSE;  // return false if the fit is rejected
  }
  return kTRUE;
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::FillAnalysisSummaryHistogram() const
{
  // Fill the analysis summary histrogram, saves all relevant analysis settigns

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
  fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
  fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
  fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
  fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
  fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
  fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
  fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
  fHistAnalysisSummary->SetBinContent(37, 1.);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
  fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
  fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
  fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
  fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
  fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fLocalJetMinEta");
  fHistAnalysisSummary->SetBinContent(45,fLocalJetMinEta );
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fLocalJetMaxEta");
  fHistAnalysisSummary->SetBinContent(46, fLocalJetMaxEta);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "fLocalJetMinPhi");
  fHistAnalysisSummary->SetBinContent(47, fLocalJetMinPhi);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "fLocalJetMaxPhi");
  fHistAnalysisSummary->SetBinContent(48, fLocalJetMaxPhi);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(49, "fSoftTrackMinPt");
  fHistAnalysisSummary->SetBinContent(49, fSoftTrackMinPt);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(50, "fSoftTrackMaxPt");
  fHistAnalysisSummary->SetBinContent(50, fSoftTrackMaxPt);
  fHistAnalysisSummary->GetXaxis()->SetBinLabel(51, "fUseScaledRho");
  fHistAnalysisSummary->SetBinContent(51, fUseScaledRho);
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::FillEventPlaneHistograms(Double_t psi2, Double_t psi3) const
{
  // Fill event plane histograms

  if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
  fHistPsi2[fInCentralitySelection]->Fill(psi2);
  fHistPsi3[fInCentralitySelection]->Fill(psi3);    
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::Terminate(Option_t *)
{
  // Terminate
}

//_____________________________________________________________________________
void AliAnalysisTaskLocalRho::SetModulationFit(TF1* fit) 
{
  // Set function to fit modulation

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