ROOT logo
/**************************************************************************
 * Copyright(c) 1998-2010, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/* $Id$ */

//***********************************************************************
// Class AliHFPtSpectrum
// Base class for feed-down corrections on heavy-flavour decays
// computes the cross-section via one of the three implemented methods:
//   0) Consider no feed-down prediction 
//   1) Subtract the feed-down with the "fc" method 
//       Yield = Reco * fc;  where fc = 1 / ( 1 + (eff_b/eff_c)*(N_b/N_c) ) ;
//   2) Subtract the feed-down with the "Nb" method
//       Yield = Reco - Feed-down (exact formula on the function implementation)
//
//  (the corrected yields per bin are divided by the bin-width)
//
//
//  In HIC you can also evaluate how the feed-down correction is influenced by an energy loss hypothesis: 
//      Raa(c-->D) / Raa(b-->D) defined here as Rcb for the "fc" method
//      Raa(b-->D) defined here as Rb for the "Nb" method
//
// Author: Z.Conesa, zconesa@in2p3.fr
//***********************************************************************

#include <Riostream.h>

#include "TMath.h"
#include "TH1.h"
#include "TH1D.h"
#include "TH2.h"
#include "TH2D.h"
#include "TNtuple.h"
#include "TGraphAsymmErrors.h"
#include "TNamed.h"
#include "TCanvas.h"
#include "TLegend.h"

//#include "AliLog.h"
#include "AliHFSystErr.h"
#include "AliHFPtSpectrum.h"

ClassImp(AliHFPtSpectrum)

//_________________________________________________________________________________________________________
AliHFPtSpectrum::AliHFPtSpectrum(const char* name, const char* title, Int_t option):
  TNamed(name,title),
  fhDirectMCpt(NULL),
  fhFeedDownMCpt(NULL),
  fhDirectMCptMax(NULL),
  fhDirectMCptMin(NULL),
  fhFeedDownMCptMax(NULL),
  fhFeedDownMCptMin(NULL),
  fhDirectEffpt(NULL),
  fhFeedDownEffpt(NULL),
  fhRECpt(NULL),
  fgRECSystematics(NULL),
  fNevts(1),
  fLuminosity(),
  fTrigEfficiency(),
  fGlobalEfficiencyUncertainties(),
  fTab(),
  fhFc(NULL),
  fhFcMax(NULL),
  fhFcMin(NULL),
  fhFcRcb(NULL),
  fgFcExtreme(NULL),
  fgFcConservative(NULL),
  fhYieldCorr(NULL),
  fhYieldCorrMax(NULL),
  fhYieldCorrMin(NULL),
  fhYieldCorrRcb(NULL),
  fgYieldCorr(NULL),
  fgYieldCorrExtreme(NULL),
  fgYieldCorrConservative(NULL),
  fhSigmaCorr(NULL),
  fhSigmaCorrMax(NULL),
  fhSigmaCorrMin(NULL),
  fhSigmaCorrDataSyst(NULL),
  fhSigmaCorrRcb(NULL),
  fgSigmaCorr(NULL),
  fgSigmaCorrExtreme(NULL),
  fgSigmaCorrConservative(NULL),
  fnSigma(NULL),
  fnHypothesis(NULL),
  fFeedDownOption(option),
  fAsymUncertainties(kTRUE),
  fPbPbElossHypothesis(kFALSE),
  fIsStatUncEff(kTRUE),
  fParticleAntiParticle(2),
  fIsEventPlane(kFALSE),
  fnPtBins(0),
  fPtBinLimits(NULL),
  fPtBinWidths(NULL),
  fhStatUncEffcSigma(NULL),
  fhStatUncEffbSigma(NULL),
  fhStatUncEffcFD(NULL),
  fhStatUncEffbFD(NULL)
{
  //
  // Default constructor
  //

  fLuminosity[0]=1.;  fLuminosity[1]=0.;  
  fTrigEfficiency[0]=1.; fTrigEfficiency[1]=0.; 
  fGlobalEfficiencyUncertainties[0]=0.; fGlobalEfficiencyUncertainties[1]=0.;
  fTab[0]=1.;  fTab[1]=0.;

}

//_________________________________________________________________________________________________________
AliHFPtSpectrum::AliHFPtSpectrum(const AliHFPtSpectrum &rhs):
  TNamed(rhs),
  fhDirectMCpt(rhs.fhDirectMCpt),
  fhFeedDownMCpt(rhs.fhFeedDownMCpt),
  fhDirectMCptMax(rhs.fhDirectMCptMax),
  fhDirectMCptMin(rhs.fhDirectMCptMin),
  fhFeedDownMCptMax(rhs.fhFeedDownMCptMax),
  fhFeedDownMCptMin(rhs.fhFeedDownMCptMin),
  fhDirectEffpt(rhs.fhDirectEffpt),
  fhFeedDownEffpt(rhs.fhFeedDownEffpt),
  fhRECpt(rhs.fhRECpt),
  fgRECSystematics(rhs.fgRECSystematics),
  fNevts(rhs.fNevts),
  fLuminosity(),
  fTrigEfficiency(),
  fGlobalEfficiencyUncertainties(),
  fTab(),
  fhFc(rhs.fhFc),
  fhFcMax(rhs.fhFcMax),
  fhFcMin(rhs.fhFcMin),
  fhFcRcb(rhs.fhFcRcb),
  fgFcExtreme(rhs.fgFcExtreme),
  fgFcConservative(rhs.fgFcConservative),
  fhYieldCorr(rhs.fhYieldCorr),
  fhYieldCorrMax(rhs.fhYieldCorrMax),
  fhYieldCorrMin(rhs.fhYieldCorrMin),
  fhYieldCorrRcb(rhs.fhYieldCorrRcb),
  fgYieldCorr(rhs.fgYieldCorr),
  fgYieldCorrExtreme(rhs.fgYieldCorrExtreme),
  fgYieldCorrConservative(rhs.fgYieldCorrConservative),
  fhSigmaCorr(rhs.fhSigmaCorr),
  fhSigmaCorrMax(rhs.fhSigmaCorrMax),
  fhSigmaCorrMin(rhs.fhSigmaCorrMin),
  fhSigmaCorrDataSyst(rhs.fhSigmaCorrDataSyst),
  fhSigmaCorrRcb(rhs.fhSigmaCorrRcb),
  fgSigmaCorr(rhs.fgSigmaCorr),
  fgSigmaCorrExtreme(rhs.fgSigmaCorrExtreme),
  fgSigmaCorrConservative(rhs.fgSigmaCorrConservative),
  fnSigma(rhs.fnSigma),
  fnHypothesis(rhs.fnHypothesis),
  fFeedDownOption(rhs.fFeedDownOption),
  fAsymUncertainties(rhs.fAsymUncertainties),
  fPbPbElossHypothesis(rhs.fPbPbElossHypothesis),
  fIsStatUncEff(rhs.fIsStatUncEff),
  fParticleAntiParticle(rhs.fParticleAntiParticle),
  fIsEventPlane(rhs.fIsEventPlane),
  fnPtBins(rhs.fnPtBins),
  fPtBinLimits(),
  fPtBinWidths(),
  fhStatUncEffcSigma(NULL),
  fhStatUncEffbSigma(NULL),
  fhStatUncEffcFD(NULL),
  fhStatUncEffbFD(NULL)
{
  //
  // Copy constructor
  //

  for(Int_t i=0; i<2; i++){
    fLuminosity[i] = rhs.fLuminosity[i];
    fTrigEfficiency[i] = rhs.fTrigEfficiency[i];
    fGlobalEfficiencyUncertainties[i] = rhs.fGlobalEfficiencyUncertainties[i];
    fTab[i] = rhs.fTab[i];
  }

  for(Int_t i=0; i<fnPtBins; i++){
    fPtBinLimits[i] = rhs.fPtBinLimits[i];
    fPtBinWidths[i] = rhs.fPtBinWidths[i];
  }
  fPtBinLimits[fnPtBins] = rhs.fPtBinLimits[fnPtBins];

}

//_________________________________________________________________________________________________________
AliHFPtSpectrum &AliHFPtSpectrum::operator=(const AliHFPtSpectrum &source){
  //
  // Assignment operator
  //

  if (&source == this) return *this;
  
  fhDirectMCpt = source.fhDirectMCpt;
  fhFeedDownMCpt = source.fhFeedDownMCpt;
  fhDirectMCptMax = source.fhDirectMCptMax;
  fhDirectMCptMin = source.fhDirectMCptMin;
  fhFeedDownMCptMax = source.fhFeedDownMCptMax;
  fhFeedDownMCptMin = source.fhFeedDownMCptMin;
  fhDirectEffpt = source.fhDirectEffpt;
  fhFeedDownEffpt = source.fhFeedDownEffpt;
  fhRECpt = source.fhRECpt;
  fgRECSystematics = source.fgRECSystematics;
  fNevts = source.fNevts;
  fhFc = source.fhFc;
  fhFcMax = source.fhFcMax;
  fhFcMin = source.fhFcMin;
  fhFcRcb = source.fhFcRcb;
  fgFcExtreme = source.fgFcExtreme;
  fgFcConservative = source.fgFcConservative;
  fhYieldCorr = source.fhYieldCorr;
  fhYieldCorrMax = source.fhYieldCorrMax;
  fhYieldCorrMin = source.fhYieldCorrMin;
  fhYieldCorrRcb = source.fhYieldCorrRcb;
  fgYieldCorr = source.fgYieldCorr;
  fgYieldCorrExtreme = source.fgYieldCorrExtreme;
  fgYieldCorrConservative = source.fgYieldCorrConservative;
  fhSigmaCorr = source.fhSigmaCorr;
  fhSigmaCorrMax = source.fhSigmaCorrMax;
  fhSigmaCorrMin = source.fhSigmaCorrMin;
  fhSigmaCorrDataSyst = source.fhSigmaCorrDataSyst;
  fhSigmaCorrRcb = source.fhSigmaCorrRcb;
  fgSigmaCorr = source.fgSigmaCorr;
  fgSigmaCorrExtreme = source.fgSigmaCorrExtreme;
  fgSigmaCorrConservative = source.fgSigmaCorrConservative;
  fnSigma = source.fnSigma;
  fnHypothesis = source.fnHypothesis;
  fFeedDownOption = source.fFeedDownOption;
  fAsymUncertainties = source.fAsymUncertainties;
  fPbPbElossHypothesis = source.fPbPbElossHypothesis;
  fIsStatUncEff = source.fIsStatUncEff;
  fParticleAntiParticle = source.fParticleAntiParticle;
  fIsEventPlane = source.fIsEventPlane;
  
  for(Int_t i=0; i<2; i++){
    fLuminosity[i] = source.fLuminosity[i];
    fTrigEfficiency[i] = source.fTrigEfficiency[i];
    fGlobalEfficiencyUncertainties[i] = source.fGlobalEfficiencyUncertainties[i];
    fTab[i] = source.fTab[i];
  }

  fnPtBins = source.fnPtBins;
  for(Int_t i=0; i<fnPtBins; i++){
    fPtBinLimits[i] = source.fPtBinLimits[i];
    fPtBinWidths[i] = source.fPtBinWidths[i];
  }
  fPtBinLimits[fnPtBins] = source.fPtBinLimits[fnPtBins];

  return *this;
}

//_________________________________________________________________________________________________________
AliHFPtSpectrum::~AliHFPtSpectrum(){
  //
  // Destructor
  //
  if (fhDirectMCpt) delete fhDirectMCpt;    
  if (fhFeedDownMCpt) delete fhFeedDownMCpt;  
  if (fhDirectMCptMax) delete fhDirectMCptMax; 
  if (fhDirectMCptMin) delete fhDirectMCptMin; 
  if (fhFeedDownMCptMax) delete fhFeedDownMCptMax;
  if (fhFeedDownMCptMin) delete fhFeedDownMCptMin;
  if (fhDirectEffpt) delete fhDirectEffpt;    
  if (fhFeedDownEffpt) delete fhFeedDownEffpt;  
  if (fhRECpt) delete fhRECpt;    
  if (fgRECSystematics) delete fgRECSystematics;
  if (fhFc) delete fhFc;
  if (fhFcMax) delete fhFcMax;
  if (fhFcMin) delete fhFcMin;
  if (fhFcRcb) delete fhFcRcb;
  if (fgFcExtreme) delete fgFcExtreme;  
  if (fgFcConservative) delete fgFcConservative;
  if (fhYieldCorr) delete fhYieldCorr;                
  if (fhYieldCorrMax) delete fhYieldCorrMax;             
  if (fhYieldCorrMin) delete fhYieldCorrMin;    
  if (fhYieldCorrRcb) delete fhYieldCorrRcb;
  if (fgYieldCorr) delete fgYieldCorr;  
  if (fgYieldCorrExtreme) delete fgYieldCorrExtreme;
  if (fgYieldCorrConservative) delete fgYieldCorrConservative;
  if (fhSigmaCorr) delete fhSigmaCorr;                  
  if (fhSigmaCorrMax) delete fhSigmaCorrMax;               
  if (fhSigmaCorrMin) delete fhSigmaCorrMin; 
  if (fhSigmaCorrDataSyst) delete fhSigmaCorrDataSyst;
  if (fgSigmaCorr) delete fgSigmaCorr;    
  if (fgSigmaCorrExtreme) delete fgSigmaCorrExtreme;
  if (fgSigmaCorrConservative) delete fgSigmaCorrConservative;
  if (fnSigma) delete fnSigma;
  if (fnHypothesis) delete fnHypothesis;
  if (fPtBinLimits) delete [] fPtBinLimits;
  if (fPtBinWidths) delete [] fPtBinWidths;
}
  

//_________________________________________________________________________________________________________
TH1D * AliHFPtSpectrum::RebinTheoreticalSpectra(TH1D *hTheory, const char *name) {
  //
  // Function to rebin the theoretical spectrum 
  //  with respect to the real-data reconstructed spectrum binning 
  //
  
  if (!hTheory || !fhRECpt) {
    AliError("Feed-down or reconstructed spectra don't exist");
    return NULL;
  }

  //
  // Get the reconstructed spectra bins & limits
  Int_t nbinsMC = hTheory->GetNbinsX();

  // Check that the reconstructed spectra binning 
  // is larger than the theoretical one
  Double_t thbinwidth = hTheory->GetBinWidth(1);
  for (Int_t i=1; i<=fnPtBins; i++) {
    if ( thbinwidth > fPtBinWidths[i-1] ) {
      AliInfo(" Beware it seems that the reconstructed spectra has a smaller binning than the theoretical predictions !! ");
    }
  }

  //
  // Define a new histogram with the real-data reconstructed spectrum binning 
  TH1D * hTheoryRebin = new TH1D(name," theoretical rebinned prediction",fnPtBins,fPtBinLimits);

  Double_t sum[fnPtBins], items[fnPtBins];
  for (Int_t ibin=0; ibin<fnPtBins; ibin++) {
    sum[ibin]=0.; items[ibin]=0.;
  }
  for (Int_t ibin=0; ibin<=nbinsMC; ibin++){
    
    for (Int_t ibinrec=0; ibinrec<fnPtBins; ibinrec++){
      if (hTheory->GetBinCenter(ibin)>fPtBinLimits[ibinrec] &&
	  hTheory->GetBinCenter(ibin)<fPtBinLimits[ibinrec+1]){
	sum[ibinrec]+=hTheory->GetBinContent(ibin);
	items[ibinrec]+=1.;
      }
    }
    
  }

  // set the theoretical rebinned spectra to ( sum-bins / n-bins ) per new bin
  for (Int_t ibinrec=0; ibinrec<fnPtBins; ibinrec++) {
    hTheoryRebin->SetBinContent(ibinrec+1,sum[ibinrec]/items[ibinrec]);
  }
  
  return (TH1D*)hTheoryRebin;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetMCptSpectra(TH1D *hDirect, TH1D *hFeedDown){
  //
  // Set the MonteCarlo or Theoretical spectra
  //  both for direct and feed-down contributions
  //
  
  if (!hDirect || !hFeedDown || !fhRECpt) {
    AliError("One or both (direct, feed-down) spectra or the reconstructed spectra don't exist");
    return;
  }

  Bool_t areconsistent = kTRUE;
  areconsistent = CheckHistosConsistency(hDirect,hFeedDown);
  if (!areconsistent) {
    AliInfo("Histograms are not consistent (bin width, bounds)"); 
    return;
  }

  //
  // Rebin the theoretical predictions to the reconstructed spectra binning 
  //
  fhDirectMCpt = RebinTheoreticalSpectra(hDirect,"fhDirectMCpt");
  fhDirectMCpt->SetNameTitle("fhDirectMCpt"," direct theoretical prediction");
  fhFeedDownMCpt = RebinTheoreticalSpectra(hFeedDown,"fhFeedDownMCpt");
  fhFeedDownMCpt->SetNameTitle("fhFeedDownMCpt"," feed-down theoretical prediction");

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetFeedDownMCptSpectra(TH1D *hFeedDown){
  //
  // Set the MonteCarlo or Theoretical spectra
  //  for feed-down contribution
  //
  
  if (!hFeedDown || !fhRECpt) {
    AliError("Feed-down or reconstructed spectra don't exist");
    return;
  }

  //
  // Rebin the theoretical predictions to the reconstructed spectra binning 
  //
  fhFeedDownMCpt = RebinTheoreticalSpectra(hFeedDown,"fhFeedDownMCpt");
  fhFeedDownMCpt->SetNameTitle("fhFeedDownMCpt"," feed-down theoretical prediction");

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetMCptDistributionsBounds(TH1D *hDirectMax, TH1D *hDirectMin, TH1D *hFeedDownMax, TH1D *hFeedDownMin){
  //
  // Set the maximum and minimum MonteCarlo or Theoretical spectra
  //  both for direct and feed-down contributions
  // used in case uncertainties are asymmetric and ca not be on the "basic histograms"
  //

  if (!hDirectMax || !hDirectMin || !hFeedDownMax|| !hFeedDownMin || !fhRECpt) {
    AliError("One or all of the max/min direct/feed-down or the reconstructed spectra don't exist");
    return;
  }

  Bool_t areconsistent = kTRUE; 
  areconsistent &= CheckHistosConsistency(hDirectMax,hDirectMin);
  areconsistent &= CheckHistosConsistency(hFeedDownMax,hFeedDownMin);
  areconsistent &= CheckHistosConsistency(hDirectMax,hFeedDownMax);
  if (!areconsistent) {
    AliInfo("Histograms are not consistent (bin width, bounds)"); 
    return;
  }


  //
  // Rebin the theoretical predictions to the reconstructed spectra binning 
  //
  fhDirectMCptMax = RebinTheoreticalSpectra(hDirectMax,"fhDirectMCptMax");
  fhDirectMCptMax->SetNameTitle("fhDirectMCptMax"," maximum direct theoretical prediction");
  fhDirectMCptMin = RebinTheoreticalSpectra(hDirectMin,"fhDirectMCptMin");
  fhDirectMCptMin->SetNameTitle("fhDirectMCptMin"," minimum direct theoretical prediction");
  fhFeedDownMCptMax = RebinTheoreticalSpectra(hFeedDownMax,"fhFeedDownMCptMax");
  fhFeedDownMCptMax->SetNameTitle("fhFeedDownMCptMax"," maximum feed-down theoretical prediction");
  fhFeedDownMCptMin = RebinTheoreticalSpectra(hFeedDownMin,"fhFeedDownMCptMin");
  fhFeedDownMCptMin->SetNameTitle("fhFeedDownMCptMin"," minimum feed-down theoretical prediction");

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetFeedDownMCptDistributionsBounds(TH1D *hFeedDownMax, TH1D *hFeedDownMin){
  //
  // Set the maximum and minimum MonteCarlo or Theoretical spectra
  //   for feed-down contributions
  // used in case uncertainties are asymmetric and can not be on the "basic histogram"
  //

  if (!hFeedDownMax || !hFeedDownMin || !fhRECpt) {
    AliError("One or all of the max/min direct/feed-down spectra don't exist");
    return;
  }

  Bool_t areconsistent = kTRUE; 
  areconsistent &= CheckHistosConsistency(hFeedDownMax,hFeedDownMin);
  if (!areconsistent) {
    AliInfo("Histograms are not consistent (bin width, bounds)"); 
    return;
  }


  //
  // Rebin the theoretical predictions to the reconstructed spectra binning 
  //
  fhFeedDownMCptMax = RebinTheoreticalSpectra(hFeedDownMax,"fhFeedDownMCptMax");
  fhFeedDownMCptMax->SetNameTitle("fhFeedDownMCptMax"," maximum feed-down theoretical prediction");
  fhFeedDownMCptMin = RebinTheoreticalSpectra(hFeedDownMin,"fhFeedDownMCptMin");
  fhFeedDownMCptMin->SetNameTitle("fhFeedDownMCptMin"," minimum feed-down theoretical prediction");

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetDirectAccEffCorrection(TH1D *hDirectEff){
  //
  // Set the Acceptance and Efficiency corrections 
  //   for the direct contribution
  //
  
  if (!hDirectEff) {
    AliError("The direct acceptance and efficiency corrections doesn't exist");
    return;
  }

  fhDirectEffpt = (TH1D*)hDirectEff->Clone();
  fhDirectEffpt->SetNameTitle("fhDirectEffpt"," direct acceptance x efficiency correction");
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetAccEffCorrection(TH1D *hDirectEff, TH1D *hFeedDownEff){
  //
  // Set the Acceptance and Efficiency corrections 
  //  both for direct and feed-down contributions
  //
  
  if (!hDirectEff || !hFeedDownEff) {
    AliError("One or both (direct, feed-down) acceptance and efficiency corrections don't exist");
    return;
  }

  Bool_t areconsistent=kTRUE;
  areconsistent = CheckHistosConsistency(hDirectEff,hFeedDownEff);
  if (!areconsistent) {
    AliInfo("Histograms are not consistent (bin width, bounds)"); 
    return;
  }

  fhDirectEffpt = (TH1D*)hDirectEff->Clone();
  fhFeedDownEffpt = (TH1D*)hFeedDownEff->Clone();
  fhDirectEffpt->SetNameTitle("fhDirectEffpt"," direct acceptance x efficiency correction");
  fhFeedDownEffpt->SetNameTitle("fhFeedDownEffpt"," feed-down acceptance x efficiency correction");
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetReconstructedSpectrum(TH1D *hRec) {
  //
  // Set the reconstructed spectrum
  //
  
  if (!hRec) {
    AliError("The reconstructed spectrum doesn't exist");
    return;
  }

  fhRECpt = (TH1D*)hRec->Clone();
  fhRECpt->SetNameTitle("fhRECpt"," reconstructed spectrum");

  //
  // Evaluate the number of intervals and limits of the spectrum
  //
  fnPtBins = fhRECpt->GetNbinsX();
  fPtBinLimits = new Double_t[fnPtBins+1];
  fPtBinWidths = new Double_t[fnPtBins];
  Double_t xlow=0., binwidth=0.;
  for (Int_t i=1; i<=fnPtBins; i++) {
    binwidth = fhRECpt->GetBinWidth(i);
    xlow = fhRECpt->GetBinLowEdge(i);
    fPtBinLimits[i-1] = xlow;
    fPtBinWidths[i-1] = binwidth;
  }
  fPtBinLimits[fnPtBins] = xlow + binwidth;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::SetReconstructedSpectrumSystematics(TGraphAsymmErrors *gRec) {
  //
  // Set the reconstructed spectrum (uncorrected yield) systematic uncertainties
  // 

  // Check the compatibility with the reconstructed spectrum
  Double_t gbinwidth = gRec->GetErrorXlow(1) + gRec->GetErrorXhigh(1) ;
  Double_t hbinwidth = fhRECpt->GetBinWidth(1);
  Double_t gxbincenter=0., gybincenter=0.;
  gRec->GetPoint(1,gxbincenter,gybincenter);
  Double_t hbincenter = fhRECpt->GetBinCenter(1);
  if ( (gbinwidth != hbinwidth) || (gxbincenter!=hbincenter) ) {
    AliError(" The reconstructed spectrum and its systematics don't seem compatible");
    return;
  }
  
  fgRECSystematics = gRec; 
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::ComputeHFPtSpectrum(Double_t deltaY, Double_t branchingRatioC, Double_t branchingRatioBintoFinalDecay) {
  //
  // Main function to compute the corrected cross-section:
  // variables : analysed delta_y, BR for the final correction,
  //             BR b --> D --> decay (relative to the input theoretical prediction)
  //
  //   Sigma = ( 1. / (lumi * delta_y * BR_c * ParticleAntiPartFactor * eff_trig * eff_c ) ) * spectra (corrected for feed-down)
  //
  // Uncertainties: (stat) delta_sigma = sigma * sqrt ( (delta_spectra/spectra)^2 )
  //  (syst but feed-down) delta_sigma = sigma * sqrt ( (delta_spectra_syst/spectra)^2 + (delta_lumi/lumi)^2 + (delta_eff_trig/eff_trig)^2 + (delta_eff/eff)^2 )
  //      (feed-down syst) delta_sigma = sigma * sqrt ( (delta_spectra_fd/spectra_fd)^2 )
  //
  //  In HIC the feed-down correction varies with an energy loss hypothesis:
  //      Raa(c-->D) / Raa(b-->D) for the "fc" method, Raa(b-->D) for the "Nb" method (see exact formulas in the functions)
  //

  //
  // First: Initialization
  //
  Bool_t areHistosOk = Initialize();
  if (!areHistosOk) {
    AliInfo(" Histos not properly initialized. Check : inconsistent binning ? missing histos ?");
    return;
  }
  // Reset the statistical uncertainties on the efficiencies if needed
  if(!fIsStatUncEff) ResetStatUncEff();

  //
  // Second: Correct for feed-down
  //
  if (fFeedDownOption==1) {
    // Compute the feed-down correction via fc-method
    CalculateFeedDownCorrectionFc(); 
    // Correct the yield for feed-down correction via fc-method
    CalculateFeedDownCorrectedSpectrumFc(); 
  }
  else if (fFeedDownOption==2) {
    // Correct the yield for feed-down correction via Nb-method
    CalculateFeedDownCorrectedSpectrumNb(deltaY,branchingRatioBintoFinalDecay); 
  }
  else if (fFeedDownOption==0) { 
    // If there is no need for feed-down correction,
    //    the "corrected" yield is equal to the raw yield
    CalculateCorrectedSpectrumNoFeedDown();
  }
  else { 
    AliInfo(" Are you sure the feed-down correction option is right ?"); 
  }


  // Print out information
  printf("\n\n     Correcting the spectra with : \n   luminosity = %2.2e +- %2.2e, trigger efficiency = %2.2e +- %2.2e, \n    delta_y = %2.2f, BR_c = %2.2e, BR_b_decay = %2.2e \n    %2.2f percent uncertainty on the efficiencies, and %2.2f percent uncertainty on the b/c efficiencies ratio \n\n",fLuminosity[0],fLuminosity[1],fTrigEfficiency[0],fTrigEfficiency[1],deltaY,branchingRatioC,branchingRatioBintoFinalDecay,fGlobalEfficiencyUncertainties[0],fGlobalEfficiencyUncertainties[1]);
  if (fPbPbElossHypothesis)  printf("\n\n     The considered Tab is  %4.2e +- %2.2e \n\n",fTab[0],fTab[1]);

  //
  // Finally: Correct from yields to cross-section
  //
  
  // declare the output histograms
  fhSigmaCorr = new TH1D("fhSigmaCorr","corrected sigma",fnPtBins,fPtBinLimits);
  fgSigmaCorr = new TGraphAsymmErrors(fnPtBins+1);

  fhSigmaCorrMax = new TH1D("fhSigmaCorrMax","max corrected sigma",fnPtBins,fPtBinLimits);
  fhSigmaCorrMin = new TH1D("fhSigmaCorrMin","min corrected sigma",fnPtBins,fPtBinLimits);
  fhSigmaCorrDataSyst = new TH1D("fhSigmaCorrDataSyst","data syst uncertainties on the corrected sigma",fnPtBins,fPtBinLimits);

  if (fPbPbElossHypothesis && fFeedDownOption==1) {
    fhSigmaCorrRcb = new TH2D("fhSigmaCorrRcb","corrected sigma vs Rcb Eloss hypothesis; p_{T} [GeV/c] ; Rcb Eloss hypothesis ; #sigma",fnPtBins,fPtBinLimits,800,0.,4.);
    fnSigma = new TNtuple("fnSigma"," Sigma ntuple calculation","pt:Signal:Rcb:fc:Yield:Sigma:SigmaStatUnc:SigmaMax:SigmaMin");
  }
  if (fPbPbElossHypothesis && fFeedDownOption==2) {
    fhSigmaCorrRcb = new TH2D("fhSigmaCorrRcb","corrected sigma vs Rb Eloss hypothesis; p_{T} [GeV/c] ; Rb Eloss hypothesis ; #sigma",fnPtBins,fPtBinLimits,800,0.,4.);
    fnSigma = new TNtuple("fnSigma"," Sigma ntuple calculation","pt:Signal:Rb:fc:Yield:Sigma:SigmaStatUnc:SigmaMax:SigmaMin");
  }
  // and the output TGraphAsymmErrors
  if (fAsymUncertainties){
    fgSigmaCorrExtreme = new TGraphAsymmErrors(fnPtBins+1);
    fgSigmaCorrConservative = new TGraphAsymmErrors(fnPtBins+1);
  }
  fhStatUncEffcSigma = new TH1D("fhStatUncEffcSigma","direct charm stat unc on the cross section",fnPtBins,fPtBinLimits);
  fhStatUncEffbSigma = new TH1D("fhStatUncEffbSigma","secondary charm stat unc on the cross section",fnPtBins,fPtBinLimits);


  // protect against null denominator
  if (deltaY==0. || fLuminosity[0]==0. || fTrigEfficiency[0]==0. || branchingRatioC==0.) {
    AliError(" Hey you ! Why luminosity or trigger-efficiency or the c-BR or delta_y are set to zero ?! ");
    return ;
  }

  Double_t value=0, errValue=0, errvalueMax=0., errvalueMin=0.;
  Double_t errvalueExtremeMax=0., errvalueExtremeMin=0.;
  Double_t errvalueConservativeMax=0., errvalueConservativeMin=0.;
  Double_t errvalueStatUncEffc=0.;
  for(Int_t ibin=1; ibin<=fnPtBins; ibin++){
    
    // Variables initialization
    value=0.; errValue=0.; errvalueMax=0.; errvalueMin=0.;
    errvalueExtremeMax=0.; errvalueExtremeMin=0.;
    errvalueConservativeMax=0.; errvalueConservativeMin=0.;
    errvalueStatUncEffc=0.;

    // Sigma calculation
    //   Sigma = ( 1. / (lumi * delta_y * BR_c * ParticleAntiPartFactor * eff_trig * eff_c ) ) * spectra (corrected for feed-down)
    value = (fhDirectEffpt->GetBinContent(ibin) && fhDirectEffpt->GetBinContent(ibin)!=0. && fhRECpt->GetBinContent(ibin)>0.) ? 
      ( fhYieldCorr->GetBinContent(ibin) / ( deltaY * branchingRatioC * fParticleAntiParticle * fLuminosity[0] * fTrigEfficiency[0] * fhDirectEffpt->GetBinContent(ibin) ) )
      : 0. ;

    // Sigma statistical uncertainty:
    //   delta_sigma = sigma * sqrt ( (delta_spectra/spectra)^2 )
    errValue = (value!=0.) ?  value * (fhYieldCorr->GetBinError(ibin)/fhYieldCorr->GetBinContent(ibin))  : 0. ;

    //    cout<< " x "<< fhRECpt->GetBinCenter(ibin) << " sigma " << value << " +- "<< errValue << " (stat)"<<endl;

    //
    // Sigma systematic uncertainties
    //
    if (fAsymUncertainties && value>0.) {
      
      //  (syst but feed-down) delta_sigma = sigma * sqrt ( (delta_spectra_syst/spectra)^2 +
      //                                     (delta_lumi/lumi)^2 + (delta_eff_trig/eff_trig)^2 + (delta_eff/eff)^2  + (global_eff)^2 )
      errvalueMax = value * TMath::Sqrt( (fgYieldCorr->GetErrorYhigh(ibin)/fhYieldCorr->GetBinContent(ibin))*(fgYieldCorr->GetErrorYhigh(ibin)/fhYieldCorr->GetBinContent(ibin)) + 
					 (fLuminosity[1]/fLuminosity[0])*(fLuminosity[1]/fLuminosity[0]) + 
					 (fTrigEfficiency[1]/fTrigEfficiency[0])*(fTrigEfficiency[1]/fTrigEfficiency[0])  +
					 (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))*(fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin)) +
					 fGlobalEfficiencyUncertainties[0]*fGlobalEfficiencyUncertainties[0] );
      errvalueMin = value * TMath::Sqrt( (fgYieldCorr->GetErrorYlow(ibin)/fhYieldCorr->GetBinContent(ibin))*(fgYieldCorr->GetErrorYlow(ibin)/fhYieldCorr->GetBinContent(ibin)) + 
					 (fLuminosity[1]/fLuminosity[0])*(fLuminosity[1]/fLuminosity[0]) + 
					 (fTrigEfficiency[1]/fTrigEfficiency[0])*(fTrigEfficiency[1]/fTrigEfficiency[0])  +
					 (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))*(fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))  +
					 fGlobalEfficiencyUncertainties[0]*fGlobalEfficiencyUncertainties[0] );

      // Uncertainties from feed-down
      //      (feed-down syst) delta_sigma = sigma * sqrt ( (delta_spectra_fd/spectra_fd)^2 )
      //   extreme case
      errvalueExtremeMax = value * (fgYieldCorrExtreme->GetErrorYhigh(ibin)/fhYieldCorr->GetBinContent(ibin));
      errvalueExtremeMin =  value * (fgYieldCorrExtreme->GetErrorYlow(ibin)/fhYieldCorr->GetBinContent(ibin));
      //
      //   conservative case
      errvalueConservativeMax = value * (fgYieldCorrConservative->GetErrorYhigh(ibin)/fhYieldCorr->GetBinContent(ibin));
      errvalueConservativeMin =  value * (fgYieldCorrConservative->GetErrorYlow(ibin)/fhYieldCorr->GetBinContent(ibin));


      // stat unc of the efficiencies, separately
      errvalueStatUncEffc = value * (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin)) ;

    }
    else {
      // protect against null denominator
      errvalueMax = (value!=0.) ?
	value * TMath::Sqrt( (fLuminosity[1]/fLuminosity[0])*(fLuminosity[1]/fLuminosity[0]) + 
			     (fTrigEfficiency[1]/fTrigEfficiency[0])*(fTrigEfficiency[1]/fTrigEfficiency[0])  +
			     (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))*(fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))  +
			     fGlobalEfficiencyUncertainties[0]*fGlobalEfficiencyUncertainties[0] )
	: 0. ;
      errvalueMin = errvalueMax;
    }
    
    //
    // Fill the histograms
    //
    fhSigmaCorr->SetBinContent(ibin,value);
    fhSigmaCorr->SetBinError(ibin,errValue);

    Double_t x = fhYieldCorr->GetBinCenter(ibin);
    fgSigmaCorr->SetPoint(ibin,x,value); // i,x,y
    fgSigmaCorr->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueMin,errvalueMax); // i,xl,xh,yl,yh
    fhSigmaCorrMax->SetBinContent(ibin,value+errvalueMax);
    fhSigmaCorrMin->SetBinContent(ibin,value-errvalueMin);

    //
    // Fill the histos and ntuple vs the Eloss hypothesis
    //
    if (fPbPbElossHypothesis) {

      // Loop over the Eloss hypothesis
      if(!fnHypothesis) { 
	AliError("Error reading the fc hypothesis no ntuple, please check !!");
	return ;
      }
      Int_t entriesHypo = fnHypothesis->GetEntries();
      Float_t pt=0., Rhy=0., fc=0., fcMin=0., fcMax=0.;
      fnHypothesis->SetBranchAddress("pt",&pt);
      if (fFeedDownOption==2) fnHypothesis->SetBranchAddress("Rb",&Rhy);
      else if (fFeedDownOption==1) fnHypothesis->SetBranchAddress("Rcb",&Rhy);
      fnHypothesis->SetBranchAddress("fc",&fc);
      fnHypothesis->SetBranchAddress("fcMax",&fcMax);
      fnHypothesis->SetBranchAddress("fcMin",&fcMin);

      for (Int_t item=0; item<entriesHypo; item++) {

	fnHypothesis->GetEntry(item);
	if ( TMath::Abs( pt - fhDirectEffpt->GetBinCenter(ibin) ) > 0.15 ) continue;

	Double_t yieldRcbvalue = (fhRECpt->GetBinContent(ibin) ) ? fhRECpt->GetBinContent(ibin) * fc : 0. ;
	yieldRcbvalue /= fhRECpt->GetBinWidth(ibin) ;
	Double_t yieldRcbvalueMax = (fhRECpt->GetBinContent(ibin) ) ? fhRECpt->GetBinContent(ibin) * fcMax : 0. ;
	yieldRcbvalueMax /= fhRECpt->GetBinWidth(ibin) ;
	Double_t yieldRcbvalueMin = (fhRECpt->GetBinContent(ibin) ) ? fhRECpt->GetBinContent(ibin) * fcMin : 0. ;
	yieldRcbvalueMin /= fhRECpt->GetBinWidth(ibin) ;

	// Sigma calculation
	//   Sigma = ( 1. / (lumi * delta_y * BR_c * ParticleAntiPartFactor * eff_trig * eff_c ) ) * spectra (corrected for feed-down)
	Double_t sigmaRcbvalue = (fhDirectEffpt->GetBinContent(ibin) && fhDirectEffpt->GetBinContent(ibin)>0.) ? 
	  ( yieldRcbvalue / ( deltaY * branchingRatioC * fParticleAntiParticle * fLuminosity[0] * fTrigEfficiency[0] * fhDirectEffpt->GetBinContent(ibin) ) )
	  : 0. ;
	Double_t sigmaRcbvalueMax = (sigmaRcbvalue!=0.) ? 
	  ( yieldRcbvalueMax / ( deltaY * branchingRatioC * fParticleAntiParticle * fLuminosity[0] * fTrigEfficiency[0] * fhDirectEffpt->GetBinContent(ibin) ) )
	  : 0. ;
	Double_t sigmaRcbvalueMin = (sigmaRcbvalue!=0.) ? 
	  ( yieldRcbvalueMin / ( deltaY * branchingRatioC * fParticleAntiParticle * fLuminosity[0] * fTrigEfficiency[0] * fhDirectEffpt->GetBinContent(ibin) ) )
	  : 0. ;
	// Sigma statistical uncertainty:
	//   delta_sigma = sigma * sqrt ( (delta_spectra/spectra)^2 ) = sigma * ( delta_spectra / (spectra-corr * binwidth) )
	Double_t sigmaRcbvalueStatUnc = (sigmaRcbvalue!=0.) ? 
	  sigmaRcbvalue * ( fhRECpt->GetBinError(ibin) / ( yieldRcbvalue * fhRECpt->GetBinWidth(ibin) ) )  : 0. ;

	fhSigmaCorrRcb->Fill( fhSigmaCorr->GetBinCenter(ibin) , Rhy, sigmaRcbvalue );
	// 	if(ibin==3) 
	// 	  cout << " pt "<< fhRECpt->GetBinCenter(ibin) <<" bin "<< ibin<<" rval="<<rval<<", rbin="<<rbin<<" fc-value="<< fhFcRcb->GetBinContent(ibin,rbin) <<", yield-fcRbvalue="<<yieldRcbvalue<<", sigma-fcRbvalue="<<sigmaRcbvalue<<endl;
	fnSigma->Fill(fhRECpt->GetBinCenter(ibin), fhRECpt->GetBinContent(ibin),
		      Rhy, fc, yieldRcbvalue, sigmaRcbvalue, sigmaRcbvalueStatUnc, sigmaRcbvalueMax, sigmaRcbvalueMin );
      }
    }
    //
    // Fill the TGraphAsymmErrors
    if (fAsymUncertainties) {
      fgSigmaCorrExtreme->SetPoint(ibin,x,value); // i,x,y
      fgSigmaCorrExtreme->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueExtremeMin,errvalueExtremeMax); // i,xl,xh,yl,yh
      fgSigmaCorrConservative->SetPoint(ibin,x,value); // i,x,y
      fgSigmaCorrConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueConservativeMin,errvalueConservativeMax); // i,xl,xh,yl,yh

      fhStatUncEffcSigma->SetBinContent(ibin,0.); 
      if(value>0.) fhStatUncEffcSigma->SetBinError(ibin,((errvalueStatUncEffc/value)*100.));
      fhStatUncEffbSigma->SetBinContent(ibin,0.); fhStatUncEffbSigma->SetBinError(ibin,0.);
      //      cout << " pt "<< fhRECpt->GetBinCenter(ibin) <<" bin "<< ibin<<" stat-unc-c-sigma "<< errvalueStatUncEffc/value << endl;
    }
    
  }

}

//_________________________________________________________________________________________________________
TH1D * AliHFPtSpectrum::EstimateEfficiencyRecoBin(TH1D *hSimu, TH1D *hReco, const char *name) {
  //
  // Function that computes the acceptance and efficiency correction
  //  based on the simulated and reconstructed spectra
  //  and using the reconstructed spectra bin width
  //
  //  eff = reco/sim ; err_eff = sqrt( eff*(1-eff) )/ sqrt( sim )
  // 

  if(!fhRECpt){
    AliInfo("Hey, the reconstructed histogram was not set yet !"); 
    return NULL;
  }

  TH1D * hEfficiency = new TH1D(name," acceptance #times efficiency",fnPtBins,fPtBinLimits);
  
  Double_t *sumSimu=new Double_t[fnPtBins];
  Double_t *sumReco=new Double_t[fnPtBins];
  for (Int_t ibin=0; ibin<fnPtBins; ibin++){
    sumSimu[ibin]=0.;  sumReco[ibin]=0.;
  }
  for (Int_t ibin=0; ibin<=hSimu->GetNbinsX(); ibin++){

    for (Int_t ibinrec=0; ibinrec<fnPtBins; ibinrec++){
      if ( hSimu->GetBinCenter(ibin)>fPtBinLimits[ibinrec] &&
	   hSimu->GetBinCenter(ibin)<fPtBinLimits[ibinrec+1] ) {
	sumSimu[ibinrec]+=hSimu->GetBinContent(ibin);
      }
      if ( hReco->GetBinCenter(ibin)>fPtBinLimits[ibinrec] &&
	   hReco->GetBinCenter(ibin)<fPtBinLimits[ibinrec+1] ) {
	sumReco[ibinrec]+=hReco->GetBinContent(ibin);
      }
    }
    
  }


  // the efficiency is computed as reco/sim (in each bin)
  //  its uncertainty is err_eff = sqrt( eff*(1-eff) )/ sqrt( sim )
  Double_t eff=0., erreff=0.;
  for (Int_t ibinrec=0; ibinrec<fnPtBins; ibinrec++) {
    if (sumSimu[ibinrec]!= 0. && sumReco[ibinrec]!=0.) {
      eff = sumReco[ibinrec] / sumSimu[ibinrec] ;
      // protection in case eff > 1.0
      // test calculation (make the argument of the sqrt positive)
      erreff = TMath::Sqrt( eff * TMath::Abs(1.0 - eff) ) / TMath::Sqrt( sumSimu[ibinrec] );
    }
    else { eff=0.0; erreff=0.; }
    hEfficiency->SetBinContent(ibinrec+1,eff);
    hEfficiency->SetBinError(ibinrec+1,erreff);
  }

  delete [] sumSimu;
  delete [] sumReco;

  return (TH1D*)hEfficiency;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::EstimateAndSetDirectEfficiencyRecoBin(TH1D *hSimu, TH1D *hReco) {
  //
  // Function that computes the Direct  acceptance and efficiency correction
  //  based on the simulated and reconstructed spectra
  //  and using the reconstructed spectra bin width
  //
  //  eff = reco/sim ; err_eff = sqrt( eff*(1-eff) )/ sqrt( sim )
  // 

  if(!fhRECpt || !hSimu || !hReco){
    AliError("Hey, the reconstructed histogram was not set yet !"); 
    return;
  }

  fhDirectEffpt = EstimateEfficiencyRecoBin(hSimu,hReco,"fhDirectEffpt");
  fhDirectEffpt->SetNameTitle("fhDirectEffpt"," direct acceptance #times efficiency");

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::EstimateAndSetFeedDownEfficiencyRecoBin(TH1D *hSimu, TH1D *hReco) {
  //
  // Function that computes the Feed-Down acceptance and efficiency correction
  //  based on the simulated and reconstructed spectra
  //  and using the reconstructed spectra bin width
  //
  //  eff = reco/sim ; err_eff = sqrt( eff*(1-eff) )/ sqrt( sim )
  // 
  
  if(!fhRECpt || !hSimu || !hReco){
    AliError("Hey, the reconstructed histogram was not set yet !"); 
    return;
  }

  fhFeedDownEffpt = EstimateEfficiencyRecoBin(hSimu,hReco,"fhFeedDownEffpt");
  fhFeedDownEffpt->SetNameTitle("fhFeedDownEffpt"," feed-down acceptance #times efficiency");

}

//_________________________________________________________________________________________________________
Bool_t AliHFPtSpectrum::Initialize(){
  //
  // Initialization of the variables (histograms)
  //

  if (fFeedDownOption==0) { 
    AliInfo("Getting ready for the corrections without feed-down consideration");
  } else if (fFeedDownOption==1) { 
    AliInfo("Getting ready for the fc feed-down correction calculation");
  } else if (fFeedDownOption==2) {
    AliInfo("Getting ready for the Nb feed-down correction calculation");
  } else { AliError("The calculation option must be <=2"); return kFALSE; }

  // Start checking the input histograms consistency
  Bool_t areconsistent=kTRUE;

  // General checks 
  if (!fhDirectEffpt || !fhRECpt) {
    AliError(" Reconstructed spectra and/or the Nc efficiency distributions are not defined");
    return kFALSE;
  }
  areconsistent &= CheckHistosConsistency(fhRECpt,fhDirectEffpt);
  if (!areconsistent) {
    AliInfo("Histograms required for Nb correction are not consistent (bin width, bounds)"); 
    return kFALSE;
  }
  if (fFeedDownOption==0) return kTRUE;

  //
  // Common checks for options 1 (fc) & 2(Nb)
  if (!fhFeedDownMCpt || !fhFeedDownEffpt) {
    AliError(" Theoretical Nb and/or the Nb efficiency distributions are not defined");
    return kFALSE;
  }
  areconsistent &= CheckHistosConsistency(fhRECpt,fhFeedDownMCpt);
  areconsistent &= CheckHistosConsistency(fhFeedDownMCpt,fhFeedDownEffpt);
  if (fAsymUncertainties) {
    if (!fhFeedDownMCptMax || !fhFeedDownMCptMin) {
      AliError(" Max/Min theoretical Nb distributions are not defined");
      return kFALSE;
    }
    areconsistent &= CheckHistosConsistency(fhFeedDownMCpt,fhFeedDownMCptMax);
  }
  if (!areconsistent) {
    AliInfo("Histograms required for Nb correction are not consistent (bin width, bounds)"); 
    return kFALSE;
  }
  if (fFeedDownOption>1) return kTRUE;  

  //
  // Now checks for option 1 (fc correction) 
  if (!fhDirectMCpt) {
    AliError("Theoretical Nc distributions is not defined");
    return kFALSE;
  }
  areconsistent &= CheckHistosConsistency(fhDirectMCpt,fhFeedDownMCpt);
  areconsistent &= CheckHistosConsistency(fhDirectMCpt,fhDirectEffpt);
  if (fAsymUncertainties) {
    if (!fhDirectMCptMax || !fhDirectMCptMin) {
      AliError(" Max/Min theoretical Nc distributions are not defined");
      return kFALSE;
    }
    areconsistent &= CheckHistosConsistency(fhDirectMCpt,fhDirectMCptMax);
  }
  if (!areconsistent) {
    AliInfo("Histograms required for fc correction are not consistent (bin width, bounds)"); 
    return kFALSE;
  }

  return kTRUE;
}

//_________________________________________________________________________________________________________
Bool_t AliHFPtSpectrum::CheckHistosConsistency(TH1D *h1, TH1D *h2){
  //
  // Check the histograms consistency (bins, limits)
  //

  if (!h1 || !h2) {
    AliError("One or both histograms don't exist");
    return kFALSE;
  }

  Double_t binwidth1 = h1->GetBinWidth(1);
  Double_t binwidth2 = h2->GetBinWidth(1);
  Double_t min1 = h1->GetBinCenter(1) - (binwidth1/2.) ; 
//   Double_t max1 = h1->GetBinCenter(nbins1) + (binwidth1/2.) ;
  Double_t min2 = h2->GetBinCenter(1) - (binwidth2/2.) ;
//   Double_t max2 = h2->GetBinCenter(nbins2) + (binwidth2/2.) ;

  if (binwidth1!=binwidth2) {
    AliInfo(" histograms with different bin width");
    return kFALSE;
  }
  if (min1!=min2) {
    AliInfo(" histograms with different minimum");
    return kFALSE;
  }
//   if (max1!=max2) {
//     AliInfo(" histograms with different maximum");
//     return kFALSE;
//   }

  return kTRUE;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::CalculateCorrectedSpectrumNoFeedDown(){
  //
  // Compute the corrected spectrum with no feed-down correction
  //


  // declare the output histograms
  fhYieldCorr = new TH1D("fhYieldCorr","corrected yield (no feed-down corr)",fnPtBins,fPtBinLimits);
  fhYieldCorrMax = new TH1D("fhYieldCorrMax","max corrected yield (no feed-down corr)",fnPtBins,fPtBinLimits);
  fhYieldCorrMin = new TH1D("fhYieldCorrMin","min corrected yield (no feed-down corr)",fnPtBins,fPtBinLimits);
  fgYieldCorr = new TGraphAsymmErrors(fnPtBins+1);

  //
  // Do the calculation
  //
  for (Int_t ibin=1; ibin<=fnPtBins; ibin++) {
    Double_t value = fhRECpt->GetBinContent(ibin) /fhRECpt->GetBinWidth(ibin);
    Double_t errvalue = fhRECpt->GetBinError(ibin) /fhRECpt->GetBinWidth(ibin);
    fhYieldCorr->SetBinContent(ibin,value);
    fhYieldCorr->SetBinError(ibin,errvalue);
    fhYieldCorrMax->SetBinContent(ibin,value);
    fhYieldCorrMax->SetBinError(ibin,errvalue);
    fhYieldCorrMin->SetBinContent(ibin,value);
    fhYieldCorrMin->SetBinError(ibin,errvalue);
  }

  fAsymUncertainties=kFALSE;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::CalculateFeedDownCorrectionFc(){ 
  //
  // Compute fc factor and its uncertainties bin by bin
  //   fc = 1 / ( 1 + (eff_b/eff_c)*(N_b/N_c) ) 
  //
  // uncertainties: (conservative) combine the upper/lower N_b & N_c predictions together
  //                (extreme) combine the upper N_b predictions with the lower N_c predictions & viceversa
  //                systematic uncertainty on the acceptance x efficiency b/c ratio are included 
  //
  //  In addition, in HIC the feed-down correction varies with an energy loss hypothesis: Raa(c-->D) / Raa(b-->D) = Rcb
  //	       fc (Rcb) = ( 1. / ( 1 + (eff_b/eff_c)*(N_b/N_c)* (1/Rcb) ) );
  //
  AliInfo("Calculating the feed-down correction factor (fc method)");
  
  // define the variables
  Double_t correction=1.;
  Double_t theoryRatio=1.;
  Double_t effRatio=1.; 
  Double_t correctionExtremeA=1., correctionExtremeB=1.;
  Double_t theoryRatioExtremeA=1., theoryRatioExtremeB=1.;
  Double_t correctionConservativeA=1., correctionConservativeB=1.;
  Double_t theoryRatioConservativeA=1., theoryRatioConservativeB=1.;
  //  Double_t correctionUnc=0.;
  Double_t correctionExtremeAUnc=0., correctionExtremeBUnc=0.;
  Double_t correctionConservativeAUnc=0., correctionConservativeBUnc=0.;

  // declare the output histograms
  fhFc = new TH1D("fhFc","fc correction factor",fnPtBins,fPtBinLimits);
  fhFcMax = new TH1D("fhFcMax","max fc correction factor",fnPtBins,fPtBinLimits);
  fhFcMin = new TH1D("fhFcMin","min fc correction factor",fnPtBins,fPtBinLimits);
  if(fPbPbElossHypothesis) {
    fhFcRcb = new TH2D("fhFcRcb","fc correction factor vs Rcb Eloss hypothesis; p_{T} [GeV/c] ; Rcb Eloss hypothesis ; fc correction",fnPtBins,fPtBinLimits,800,0.,4.);
    fnHypothesis = new TNtuple("fnHypothesis"," Feed-down correction vs hypothesis (fc)","pt:Rcb:fc:fcMax:fcMin");
  }
  // two local control histograms
  TH1D *hTheoryRatio = new TH1D("hTheoryRatio","Theoretical B-->D over c-->D (feed-down/direct) ratio",fnPtBins,fPtBinLimits);
  TH1D *hEffRatio = new TH1D("hEffRatio","Efficiency B-->D over c-->D (feed-down/direct) ratio",fnPtBins,fPtBinLimits);
  // and the output TGraphAsymmErrors
  if (fAsymUncertainties) {
    fgFcExtreme = new TGraphAsymmErrors(fnPtBins+1);
    fgFcExtreme->SetNameTitle("fgFcExtreme","fgFcExtreme");
    fgFcConservative = new TGraphAsymmErrors(fnPtBins+1);
    fgFcConservative->SetNameTitle("fgFcConservative","fgFcConservative");
  }

  fhStatUncEffcFD = new TH1D("fhStatUncEffcFD","direct charm stat unc on the feed-down correction",fnPtBins,fPtBinLimits);
  fhStatUncEffbFD = new TH1D("fhStatUncEffbFD","secondary charm stat unc on the feed-down correction",fnPtBins,fPtBinLimits);
  Double_t correctionConservativeAUncStatEffc=0., correctionConservativeBUncStatEffc=0.;
  Double_t correctionConservativeAUncStatEffb=0., correctionConservativeBUncStatEffb=0.;

  //
  // Compute fc
  //
  for (Int_t ibin=1; ibin<=fnPtBins; ibin++) {

    // Variables initialization
    correction=1.; theoryRatio=1.; effRatio=1.; 
    correctionExtremeA=1.; correctionExtremeB=1.;
    theoryRatioExtremeA=1.; theoryRatioExtremeB=1.;
    correctionConservativeA=1.; correctionConservativeB=1.;
    theoryRatioConservativeA=1.; theoryRatioConservativeB=1.;
    //    correctionUnc=0.;
    correctionExtremeAUnc=0.; correctionExtremeBUnc=0.;
    correctionConservativeAUnc=0.; correctionConservativeBUnc=0.;
    correctionConservativeAUncStatEffc=0.; correctionConservativeBUncStatEffc=0.;
    correctionConservativeAUncStatEffb=0.; correctionConservativeBUncStatEffb=0.;

    //  theory_ratio = (N_b/N_c) 
    theoryRatio = (fhDirectMCpt->GetBinContent(ibin)>0. && fhFeedDownMCpt->GetBinContent(ibin)>0.) ? 
      fhFeedDownMCpt->GetBinContent(ibin) / fhDirectMCpt->GetBinContent(ibin) : 1.0 ;

    //
    // Calculate the uncertainty [ considering only the theoretical uncertainties on Nb & Nc for now !!! ]
    //
    // extreme A = direct-max, feed-down-min
    theoryRatioExtremeA = (fhDirectMCptMax->GetBinContent(ibin)>0. && fhFeedDownMCptMin->GetBinContent(ibin)>0.) ? 
      fhFeedDownMCptMin->GetBinContent(ibin) / fhDirectMCptMax->GetBinContent(ibin) : 1.0 ;
    // extreme B = direct-min, feed-down-max
    theoryRatioExtremeB = (fhDirectMCptMin->GetBinContent(ibin)>0. && fhDirectMCptMax->GetBinContent(ibin)>0.) ? 
      fhFeedDownMCptMax->GetBinContent(ibin) / fhDirectMCptMin->GetBinContent(ibin) : 1.0 ;
    // conservative A = direct-max, feed-down-max
    theoryRatioConservativeA = (fhDirectMCptMax->GetBinContent(ibin)>0. && fhFeedDownMCptMin->GetBinContent(ibin)>0.) ? 
      fhFeedDownMCptMax->GetBinContent(ibin) / fhDirectMCptMax->GetBinContent(ibin) : 1.0 ;
    // conservative B = direct-min, feed-down-min
    theoryRatioConservativeB = (fhDirectMCptMin->GetBinContent(ibin)>0. && fhDirectMCptMax->GetBinContent(ibin)>0.) ? 
      fhFeedDownMCptMin->GetBinContent(ibin) / fhDirectMCptMin->GetBinContent(ibin) : 1.0 ;

    //  eff_ratio = (eff_b/eff_c)
    effRatio = (fhDirectEffpt->GetBinContent(ibin) && fhDirectEffpt->GetBinContent(ibin)!=0.) ? 
      fhFeedDownEffpt->GetBinContent(ibin) / fhDirectEffpt->GetBinContent(ibin) : 1.0 ;

    //   fc = 1 / ( 1 + (eff_b/eff_c)*(N_b/N_c) ) 
    if( TMath::Abs(effRatio - 1.0)<0.0001 || TMath::Abs(theoryRatio - 1.0)<0.0001 ) {
      correction = 1.0;
      correctionExtremeA = 1.0;
      correctionExtremeB = 1.0;
      correctionConservativeA = 1.0; 
      correctionConservativeB = 1.0;
    }
    else {
      correction = ( 1. / ( 1 + ( effRatio * theoryRatio ) ) );
      correctionExtremeA = ( 1. / ( 1 + ( effRatio * theoryRatioExtremeA ) ) );
      correctionExtremeB = ( 1. / ( 1 + ( effRatio * theoryRatioExtremeB ) ) );
      correctionConservativeA = ( 1. / ( 1 + ( effRatio * theoryRatioConservativeA ) ) );
      correctionConservativeB = ( 1. / ( 1 + ( effRatio * theoryRatioConservativeB ) ) );
    }


    // fc uncertainty from (eff_b/eff_c) = fc^2 * (N_b/N_c) * delta(eff_b/eff_c)
    //  delta(eff_b/eff_c) is a percentage = effRatio * sqrt( fGlobalEfficiencyUncertainties[1]^2 + unc_eff_c ^2 + unc_eff_b ^2 ) 
    Double_t relEffUnc =  TMath::Sqrt( fGlobalEfficiencyUncertainties[1]*fGlobalEfficiencyUncertainties[1] + 
				       (fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))*(fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin)) +
				       (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin))*(fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin)) 
				       );

    //    correctionUnc = correction*correction * theoryRatio * effRatio * relEffUnc;
    correctionExtremeAUnc = correctionExtremeA*correctionExtremeA * theoryRatioExtremeA  * effRatio * relEffUnc;
    correctionExtremeBUnc = correctionExtremeB*correctionExtremeB * theoryRatioExtremeB  * effRatio * relEffUnc;

    correctionConservativeAUnc = correctionConservativeA*correctionConservativeA * theoryRatioConservativeA  *effRatio * relEffUnc;
    //
    correctionConservativeAUncStatEffc = correctionConservativeA*correctionConservativeA * theoryRatioConservativeA  *effRatio * 
      (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin));
    correctionConservativeAUncStatEffb = correctionConservativeA*correctionConservativeA * theoryRatioConservativeA  *effRatio * 
      (fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin));

    correctionConservativeBUnc = correctionConservativeB*correctionConservativeB * theoryRatioConservativeB  *effRatio * relEffUnc;

    correctionConservativeBUncStatEffb = correctionConservativeB*correctionConservativeB * theoryRatioConservativeB  *effRatio * 
      (fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin));
    correctionConservativeBUncStatEffc = correctionConservativeB*correctionConservativeB * theoryRatioConservativeB  *effRatio * 
      (fhDirectEffpt->GetBinError(ibin)/fhDirectEffpt->GetBinContent(ibin));


    // Fill in the histograms
    hTheoryRatio->SetBinContent(ibin,theoryRatio);
    hEffRatio->SetBinContent(ibin,effRatio);
    fhFc->SetBinContent(ibin,correction);
    //
    // Estimate how the result varies vs charm/beauty Eloss hypothesis
    //
    if ( correction>1.0e-16 && fPbPbElossHypothesis){
      // Loop over the Eloss hypothesis
      //      Int_t rbin=0;
      for (Float_t rval=0.0025; rval<4.0; rval+=0.005){
	// Central fc with Eloss hypothesis
	Double_t correctionRcb = ( 1. / ( 1 + ( effRatio * theoryRatio * (1/rval) ) ) );
	fhFcRcb->Fill( fhFc->GetBinCenter(ibin) , rval, correctionRcb );
	// 	if(ibin==3){
	// 	  cout << " pt "<< fhFc->GetBinCenter(ibin) <<" bin "<< ibin<<" rval="<<rval<<", rbin="<<rbin<<", fc-Rcb-value="<<correctionRcb<<endl;
	// 	  rbin++;
	// 	}
	// Upper / lower fc with up / low FONLL bands and Eloss hypothesis
	Double_t correctionConservativeARcb = ( 1. / ( 1 + ( effRatio * theoryRatioConservativeA * (1/rval) ) ) );
	Double_t correctionConservativeBRcb = ( 1. / ( 1 + ( effRatio * theoryRatioConservativeB * (1/rval) ) ) );
	Double_t correctionConservativeARcbUnc = correctionConservativeARcb*correctionConservativeARcb * theoryRatioConservativeA * (1/rval) *effRatio * relEffUnc;
	Double_t correctionConservativeBRcbUnc = correctionConservativeBRcb*correctionConservativeBRcb * theoryRatioConservativeB * (1/rval) *effRatio * relEffUnc;
	Double_t consvalRcb[4] = { correctionConservativeARcb - correctionConservativeARcbUnc, correctionConservativeARcb + correctionConservativeARcbUnc,
				   correctionConservativeBRcb - correctionConservativeBRcbUnc, correctionConservativeBRcb + correctionConservativeBRcbUnc};
	Double_t uncConservativeRcbMin = correctionRcb - TMath::MinElement(4,consvalRcb);
	Double_t uncConservativeRcbMax =  TMath::MaxElement(4,consvalRcb) - correctionRcb;
// 	if(ibin==3)
// 	  cout << " pt="<<fhDirectEffpt->GetBinCenter(ibin)<<", hypo="<<rval<<", fc="<<correctionRcb<<" +"<<uncConservativeRcbMax<<" -"<<uncConservativeRcbMin<<endl;
	fnHypothesis->Fill( fhDirectEffpt->GetBinCenter(ibin), rval, correctionRcb, correctionRcb+uncConservativeRcbMax, correctionRcb-uncConservativeRcbMin);
      }
    }
    //
    // Fill the rest of (asymmetric) histograms
    //
    if (fAsymUncertainties) {
      Double_t x = fhDirectMCpt->GetBinCenter(ibin);
      Double_t val[4] = { correctionExtremeA + correctionExtremeAUnc, correctionExtremeA - correctionExtremeAUnc, 
			  correctionExtremeB + correctionExtremeBUnc, correctionExtremeB - correctionExtremeBUnc };
      Double_t uncExtremeMin = correction - TMath::MinElement(4,val);
      Double_t uncExtremeMax = TMath::MaxElement(4,val) - correction;
      fgFcExtreme->SetPoint(ibin,x,correction); // i,x,y
      fgFcExtreme->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),uncExtremeMin,uncExtremeMax); // i,xl,xh,yl,yh
      fhFcMax->SetBinContent(ibin,correction+uncExtremeMax);
      fhFcMin->SetBinContent(ibin,correction-uncExtremeMin);
      Double_t consval[4] = { correctionConservativeA - correctionConservativeAUnc, correctionConservativeA + correctionConservativeAUnc,
			      correctionConservativeB - correctionConservativeBUnc, correctionConservativeB + correctionConservativeBUnc};
      Double_t uncConservativeMin = correction - TMath::MinElement(4,consval);
      Double_t uncConservativeMax =  TMath::MaxElement(4,consval) - correction;
      fgFcConservative->SetPoint(ibin,x,correction); // i,x,y
      fgFcConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),uncConservativeMin,uncConservativeMax); // i,xl,xh,yl,yh
      if( !(correction>0.) ){
	fgFcExtreme->SetPoint(ibin,x,0.); // i,x,y
	fgFcExtreme->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),0.,0.); // i,xl,xh,yl,yh
	fgFcConservative->SetPoint(ibin,x,0.); // i,x,y
	fgFcConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),0.,0.); // i,xl,xh,yl,yh
      }

      Double_t valStatEffc[2] = { correctionConservativeAUncStatEffc/correctionConservativeA, 
				  correctionConservativeBUncStatEffc/correctionConservativeB };
      Double_t valStatEffb[2] = { correctionConservativeAUncStatEffb/correctionConservativeA, 
				  correctionConservativeBUncStatEffb/correctionConservativeB };
      Double_t uncConservativeStatEffc = TMath::MaxElement(2,valStatEffc);
      Double_t uncConservativeStatEffb = TMath::MaxElement(2,valStatEffb);
      fhStatUncEffcFD->SetBinContent(ibin,0.); fhStatUncEffcFD->SetBinError(ibin,uncConservativeStatEffc*100.);
      fhStatUncEffbFD->SetBinContent(ibin,0.); fhStatUncEffbFD->SetBinError(ibin,uncConservativeStatEffb*100.);
      //      cout << " pt "<< fhStatUncEffcFD->GetBinCenter(ibin) <<" bin "<< ibin<<" fc-stat-c ="<<uncConservativeStatEffc<<" fc-stat-b ="<<uncConservativeStatEffb<<endl;
    }

  }

}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::CalculateFeedDownCorrectedSpectrumFc(){
  //
  // Compute the feed-down corrected spectrum if feed-down correction is done via fc factor (bin by bin)
  //    physics = reco * fc / bin-width
  //
  //    uncertainty:             (stat) delta_physics = physics * sqrt ( (delta_reco/reco)^2 )
  //               (syst but feed-down) delta_physics = physics * sqrt ( (delta_reco_syst/reco)^2 )
  //                   (feed-down syst) delta_physics = physics * sqrt ( (delta_fc/fc)^2 )
  //
  //    ( Calculation done bin by bin )
  //
  //  In addition, in HIC the feed-down correction varies with an energy loss hypothesis: Raa(c-->D) / Raa(b-->D) = Rcb

  AliInfo(" Calculating the feed-down corrected spectrum (fc method)");

  if (!fhFc || !fhRECpt) {
    AliError(" Reconstructed or fc distributions are not defined");
    return;
  }

  Double_t value = 0., errvalue = 0., errvalueMax= 0., errvalueMin= 0.;
  Double_t valueExtremeMax= 0., valueExtremeMin= 0.;
  Double_t valueConservativeMax= 0., valueConservativeMin= 0.;
  
  // declare the output histograms
  fhYieldCorr = new TH1D("fhYieldCorr","corrected yield (by fc)",fnPtBins,fPtBinLimits);
  fhYieldCorrMax = new TH1D("fhYieldCorrMax","max corrected yield (by fc)",fnPtBins,fPtBinLimits);
  fhYieldCorrMin = new TH1D("fhYieldCorrMin","min corrected yield (by fc)",fnPtBins,fPtBinLimits);
  if(fPbPbElossHypothesis) fhYieldCorrRcb = new TH2D("fhYieldCorrRcb","corrected yield (by fc) vs Rcb Eloss hypothesis; p_{T} [GeV/c] ; Rcb Eloss hypothesis ; corrected yield",fnPtBins,fPtBinLimits,800,0.,4.);
  // and the output TGraphAsymmErrors
  fgYieldCorr = new TGraphAsymmErrors(fnPtBins+1);
  if (fAsymUncertainties){
    fgYieldCorrExtreme = new TGraphAsymmErrors(fnPtBins+1);
    fgYieldCorrConservative = new TGraphAsymmErrors(fnPtBins+1);
  }
  
  //
  // Do the calculation
  // 
  for (Int_t ibin=1; ibin<=fnPtBins; ibin++) {

    // Variables initialization
    value = 0.; errvalue = 0.; errvalueMax= 0.; errvalueMin= 0.;
    valueExtremeMax= 0.; valueExtremeMin= 0.;
    valueConservativeMax= 0.; valueConservativeMin= 0.;


    // calculate the value 
    //    physics = reco * fc / bin-width
    value = (fhRECpt->GetBinContent(ibin) && fhFc->GetBinContent(ibin)) ? 
      fhRECpt->GetBinContent(ibin) * fhFc->GetBinContent(ibin) : 0. ;
    value /= fhRECpt->GetBinWidth(ibin) ;
    
    // Statistical uncertainty 
    //    (stat) delta_physics = physics * sqrt ( (delta_reco/reco)^2 )
    errvalue = (value!=0. && fhRECpt->GetBinContent(ibin) && fhRECpt->GetBinContent(ibin)!=0.) ?
      value * (fhRECpt->GetBinError(ibin)/fhRECpt->GetBinContent(ibin)) : 0. ; 

    // Calculate the systematic uncertainties
    //    (syst but feed-down) delta_physics = physics * sqrt ( (delta_reco_syst/reco)^2 )
    //        (feed-down syst) delta_physics = physics * sqrt ( (delta_fc/fc)^2 )
    //
    // Protect against null denominator. If so, define uncertainty as null
    if (fhRECpt->GetBinContent(ibin) && fhRECpt->GetBinContent(ibin)!=0.) {

      if (fAsymUncertainties) {

	// Systematics but feed-down
	if (fgRECSystematics) { 
	  errvalueMax = value * ( fgRECSystematics->GetErrorYhigh(ibin) / fhRECpt->GetBinContent(ibin) );
	  errvalueMin = value * ( fgRECSystematics->GetErrorYlow(ibin) / fhRECpt->GetBinContent(ibin) );
	}
	else { errvalueMax = 0.; errvalueMin = 0.; }
	
	// Extreme feed-down systematics
	valueExtremeMax = fhRECpt->GetBinContent(ibin) * ( fhFc->GetBinContent(ibin) + fgFcExtreme->GetErrorYhigh(ibin) ) / fhRECpt->GetBinWidth(ibin) ;
	valueExtremeMin = fhRECpt->GetBinContent(ibin) * ( fhFc->GetBinContent(ibin) - fgFcExtreme->GetErrorYlow(ibin) ) / fhRECpt->GetBinWidth(ibin) ;
	
	// Conservative feed-down systematics
	valueConservativeMax = fhRECpt->GetBinContent(ibin) * ( fhFc->GetBinContent(ibin) + fgFcConservative->GetErrorYhigh(ibin) ) / fhRECpt->GetBinWidth(ibin) ;
	valueConservativeMin = fhRECpt->GetBinContent(ibin) * ( fhFc->GetBinContent(ibin) - fgFcConservative->GetErrorYlow(ibin) ) / fhRECpt->GetBinWidth(ibin) ;

      }

    }
    else { errvalueMax = 0.; errvalueMin = 0.; }
    
    //
    // Fill in the histograms
    //
    fhYieldCorr->SetBinContent(ibin,value);
    fhYieldCorr->SetBinError(ibin,errvalue);
    //
    // Fill the histos and ntuple vs the Eloss hypothesis
    //
    if (fPbPbElossHypothesis) {
      // Loop over the Eloss hypothesis
      for (Float_t rval=0.0025; rval<4.0; rval+=0.005){
	Int_t rbin = FindTH2YBin(fhYieldCorrRcb,rval);
	Double_t fcRcbvalue = fhFcRcb->GetBinContent(ibin,rbin);
	//    physics = reco * fcRcb / bin-width
	Double_t Rcbvalue = (fhRECpt->GetBinContent(ibin) && fcRcbvalue) ? 
	  fhRECpt->GetBinContent(ibin) * fcRcbvalue : 0. ;
	Rcbvalue /= fhRECpt->GetBinWidth(ibin) ;
	fhYieldCorrRcb->Fill( fhYieldCorr->GetBinCenter(ibin) , rval, Rcbvalue );
	// 	  cout << " pt "<< fhRECpt->GetBinCenter(ibin) <<" bin "<< ibin<<" rval="<<rval<<", rbin="<<rbin<<" fc-fcRbvalue="<<fcRcbvalue<<", yield="<<Rcbvalue<<endl;
      }
    }
    if (fAsymUncertainties) {
      Double_t center = fhYieldCorr->GetBinCenter(ibin);
      fgYieldCorr->SetPoint(ibin,center,value); // i,x,y
      fgYieldCorr->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueMin,errvalueMax); // i,xl,xh,yl,yh
      fhYieldCorrMax->SetBinContent(ibin,value+errvalueMax); 
      fhYieldCorrMin->SetBinContent(ibin,value-errvalueMin);
      fgYieldCorrExtreme->SetPoint(ibin,center,value);
      fgYieldCorrExtreme->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),value-valueExtremeMin,valueExtremeMax-value);
      fgYieldCorrConservative->SetPoint(ibin,center,value);
      fgYieldCorrConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),value-valueConservativeMin,valueConservativeMax-value);
    }

  }
  
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::CalculateFeedDownCorrectedSpectrumNb(Double_t deltaY, Double_t branchingRatioBintoFinalDecay) {
  //
  // Compute the feed-down corrected spectrum if feed-down correction is done via Nb (bin by bin)
  //    physics =  [ reco  - (lumi * delta_y * BR_b * eff_trig * eff_b * Nb_th) ] / bin-width
  //
  //    uncertainty:   (stat)  delta_physics = sqrt ( (delta_reco)^2 )  / bin-width
  //     (syst but feed-down)  delta_physics = sqrt ( (delta_reco_syst)^2 )  / bin-width
  //         (feed-down syst)  delta_physics = sqrt ( (k*delta_lumi/lumi)^2 + (k*delta_eff_trig/eff_trig)^2 
  //                                                   + (k*delta_Nb/Nb)^2 + (k*delta_eff/eff)^2  + (k*global_eff_ratio)^2 ) / bin-width
  //                    where k = lumi * delta_y * BR_b * eff_trig * eff_b * Nb_th
  //
  //  In addition, in HIC the feed-down correction varies with an energy loss hypothesis: Raa(b-->D) = Rb
  //    physics =  [ reco  - ( Tab * Nevt * delta_y * BR_b * eff_trig * eff_b * Nb_th * Rb ) ] / bin-width
  //
  AliInfo("Calculating the feed-down correction factor and spectrum (Nb method)");

  Double_t value = 0., errvalue = 0., errvalueMax = 0., errvalueMin = 0., kfactor = 0.;
  Double_t errvalueExtremeMax = 0., errvalueExtremeMin = 0.;
  
  // declare the output histograms
  fhYieldCorr = new TH1D("fhYieldCorr","corrected yield (by Nb)",fnPtBins,fPtBinLimits);
  fhYieldCorrMax = new TH1D("fhYieldCorrMax","max corrected yield (by Nb)",fnPtBins,fPtBinLimits);
  fhYieldCorrMin = new TH1D("fhYieldCorrMin","min corrected yield (by Nb)",fnPtBins,fPtBinLimits);
  if(fPbPbElossHypothesis) {
    fhFcRcb = new TH2D("fhFcRcb","fc correction factor (Nb method) vs Rb Eloss hypothesis; p_{T} [GeV/c] ; Rb Eloss hypothesis ; fc correction",fnPtBins,fPtBinLimits,800,0.,4.);
    fhYieldCorrRcb = new TH2D("fhYieldCorrRcb","corrected yield (by Nb) vs Rb Eloss hypothesis; p_{T} [GeV/c] ; Rb Eloss hypothesis ; corrected yield",fnPtBins,fPtBinLimits,800,0.,4.);
    fnHypothesis = new TNtuple("fnHypothesis"," Feed-down correction vs hypothesis (Nb)","pt:Rb:fc:fcMax:fcMin");
  }
  // and the output TGraphAsymmErrors
  fgYieldCorr = new TGraphAsymmErrors(fnPtBins+1);
  if (fAsymUncertainties){
    fgYieldCorrExtreme = new TGraphAsymmErrors(fnPtBins+1);
    fgYieldCorrConservative = new TGraphAsymmErrors(fnPtBins+1);
    // Define fc-conservative 
    fgFcConservative = new TGraphAsymmErrors(fnPtBins+1);
    AliInfo(" Beware the conservative & extreme uncertainties are equal by definition !");
  }

  // variables to define fc-conservative 
  double correction=0, correctionMax=0., correctionMin=0.;

  fhStatUncEffcFD = new TH1D("fhStatUncEffcFD","direct charm stat unc on the feed-down correction",fnPtBins,fPtBinLimits);
  fhStatUncEffbFD = new TH1D("fhStatUncEffbFD","secondary charm stat unc on the feed-down correction",fnPtBins,fPtBinLimits);
  Double_t correctionUncStatEffc=0.;
  Double_t correctionUncStatEffb=0.;


  //
  // Do the calculation
  // 
  for (Int_t ibin=1; ibin<=fnPtBins; ibin++) {

    // Calculate the value
    //    physics =  [ reco  - (lumi * delta_y * BR_b * eff_trig * eff_b * Nb_th) ] / bin-width
    //  In HIC :   physics =  [ reco  - ( Tab * Nevt * delta_y * BR_b * eff_trig * eff_b * Nb_th * Rb ) ] / bin-width
    //
    //
    Double_t frac = 1.0, errfrac =0.;

    // Variables initialization
    value = 0.; errvalue = 0.; errvalueMax = 0.; errvalueMin = 0.; kfactor = 0.;
    errvalueExtremeMax = 0.; errvalueExtremeMin = 0.;
    correction=0; correctionMax=0.; correctionMin=0.;
    correctionUncStatEffc=0.; correctionUncStatEffb=0.;

    if(fPbPbElossHypothesis) {
      frac = fTab[0]*fNevts;
      if(fIsEventPlane) frac = frac/2.0;
      errfrac = frac * TMath::Sqrt( (fTab[1]/fTab[0])*(fTab[1]/fTab[0]) + (1/fNevts) );
    } else {
      frac = fLuminosity[0]; 
      errfrac = fLuminosity[1];
    }
    
    value = ( fhRECpt->GetBinContent(ibin)>0. && fhRECpt->GetBinContent(ibin)!=0. && 
	      fhFeedDownMCpt->GetBinContent(ibin)>0. && fhFeedDownEffpt->GetBinContent(ibin)>0. ) ?
      fhRECpt->GetBinContent(ibin) - frac*(deltaY*branchingRatioBintoFinalDecay*fParticleAntiParticle*fTrigEfficiency[0]*fhFeedDownEffpt->GetBinContent(ibin)*fhFeedDownMCpt->GetBinContent(ibin) * fhRECpt->GetBinWidth(ibin) )
      : 0. ;
    value /= fhRECpt->GetBinWidth(ibin);
    if (value<0.) value =0.;

    //  Statistical uncertainty:   delta_physics = sqrt ( (delta_reco)^2 )  / bin-width
    errvalue = (value!=0. && fhRECpt->GetBinError(ibin) && fhRECpt->GetBinError(ibin)!=0.)  ? 
      fhRECpt->GetBinError(ibin) : 0.;
    errvalue /= fhRECpt->GetBinWidth(ibin);

    // Correction (fc) : Estimate of the relative amount feed-down subtracted
    // correction =  [ 1  - (lumi * delta_y * BR_b * eff_trig * eff_b * Nb_th)/reco ] 
    // in HIC: correction =  [ 1  - ( Tab * Nevt * delta_y * BR_b * eff_trig * eff_b * Nb_th)/reco ]
    correction = (value>0.) ? 
      1 - (frac*deltaY*branchingRatioBintoFinalDecay*fParticleAntiParticle*fTrigEfficiency[0]*fhFeedDownEffpt->GetBinContent(ibin)*fhFeedDownMCpt->GetBinContent(ibin) * fhRECpt->GetBinWidth(ibin) ) / fhRECpt->GetBinContent(ibin)  : 0. ;
    if (correction<0.) correction = 0.;

    //    cout << " pt="<< fhRECpt->GetBinCenter(ibin) << " rec="<< fhRECpt->GetBinContent(ibin) << ", corr="<<correction<<" = 1 - "<<  (frac*deltaY*branchingRatioBintoFinalDecay*fParticleAntiParticle*fTrigEfficiency[0]*fhFeedDownEffpt->GetBinContent(ibin)*fhFeedDownMCpt->GetBinContent(ibin) * fhRECpt->GetBinWidth(ibin) ) / fhRECpt->GetBinContent(ibin) << endl;

    // Systematic uncertainties
    //     (syst but feed-down)  delta_physics = sqrt ( (delta_reco_syst)^2 )  / bin-width
    //         (feed-down syst)  delta_physics = sqrt ( (k*delta_lumi/lumi)^2 + (k*delta_eff_trig/eff_trig)^2 
    //                                                   + (k*delta_Nb/Nb)^2 + (k*delta_eff/eff)^2 + (k*global_eff_ratio)^2 ) / bin-width
    //                    where k = lumi * delta_y * BR_b * eff_trig * eff_b * Nb_th * bin-width
    kfactor = frac*deltaY*branchingRatioBintoFinalDecay*fParticleAntiParticle*fTrigEfficiency[0]*fhFeedDownEffpt->GetBinContent(ibin)*fhFeedDownMCpt->GetBinContent(ibin) * fhRECpt->GetBinWidth(ibin) ;
    //
    if (fAsymUncertainties && value>0.) {
      Double_t nb =  fhFeedDownMCpt->GetBinContent(ibin);
      Double_t nbDmax = fhFeedDownMCptMax->GetBinContent(ibin) - fhFeedDownMCpt->GetBinContent(ibin);
      Double_t nbDmin = fhFeedDownMCpt->GetBinContent(ibin) - fhFeedDownMCptMin->GetBinContent(ibin);

      // Systematics but feed-down
      if (fgRECSystematics){
	errvalueMax = fgRECSystematics->GetErrorYhigh(ibin) / fhRECpt->GetBinWidth(ibin) ;
	errvalueMin = fgRECSystematics->GetErrorYlow(ibin) / fhRECpt->GetBinWidth(ibin);
      }
      else { errvalueMax = 0.; errvalueMin = 0.; }
  
      // Feed-down systematics
      // min value with the maximum Nb
      Double_t errCom =  ( (kfactor*errfrac/frac)*(kfactor*errfrac/frac) ) +
	( (kfactor*fTrigEfficiency[1]/fTrigEfficiency[0])*(kfactor*fTrigEfficiency[1]/fTrigEfficiency[0]) ) +
	( (kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))*(kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin)) ) +
	( (kfactor*fGlobalEfficiencyUncertainties[1])*(kfactor*fGlobalEfficiencyUncertainties[1]) ) ;
      errvalueExtremeMin = TMath::Sqrt( errCom + ( (kfactor*nbDmax/nb)*(kfactor*nbDmax/nb) ) ) / fhRECpt->GetBinWidth(ibin);
      // max value with the minimum Nb
      errvalueExtremeMax =  TMath::Sqrt( errCom + ( (kfactor*nbDmin/nb)*(kfactor*nbDmin/nb) ) ) / fhRECpt->GetBinWidth(ibin);

      // Correction systematics (fc)
      // min value with the maximum Nb
      correctionMin = TMath::Sqrt( errCom + ( (kfactor*nbDmax/nb)*(kfactor*nbDmax/nb) ) ) / fhRECpt->GetBinContent(ibin) ;
      // max value with the minimum Nb
      correctionMax =  TMath::Sqrt( errCom + ( (kfactor*nbDmin/nb)*(kfactor*nbDmin/nb) ) ) / fhRECpt->GetBinContent(ibin) ;
      //
      correctionUncStatEffb = TMath::Sqrt(  ( (kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))*(kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))	)
					    ) / fhRECpt->GetBinContent(ibin) ;
      correctionUncStatEffc = 0.;
    }
    else{ // Don't consider Nb uncertainty in this case [ to be tested!!! ]
      errvalueExtremeMax =  TMath::Sqrt( ( (kfactor*errfrac/frac)*(kfactor*errfrac/frac) ) +
					 ( (kfactor*fTrigEfficiency[1]/fTrigEfficiency[0])*(kfactor*fTrigEfficiency[1]/fTrigEfficiency[0]) )  +
					 ( (kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))*(kfactor*fhFeedDownEffpt->GetBinError(ibin)/fhFeedDownEffpt->GetBinContent(ibin))	)  +
					 ( (kfactor*fGlobalEfficiencyUncertainties[1])*(kfactor*fGlobalEfficiencyUncertainties[1]) )
					 ) / fhRECpt->GetBinWidth(ibin);
      errvalueExtremeMin =  errvalueExtremeMax ;
    }
    

    // fill in histograms
    fhYieldCorr->SetBinContent(ibin,value);
    fhYieldCorr->SetBinError(ibin,errvalue);    
    //
    // Estimate how the result varies vs charm/beauty Eloss hypothesis
    //
    if ( correction>1.0e-16 && fPbPbElossHypothesis){
      // Loop over the Eloss hypothesis
      //      Int_t rbin=0;
      for (Float_t rval=0.0025; rval<4.0; rval+=0.005){
	// correction =  [ 1  - (Tab *Nevt * delta_y * BR_b * eff_trig * eff_b * Nb_th *binwidth )* (rval) /reco ] 
	Double_t fcRcbvalue = 1 - (frac*deltaY*branchingRatioBintoFinalDecay*fParticleAntiParticle*fTrigEfficiency[0]*fhFeedDownEffpt->GetBinContent(ibin)*fhFeedDownMCpt->GetBinContent(ibin)*fhRECpt->GetBinWidth(ibin) * rval ) / fhRECpt->GetBinContent(ibin) ;
	//	cout << " rval="<<rval <<" , fc="<<fcRcbvalue<<"    ";
	if(fcRcbvalue<1.e-16) fcRcbvalue=0.;
	fhFcRcb->Fill( fhRECpt->GetBinCenter(ibin) , rval, fcRcbvalue );
	//    physics = reco * fcRcb / bin-width
	Double_t Rcbvalue = (fhRECpt->GetBinContent(ibin) && fcRcbvalue) ? 
	  fhRECpt->GetBinContent(ibin) * fcRcbvalue : 0. ;
	Rcbvalue /= fhRECpt->GetBinWidth(ibin) ;
	fhYieldCorrRcb->Fill( fhYieldCorr->GetBinCenter(ibin) , rval, Rcbvalue );
	// 	if(ibin==3){
	// 	  cout << " pt "<< fhFcRcb->GetBinCenter(ibin) <<" bin "<< ibin<<" rval="<<rval<<", rbin="<<rbin<<", fc-Rb-value="<< fcRcbvalue << ", yield-Rb-value="<< Rcbvalue <<endl;
	//	cout << " pt "<< fhFcRcb->GetBinCenter(ibin) <<" bin "<< ibin<<" rval="<<rval<<", fc-Rb-value="<< fcRcbvalue << ", yield-Rb-value="<< Rcbvalue <<endl;
	// 	  rbin++;
	// 	}
	Double_t correctionMaxRcb = correctionMax*rval;
	Double_t correctionMinRcb = correctionMin*rval;
	fnHypothesis->Fill( fhYieldCorr->GetBinCenter(ibin), rval, fcRcbvalue, fcRcbvalue + correctionMaxRcb, fcRcbvalue - correctionMinRcb);
// 	if(ibin==3){
// 	  cout << " pt="<< fhFcRcb->GetBinCenter(ibin) <<", rval="<<rval<<", fc="<<fcRcbvalue<<" +"<<correctionMaxRcb<<" -"<<correctionMinRcb<<endl;}
      }
    }
    //
    // Fill the rest of (asymmetric) histograms
    //
    if (fAsymUncertainties) {
      Double_t x = fhYieldCorr->GetBinCenter(ibin);
      fgYieldCorr->SetPoint(ibin,x,value); // i,x,y
      fgYieldCorr->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueMin,errvalueMax); // i,xl,xh,yl,yh
      fhYieldCorrMax->SetBinContent(ibin,value+errvalueMax);
      fhYieldCorrMin->SetBinContent(ibin,value-errvalueMin);
      fgYieldCorrExtreme->SetPoint(ibin,x,value); // i,x,y
      fgYieldCorrExtreme->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueExtremeMin,errvalueExtremeMax); // i,xl,xh,yl,yh
      fgYieldCorrConservative->SetPoint(ibin,x,value); // i,x,y
      fgYieldCorrConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),errvalueExtremeMin,errvalueExtremeMax); // i,xl,xh,yl,yh
      //      cout << " bin " << ibin << ", correction " << correction << ", min correction unc " << correctionMin << ", max correction unc " << correctionMax << endl;
      if(correction>0.){
	fgFcConservative->SetPoint(ibin,x,correction);
	fgFcConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),correctionMin,correctionMax);
	
	fhStatUncEffbFD->SetBinContent(ibin,0.); fhStatUncEffbFD->SetBinError(ibin,correctionUncStatEffb/correction*100.);
	fhStatUncEffcFD->SetBinContent(ibin,0.); fhStatUncEffcFD->SetBinError(ibin,correctionUncStatEffc/correction*100.);
	//	cout << " pt "<< fhStatUncEffcFD->GetBinCenter(ibin) <<" bin "<< ibin<<" fc-stat-c ="<< correctionUncStatEffc/correction <<" fc-stat-b ="<< correctionUncStatEffb/correction <<endl;
      }
      else{
	fgFcConservative->SetPoint(ibin,x,0.);
	fgFcConservative->SetPointError(ibin,(fPtBinWidths[ibin-1]/2.),(fPtBinWidths[ibin-1]/2.),0.,0.);
      }
    }

  }
  
}


//_________________________________________________________________________________________________________
void AliHFPtSpectrum::ComputeSystUncertainties(AliHFSystErr *systematics, Bool_t combineFeedDown) {
  //
  // Function that re-calculates the global systematic uncertainties
  //   by calling the class AliHFSystErr and combining those
  //   (in quadrature) with the feed-down subtraction uncertainties
  //

  // Estimate the feed-down uncertainty in percentage
  Int_t nentries = 0;
  TGraphAsymmErrors *grErrFeeddown = 0;
  Double_t x=0., y=0., errx=0., erryl=0., erryh=0;
  if(fFeedDownOption!=0) {
    nentries = fgSigmaCorrConservative->GetN();
    grErrFeeddown = new TGraphAsymmErrors(nentries);
    for(Int_t i=0; i<nentries; i++) {
      x=0.; y=0.; errx=0.; erryl=0.; erryh=0.;
      fgSigmaCorrConservative->GetPoint(i,x,y);
      if(y>0.){
	errx = fgSigmaCorrConservative->GetErrorXlow(i) ;
	erryl = fgSigmaCorrConservative->GetErrorYlow(i) / y ;
	erryh = fgSigmaCorrConservative->GetErrorYhigh(i) / y ;
      }
      //    cout << " x "<< x << " +- "<<errx<<" , y "<<y<<" + "<<erryh<<" - "<<erryl<<endl; 
      grErrFeeddown->SetPoint(i,x,0.);
      grErrFeeddown->SetPointError(i,errx,errx,erryl,erryh); //i, xl, xh, yl, yh
    }
  }

  // Draw all the systematics independently
  systematics->DrawErrors(grErrFeeddown);

  // Set the sigma systematic uncertainties
  // possibly combine with the feed-down uncertainties 
  Double_t errylcomb=0., erryhcomb=0;
  for(Int_t i=1; i<nentries; i++) {
    fgSigmaCorr->GetPoint(i,x,y);
    errx = grErrFeeddown->GetErrorXlow(i) ;
    erryl = grErrFeeddown->GetErrorYlow(i);
    erryh = grErrFeeddown->GetErrorYhigh(i);
    if (combineFeedDown) {
      errylcomb = systematics->GetTotalSystErr(x,erryl) * y ;
      erryhcomb = systematics->GetTotalSystErr(x,erryh) * y ;
    } else {
      errylcomb = systematics->GetTotalSystErr(x) * y ;
      erryhcomb = systematics->GetTotalSystErr(x) * y ;
    }
    fgSigmaCorr->SetPointError(i,errx,errx,errylcomb,erryhcomb);
    //
    fhSigmaCorrDataSyst->SetBinContent(i,y);
    erryl = systematics->GetTotalSystErr(x) * y ;
    fhSigmaCorrDataSyst->SetBinError(i,erryl);
  }

}


//_________________________________________________________________________________________________________
void AliHFPtSpectrum::DrawSpectrum(TGraphAsymmErrors *gPrediction) {
  //
  // Example method to draw the corrected spectrum & the theoretical prediction
  //

  TCanvas *csigma = new TCanvas("csigma","Draw the corrected cross-section & the prediction");
  csigma->SetFillColor(0);
  gPrediction->GetXaxis()->SetTitleSize(0.05);
  gPrediction->GetXaxis()->SetTitleOffset(0.95);
  gPrediction->GetYaxis()->SetTitleSize(0.05);
  gPrediction->GetYaxis()->SetTitleOffset(0.95);
  gPrediction->GetXaxis()->SetTitle("p_{T}  [GeV]");
  gPrediction->GetYaxis()->SetTitle("BR #times #frac{d#sigma}{dp_{T}} |_{|y|<0.5}   [pb/GeV]");
  gPrediction->SetLineColor(kGreen+2);
  gPrediction->SetLineWidth(3);
  gPrediction->SetFillColor(kGreen+1);
  gPrediction->Draw("3CA");
  fgSigmaCorr->SetLineColor(kRed);
  fgSigmaCorr->SetLineWidth(1);
  fgSigmaCorr->SetFillColor(kRed);
  fgSigmaCorr->SetFillStyle(0);
  fgSigmaCorr->Draw("2");
  fhSigmaCorr->SetMarkerColor(kRed);
  fhSigmaCorr->Draw("esame");
  csigma->SetLogy();
  TLegend * leg = new TLegend(0.7,0.75,0.87,0.5);
  leg->SetBorderSize(0);
  leg->SetLineColor(0);
  leg->SetFillColor(0);
  leg->SetTextFont(42);
  leg->AddEntry(gPrediction,"FONLL ","fl");
  leg->AddEntry(fhSigmaCorr,"data stat. unc.","pl");
  leg->AddEntry(fgSigmaCorr,"data syst. unc.","f");
  leg->Draw();
  csigma->Draw();

}

//_________________________________________________________________________________________________________
TH1D * AliHFPtSpectrum::ReweightHisto(TH1D *hToReweight, TH1D *hReference){
  //
  // Function to  reweight histograms for testing purposes: 
  // This function takes the histo hToReweight and reweights 
  //  it (its pt shape) with respect to hReference 
  // 

  // check histograms consistency
  Bool_t areconsistent=kTRUE;
  areconsistent &= CheckHistosConsistency(hToReweight,hReference);
  if (!areconsistent) {
    AliInfo("the histograms to reweight are not consistent (bin width, bounds)"); 
    return NULL;
  }

  // define a new empty histogram
  TH1D *hReweighted = (TH1D*)hToReweight->Clone("hReweighted");
  hReweighted->Reset();
  Double_t weight=1.0;
  Double_t yvalue=1.0; 
  Double_t integralRef = hReference->Integral();
  Double_t integralH = hToReweight->Integral();

  // now reweight the spectra
  //
  // the weight is the relative probability of the given pt bin in the reference histo
  //  divided by its relative probability (to normalize it) on the histo to re-weight
  for (Int_t i=0; i<=hToReweight->GetNbinsX(); i++) {
    weight = (hReference->GetBinContent(i)/integralRef) / (hToReweight->GetBinContent(i)/integralH) ;
    yvalue = hToReweight->GetBinContent(i);
    hReweighted->SetBinContent(i,yvalue*weight);
  }

  return (TH1D*)hReweighted;
}

//_________________________________________________________________________________________________________
TH1D * AliHFPtSpectrum::ReweightRecHisto(TH1D *hRecToReweight, TH1D *hMCToReweight, TH1D *hMCReference){
  //
  // Function to  reweight histograms for testing purposes: 
  // This function takes the histo hToReweight and reweights 
  //  it (its pt shape) with respect to hReference /hMCToReweight
  // 

  // check histograms consistency
  Bool_t areconsistent=kTRUE;
  areconsistent &= CheckHistosConsistency(hMCToReweight,hMCReference);
  areconsistent &= CheckHistosConsistency(hRecToReweight,hMCReference);
  if (!areconsistent) {
    AliInfo("the histograms to reweight are not consistent (bin width, bounds)"); 
    return NULL;
  }

  // define a new empty histogram
  TH1D *hReweighted = (TH1D*)hMCToReweight->Clone("hReweighted");
  hReweighted->Reset();
  TH1D *hRecReweighted = (TH1D*)hRecToReweight->Clone("hRecReweighted");
  hRecReweighted->Reset();
  Double_t weight=1.0;
  Double_t yvalue=1.0, yrecvalue=1.0; 
  Double_t integralRef = hMCReference->Integral();
  Double_t integralH = hMCToReweight->Integral();

  // now reweight the spectra
  //
  // the weight is the relative probability of the given pt bin 
  //  that should be applied in the MC histo to get the reference histo shape
  //  Probabilities are properly normalized.
  for (Int_t i=0; i<=hMCToReweight->GetNbinsX(); i++) {
    weight = (hMCReference->GetBinContent(i)/integralRef) / (hMCToReweight->GetBinContent(i)/integralH) ;
    yvalue = hMCToReweight->GetBinContent(i);
    hReweighted->SetBinContent(i,yvalue*weight);
    yrecvalue = hRecToReweight->GetBinContent(i);
    hRecReweighted->SetBinContent(i,yrecvalue*weight);
  }

  return (TH1D*)hRecReweighted;
}



//_________________________________________________________________________________________________________
Int_t AliHFPtSpectrum::FindTH2YBin(TH2D *histo, Float_t yvalue){
  //
  // Function to find the y-axis bin of a TH2 for a given y-value
  //
  
  Int_t nbins = histo->GetNbinsY();
  Int_t ybin=0;
  for (int j=0; j<=nbins; j++) {
    Float_t value = histo->GetYaxis()->GetBinCenter(j);
    Float_t width = histo->GetYaxis()->GetBinWidth(j);
    //    if( TMath::Abs(yvalue-value)<= width/2. ) {
    if( TMath::Abs(yvalue-value)<= width ) {
      ybin =j;
      //      cout <<" value "<<value << ", yval "<< yvalue<<", bin width "<<width/2.<< " y ="<<ybin<<endl;
      break;
    }
  }
  
  return ybin;
}

//_________________________________________________________________________________________________________
void AliHFPtSpectrum::ResetStatUncEff(){

  Int_t entries = fhDirectEffpt->GetNbinsX();
  for(Int_t i=0; i<=entries; i++){
    fhDirectEffpt->SetBinError(i,0.);
    fhFeedDownEffpt->SetBinError(i,0.);
  }

}
 AliHFPtSpectrum.cxx:1
 AliHFPtSpectrum.cxx:2
 AliHFPtSpectrum.cxx:3
 AliHFPtSpectrum.cxx:4
 AliHFPtSpectrum.cxx:5
 AliHFPtSpectrum.cxx:6
 AliHFPtSpectrum.cxx:7
 AliHFPtSpectrum.cxx:8
 AliHFPtSpectrum.cxx:9
 AliHFPtSpectrum.cxx:10
 AliHFPtSpectrum.cxx:11
 AliHFPtSpectrum.cxx:12
 AliHFPtSpectrum.cxx:13
 AliHFPtSpectrum.cxx:14
 AliHFPtSpectrum.cxx:15
 AliHFPtSpectrum.cxx:16
 AliHFPtSpectrum.cxx:17
 AliHFPtSpectrum.cxx:18
 AliHFPtSpectrum.cxx:19
 AliHFPtSpectrum.cxx:20
 AliHFPtSpectrum.cxx:21
 AliHFPtSpectrum.cxx:22
 AliHFPtSpectrum.cxx:23
 AliHFPtSpectrum.cxx:24
 AliHFPtSpectrum.cxx:25
 AliHFPtSpectrum.cxx:26
 AliHFPtSpectrum.cxx:27
 AliHFPtSpectrum.cxx:28
 AliHFPtSpectrum.cxx:29
 AliHFPtSpectrum.cxx:30
 AliHFPtSpectrum.cxx:31
 AliHFPtSpectrum.cxx:32
 AliHFPtSpectrum.cxx:33
 AliHFPtSpectrum.cxx:34
 AliHFPtSpectrum.cxx:35
 AliHFPtSpectrum.cxx:36
 AliHFPtSpectrum.cxx:37
 AliHFPtSpectrum.cxx:38
 AliHFPtSpectrum.cxx:39
 AliHFPtSpectrum.cxx:40
 AliHFPtSpectrum.cxx:41
 AliHFPtSpectrum.cxx:42
 AliHFPtSpectrum.cxx:43
 AliHFPtSpectrum.cxx:44
 AliHFPtSpectrum.cxx:45
 AliHFPtSpectrum.cxx:46
 AliHFPtSpectrum.cxx:47
 AliHFPtSpectrum.cxx:48
 AliHFPtSpectrum.cxx:49
 AliHFPtSpectrum.cxx:50
 AliHFPtSpectrum.cxx:51
 AliHFPtSpectrum.cxx:52
 AliHFPtSpectrum.cxx:53
 AliHFPtSpectrum.cxx:54
 AliHFPtSpectrum.cxx:55
 AliHFPtSpectrum.cxx:56
 AliHFPtSpectrum.cxx:57
 AliHFPtSpectrum.cxx:58
 AliHFPtSpectrum.cxx:59
 AliHFPtSpectrum.cxx:60
 AliHFPtSpectrum.cxx:61
 AliHFPtSpectrum.cxx:62
 AliHFPtSpectrum.cxx:63
 AliHFPtSpectrum.cxx:64
 AliHFPtSpectrum.cxx:65
 AliHFPtSpectrum.cxx:66
 AliHFPtSpectrum.cxx:67
 AliHFPtSpectrum.cxx:68
 AliHFPtSpectrum.cxx:69
 AliHFPtSpectrum.cxx:70
 AliHFPtSpectrum.cxx:71
 AliHFPtSpectrum.cxx:72
 AliHFPtSpectrum.cxx:73
 AliHFPtSpectrum.cxx:74
 AliHFPtSpectrum.cxx:75
 AliHFPtSpectrum.cxx:76
 AliHFPtSpectrum.cxx:77
 AliHFPtSpectrum.cxx:78
 AliHFPtSpectrum.cxx:79
 AliHFPtSpectrum.cxx:80
 AliHFPtSpectrum.cxx:81
 AliHFPtSpectrum.cxx:82
 AliHFPtSpectrum.cxx:83
 AliHFPtSpectrum.cxx:84
 AliHFPtSpectrum.cxx:85
 AliHFPtSpectrum.cxx:86
 AliHFPtSpectrum.cxx:87
 AliHFPtSpectrum.cxx:88
 AliHFPtSpectrum.cxx:89
 AliHFPtSpectrum.cxx:90
 AliHFPtSpectrum.cxx:91
 AliHFPtSpectrum.cxx:92
 AliHFPtSpectrum.cxx:93
 AliHFPtSpectrum.cxx:94
 AliHFPtSpectrum.cxx:95
 AliHFPtSpectrum.cxx:96
 AliHFPtSpectrum.cxx:97
 AliHFPtSpectrum.cxx:98
 AliHFPtSpectrum.cxx:99
 AliHFPtSpectrum.cxx:100
 AliHFPtSpectrum.cxx:101
 AliHFPtSpectrum.cxx:102
 AliHFPtSpectrum.cxx:103
 AliHFPtSpectrum.cxx:104
 AliHFPtSpectrum.cxx:105
 AliHFPtSpectrum.cxx:106
 AliHFPtSpectrum.cxx:107
 AliHFPtSpectrum.cxx:108
 AliHFPtSpectrum.cxx:109
 AliHFPtSpectrum.cxx:110
 AliHFPtSpectrum.cxx:111
 AliHFPtSpectrum.cxx:112
 AliHFPtSpectrum.cxx:113
 AliHFPtSpectrum.cxx:114
 AliHFPtSpectrum.cxx:115
 AliHFPtSpectrum.cxx:116
 AliHFPtSpectrum.cxx:117
 AliHFPtSpectrum.cxx:118
 AliHFPtSpectrum.cxx:119
 AliHFPtSpectrum.cxx:120
 AliHFPtSpectrum.cxx:121
 AliHFPtSpectrum.cxx:122
 AliHFPtSpectrum.cxx:123
 AliHFPtSpectrum.cxx:124
 AliHFPtSpectrum.cxx:125
 AliHFPtSpectrum.cxx:126
 AliHFPtSpectrum.cxx:127
 AliHFPtSpectrum.cxx:128
 AliHFPtSpectrum.cxx:129
 AliHFPtSpectrum.cxx:130
 AliHFPtSpectrum.cxx:131
 AliHFPtSpectrum.cxx:132
 AliHFPtSpectrum.cxx:133
 AliHFPtSpectrum.cxx:134
 AliHFPtSpectrum.cxx:135
 AliHFPtSpectrum.cxx:136
 AliHFPtSpectrum.cxx:137
 AliHFPtSpectrum.cxx:138
 AliHFPtSpectrum.cxx:139
 AliHFPtSpectrum.cxx:140
 AliHFPtSpectrum.cxx:141
 AliHFPtSpectrum.cxx:142
 AliHFPtSpectrum.cxx:143
 AliHFPtSpectrum.cxx:144
 AliHFPtSpectrum.cxx:145
 AliHFPtSpectrum.cxx:146
 AliHFPtSpectrum.cxx:147
 AliHFPtSpectrum.cxx:148
 AliHFPtSpectrum.cxx:149
 AliHFPtSpectrum.cxx:150
 AliHFPtSpectrum.cxx:151
 AliHFPtSpectrum.cxx:152
 AliHFPtSpectrum.cxx:153
 AliHFPtSpectrum.cxx:154
 AliHFPtSpectrum.cxx:155
 AliHFPtSpectrum.cxx:156
 AliHFPtSpectrum.cxx:157
 AliHFPtSpectrum.cxx:158
 AliHFPtSpectrum.cxx:159
 AliHFPtSpectrum.cxx:160
 AliHFPtSpectrum.cxx:161
 AliHFPtSpectrum.cxx:162
 AliHFPtSpectrum.cxx:163
 AliHFPtSpectrum.cxx:164
 AliHFPtSpectrum.cxx:165
 AliHFPtSpectrum.cxx:166
 AliHFPtSpectrum.cxx:167
 AliHFPtSpectrum.cxx:168
 AliHFPtSpectrum.cxx:169
 AliHFPtSpectrum.cxx:170
 AliHFPtSpectrum.cxx:171
 AliHFPtSpectrum.cxx:172
 AliHFPtSpectrum.cxx:173
 AliHFPtSpectrum.cxx:174
 AliHFPtSpectrum.cxx:175
 AliHFPtSpectrum.cxx:176
 AliHFPtSpectrum.cxx:177
 AliHFPtSpectrum.cxx:178
 AliHFPtSpectrum.cxx:179
 AliHFPtSpectrum.cxx:180
 AliHFPtSpectrum.cxx:181
 AliHFPtSpectrum.cxx:182
 AliHFPtSpectrum.cxx:183
 AliHFPtSpectrum.cxx:184
 AliHFPtSpectrum.cxx:185
 AliHFPtSpectrum.cxx:186
 AliHFPtSpectrum.cxx:187
 AliHFPtSpectrum.cxx:188
 AliHFPtSpectrum.cxx:189
 AliHFPtSpectrum.cxx:190
 AliHFPtSpectrum.cxx:191
 AliHFPtSpectrum.cxx:192
 AliHFPtSpectrum.cxx:193
 AliHFPtSpectrum.cxx:194
 AliHFPtSpectrum.cxx:195
 AliHFPtSpectrum.cxx:196
 AliHFPtSpectrum.cxx:197
 AliHFPtSpectrum.cxx:198
 AliHFPtSpectrum.cxx:199
 AliHFPtSpectrum.cxx:200
 AliHFPtSpectrum.cxx:201
 AliHFPtSpectrum.cxx:202
 AliHFPtSpectrum.cxx:203
 AliHFPtSpectrum.cxx:204
 AliHFPtSpectrum.cxx:205
 AliHFPtSpectrum.cxx:206
 AliHFPtSpectrum.cxx:207
 AliHFPtSpectrum.cxx:208
 AliHFPtSpectrum.cxx:209
 AliHFPtSpectrum.cxx:210
 AliHFPtSpectrum.cxx:211
 AliHFPtSpectrum.cxx:212
 AliHFPtSpectrum.cxx:213
 AliHFPtSpectrum.cxx:214
 AliHFPtSpectrum.cxx:215
 AliHFPtSpectrum.cxx:216
 AliHFPtSpectrum.cxx:217
 AliHFPtSpectrum.cxx:218
 AliHFPtSpectrum.cxx:219
 AliHFPtSpectrum.cxx:220
 AliHFPtSpectrum.cxx:221
 AliHFPtSpectrum.cxx:222
 AliHFPtSpectrum.cxx:223
 AliHFPtSpectrum.cxx:224
 AliHFPtSpectrum.cxx:225
 AliHFPtSpectrum.cxx:226
 AliHFPtSpectrum.cxx:227
 AliHFPtSpectrum.cxx:228
 AliHFPtSpectrum.cxx:229
 AliHFPtSpectrum.cxx:230
 AliHFPtSpectrum.cxx:231
 AliHFPtSpectrum.cxx:232
 AliHFPtSpectrum.cxx:233
 AliHFPtSpectrum.cxx:234
 AliHFPtSpectrum.cxx:235
 AliHFPtSpectrum.cxx:236
 AliHFPtSpectrum.cxx:237
 AliHFPtSpectrum.cxx:238
 AliHFPtSpectrum.cxx:239
 AliHFPtSpectrum.cxx:240
 AliHFPtSpectrum.cxx:241
 AliHFPtSpectrum.cxx:242
 AliHFPtSpectrum.cxx:243
 AliHFPtSpectrum.cxx:244
 AliHFPtSpectrum.cxx:245
 AliHFPtSpectrum.cxx:246
 AliHFPtSpectrum.cxx:247
 AliHFPtSpectrum.cxx:248
 AliHFPtSpectrum.cxx:249
 AliHFPtSpectrum.cxx:250
 AliHFPtSpectrum.cxx:251
 AliHFPtSpectrum.cxx:252
 AliHFPtSpectrum.cxx:253
 AliHFPtSpectrum.cxx:254
 AliHFPtSpectrum.cxx:255
 AliHFPtSpectrum.cxx:256
 AliHFPtSpectrum.cxx:257
 AliHFPtSpectrum.cxx:258
 AliHFPtSpectrum.cxx:259
 AliHFPtSpectrum.cxx:260
 AliHFPtSpectrum.cxx:261
 AliHFPtSpectrum.cxx:262
 AliHFPtSpectrum.cxx:263
 AliHFPtSpectrum.cxx:264
 AliHFPtSpectrum.cxx:265
 AliHFPtSpectrum.cxx:266
 AliHFPtSpectrum.cxx:267
 AliHFPtSpectrum.cxx:268
 AliHFPtSpectrum.cxx:269
 AliHFPtSpectrum.cxx:270
 AliHFPtSpectrum.cxx:271
 AliHFPtSpectrum.cxx:272
 AliHFPtSpectrum.cxx:273
 AliHFPtSpectrum.cxx:274
 AliHFPtSpectrum.cxx:275
 AliHFPtSpectrum.cxx:276
 AliHFPtSpectrum.cxx:277
 AliHFPtSpectrum.cxx:278
 AliHFPtSpectrum.cxx:279
 AliHFPtSpectrum.cxx:280
 AliHFPtSpectrum.cxx:281
 AliHFPtSpectrum.cxx:282
 AliHFPtSpectrum.cxx:283
 AliHFPtSpectrum.cxx:284
 AliHFPtSpectrum.cxx:285
 AliHFPtSpectrum.cxx:286
 AliHFPtSpectrum.cxx:287
 AliHFPtSpectrum.cxx:288
 AliHFPtSpectrum.cxx:289
 AliHFPtSpectrum.cxx:290
 AliHFPtSpectrum.cxx:291
 AliHFPtSpectrum.cxx:292
 AliHFPtSpectrum.cxx:293
 AliHFPtSpectrum.cxx:294
 AliHFPtSpectrum.cxx:295
 AliHFPtSpectrum.cxx:296
 AliHFPtSpectrum.cxx:297
 AliHFPtSpectrum.cxx:298
 AliHFPtSpectrum.cxx:299
 AliHFPtSpectrum.cxx:300
 AliHFPtSpectrum.cxx:301
 AliHFPtSpectrum.cxx:302
 AliHFPtSpectrum.cxx:303
 AliHFPtSpectrum.cxx:304
 AliHFPtSpectrum.cxx:305
 AliHFPtSpectrum.cxx:306
 AliHFPtSpectrum.cxx:307
 AliHFPtSpectrum.cxx:308
 AliHFPtSpectrum.cxx:309
 AliHFPtSpectrum.cxx:310
 AliHFPtSpectrum.cxx:311
 AliHFPtSpectrum.cxx:312
 AliHFPtSpectrum.cxx:313
 AliHFPtSpectrum.cxx:314
 AliHFPtSpectrum.cxx:315
 AliHFPtSpectrum.cxx:316
 AliHFPtSpectrum.cxx:317
 AliHFPtSpectrum.cxx:318
 AliHFPtSpectrum.cxx:319
 AliHFPtSpectrum.cxx:320
 AliHFPtSpectrum.cxx:321
 AliHFPtSpectrum.cxx:322
 AliHFPtSpectrum.cxx:323
 AliHFPtSpectrum.cxx:324
 AliHFPtSpectrum.cxx:325
 AliHFPtSpectrum.cxx:326
 AliHFPtSpectrum.cxx:327
 AliHFPtSpectrum.cxx:328
 AliHFPtSpectrum.cxx:329
 AliHFPtSpectrum.cxx:330
 AliHFPtSpectrum.cxx:331
 AliHFPtSpectrum.cxx:332
 AliHFPtSpectrum.cxx:333
 AliHFPtSpectrum.cxx:334
 AliHFPtSpectrum.cxx:335
 AliHFPtSpectrum.cxx:336
 AliHFPtSpectrum.cxx:337
 AliHFPtSpectrum.cxx:338
 AliHFPtSpectrum.cxx:339
 AliHFPtSpectrum.cxx:340
 AliHFPtSpectrum.cxx:341
 AliHFPtSpectrum.cxx:342
 AliHFPtSpectrum.cxx:343
 AliHFPtSpectrum.cxx:344
 AliHFPtSpectrum.cxx:345
 AliHFPtSpectrum.cxx:346
 AliHFPtSpectrum.cxx:347
 AliHFPtSpectrum.cxx:348
 AliHFPtSpectrum.cxx:349
 AliHFPtSpectrum.cxx:350
 AliHFPtSpectrum.cxx:351
 AliHFPtSpectrum.cxx:352
 AliHFPtSpectrum.cxx:353
 AliHFPtSpectrum.cxx:354
 AliHFPtSpectrum.cxx:355
 AliHFPtSpectrum.cxx:356
 AliHFPtSpectrum.cxx:357
 AliHFPtSpectrum.cxx:358
 AliHFPtSpectrum.cxx:359
 AliHFPtSpectrum.cxx:360
 AliHFPtSpectrum.cxx:361
 AliHFPtSpectrum.cxx:362
 AliHFPtSpectrum.cxx:363
 AliHFPtSpectrum.cxx:364
 AliHFPtSpectrum.cxx:365
 AliHFPtSpectrum.cxx:366
 AliHFPtSpectrum.cxx:367
 AliHFPtSpectrum.cxx:368
 AliHFPtSpectrum.cxx:369
 AliHFPtSpectrum.cxx:370
 AliHFPtSpectrum.cxx:371
 AliHFPtSpectrum.cxx:372
 AliHFPtSpectrum.cxx:373
 AliHFPtSpectrum.cxx:374
 AliHFPtSpectrum.cxx:375
 AliHFPtSpectrum.cxx:376
 AliHFPtSpectrum.cxx:377
 AliHFPtSpectrum.cxx:378
 AliHFPtSpectrum.cxx:379
 AliHFPtSpectrum.cxx:380
 AliHFPtSpectrum.cxx:381
 AliHFPtSpectrum.cxx:382
 AliHFPtSpectrum.cxx:383
 AliHFPtSpectrum.cxx:384
 AliHFPtSpectrum.cxx:385
 AliHFPtSpectrum.cxx:386
 AliHFPtSpectrum.cxx:387
 AliHFPtSpectrum.cxx:388
 AliHFPtSpectrum.cxx:389
 AliHFPtSpectrum.cxx:390
 AliHFPtSpectrum.cxx:391
 AliHFPtSpectrum.cxx:392
 AliHFPtSpectrum.cxx:393
 AliHFPtSpectrum.cxx:394
 AliHFPtSpectrum.cxx:395
 AliHFPtSpectrum.cxx:396
 AliHFPtSpectrum.cxx:397
 AliHFPtSpectrum.cxx:398
 AliHFPtSpectrum.cxx:399
 AliHFPtSpectrum.cxx:400
 AliHFPtSpectrum.cxx:401
 AliHFPtSpectrum.cxx:402
 AliHFPtSpectrum.cxx:403
 AliHFPtSpectrum.cxx:404
 AliHFPtSpectrum.cxx:405
 AliHFPtSpectrum.cxx:406
 AliHFPtSpectrum.cxx:407
 AliHFPtSpectrum.cxx:408
 AliHFPtSpectrum.cxx:409
 AliHFPtSpectrum.cxx:410
 AliHFPtSpectrum.cxx:411
 AliHFPtSpectrum.cxx:412
 AliHFPtSpectrum.cxx:413
 AliHFPtSpectrum.cxx:414
 AliHFPtSpectrum.cxx:415
 AliHFPtSpectrum.cxx:416
 AliHFPtSpectrum.cxx:417
 AliHFPtSpectrum.cxx:418
 AliHFPtSpectrum.cxx:419
 AliHFPtSpectrum.cxx:420
 AliHFPtSpectrum.cxx:421
 AliHFPtSpectrum.cxx:422
 AliHFPtSpectrum.cxx:423
 AliHFPtSpectrum.cxx:424
 AliHFPtSpectrum.cxx:425
 AliHFPtSpectrum.cxx:426
 AliHFPtSpectrum.cxx:427
 AliHFPtSpectrum.cxx:428
 AliHFPtSpectrum.cxx:429
 AliHFPtSpectrum.cxx:430
 AliHFPtSpectrum.cxx:431
 AliHFPtSpectrum.cxx:432
 AliHFPtSpectrum.cxx:433
 AliHFPtSpectrum.cxx:434
 AliHFPtSpectrum.cxx:435
 AliHFPtSpectrum.cxx:436
 AliHFPtSpectrum.cxx:437
 AliHFPtSpectrum.cxx:438
 AliHFPtSpectrum.cxx:439
 AliHFPtSpectrum.cxx:440
 AliHFPtSpectrum.cxx:441
 AliHFPtSpectrum.cxx:442
 AliHFPtSpectrum.cxx:443
 AliHFPtSpectrum.cxx:444
 AliHFPtSpectrum.cxx:445
 AliHFPtSpectrum.cxx:446
 AliHFPtSpectrum.cxx:447
 AliHFPtSpectrum.cxx:448
 AliHFPtSpectrum.cxx:449
 AliHFPtSpectrum.cxx:450
 AliHFPtSpectrum.cxx:451
 AliHFPtSpectrum.cxx:452
 AliHFPtSpectrum.cxx:453
 AliHFPtSpectrum.cxx:454
 AliHFPtSpectrum.cxx:455
 AliHFPtSpectrum.cxx:456
 AliHFPtSpectrum.cxx:457
 AliHFPtSpectrum.cxx:458
 AliHFPtSpectrum.cxx:459
 AliHFPtSpectrum.cxx:460
 AliHFPtSpectrum.cxx:461
 AliHFPtSpectrum.cxx:462
 AliHFPtSpectrum.cxx:463
 AliHFPtSpectrum.cxx:464
 AliHFPtSpectrum.cxx:465
 AliHFPtSpectrum.cxx:466
 AliHFPtSpectrum.cxx:467
 AliHFPtSpectrum.cxx:468
 AliHFPtSpectrum.cxx:469
 AliHFPtSpectrum.cxx:470
 AliHFPtSpectrum.cxx:471
 AliHFPtSpectrum.cxx:472
 AliHFPtSpectrum.cxx:473
 AliHFPtSpectrum.cxx:474
 AliHFPtSpectrum.cxx:475
 AliHFPtSpectrum.cxx:476
 AliHFPtSpectrum.cxx:477
 AliHFPtSpectrum.cxx:478
 AliHFPtSpectrum.cxx:479
 AliHFPtSpectrum.cxx:480
 AliHFPtSpectrum.cxx:481
 AliHFPtSpectrum.cxx:482
 AliHFPtSpectrum.cxx:483
 AliHFPtSpectrum.cxx:484
 AliHFPtSpectrum.cxx:485
 AliHFPtSpectrum.cxx:486
 AliHFPtSpectrum.cxx:487
 AliHFPtSpectrum.cxx:488
 AliHFPtSpectrum.cxx:489
 AliHFPtSpectrum.cxx:490
 AliHFPtSpectrum.cxx:491
 AliHFPtSpectrum.cxx:492
 AliHFPtSpectrum.cxx:493
 AliHFPtSpectrum.cxx:494
 AliHFPtSpectrum.cxx:495
 AliHFPtSpectrum.cxx:496
 AliHFPtSpectrum.cxx:497
 AliHFPtSpectrum.cxx:498
 AliHFPtSpectrum.cxx:499
 AliHFPtSpectrum.cxx:500
 AliHFPtSpectrum.cxx:501
 AliHFPtSpectrum.cxx:502
 AliHFPtSpectrum.cxx:503
 AliHFPtSpectrum.cxx:504
 AliHFPtSpectrum.cxx:505
 AliHFPtSpectrum.cxx:506
 AliHFPtSpectrum.cxx:507
 AliHFPtSpectrum.cxx:508
 AliHFPtSpectrum.cxx:509
 AliHFPtSpectrum.cxx:510
 AliHFPtSpectrum.cxx:511
 AliHFPtSpectrum.cxx:512
 AliHFPtSpectrum.cxx:513
 AliHFPtSpectrum.cxx:514
 AliHFPtSpectrum.cxx:515
 AliHFPtSpectrum.cxx:516
 AliHFPtSpectrum.cxx:517
 AliHFPtSpectrum.cxx:518
 AliHFPtSpectrum.cxx:519
 AliHFPtSpectrum.cxx:520
 AliHFPtSpectrum.cxx:521
 AliHFPtSpectrum.cxx:522
 AliHFPtSpectrum.cxx:523
 AliHFPtSpectrum.cxx:524
 AliHFPtSpectrum.cxx:525
 AliHFPtSpectrum.cxx:526
 AliHFPtSpectrum.cxx:527
 AliHFPtSpectrum.cxx:528
 AliHFPtSpectrum.cxx:529
 AliHFPtSpectrum.cxx:530
 AliHFPtSpectrum.cxx:531
 AliHFPtSpectrum.cxx:532
 AliHFPtSpectrum.cxx:533
 AliHFPtSpectrum.cxx:534
 AliHFPtSpectrum.cxx:535
 AliHFPtSpectrum.cxx:536
 AliHFPtSpectrum.cxx:537
 AliHFPtSpectrum.cxx:538
 AliHFPtSpectrum.cxx:539
 AliHFPtSpectrum.cxx:540
 AliHFPtSpectrum.cxx:541
 AliHFPtSpectrum.cxx:542
 AliHFPtSpectrum.cxx:543
 AliHFPtSpectrum.cxx:544
 AliHFPtSpectrum.cxx:545
 AliHFPtSpectrum.cxx:546
 AliHFPtSpectrum.cxx:547
 AliHFPtSpectrum.cxx:548
 AliHFPtSpectrum.cxx:549
 AliHFPtSpectrum.cxx:550
 AliHFPtSpectrum.cxx:551
 AliHFPtSpectrum.cxx:552
 AliHFPtSpectrum.cxx:553
 AliHFPtSpectrum.cxx:554
 AliHFPtSpectrum.cxx:555
 AliHFPtSpectrum.cxx:556
 AliHFPtSpectrum.cxx:557
 AliHFPtSpectrum.cxx:558
 AliHFPtSpectrum.cxx:559
 AliHFPtSpectrum.cxx:560
 AliHFPtSpectrum.cxx:561
 AliHFPtSpectrum.cxx:562
 AliHFPtSpectrum.cxx:563
 AliHFPtSpectrum.cxx:564
 AliHFPtSpectrum.cxx:565
 AliHFPtSpectrum.cxx:566
 AliHFPtSpectrum.cxx:567
 AliHFPtSpectrum.cxx:568
 AliHFPtSpectrum.cxx:569
 AliHFPtSpectrum.cxx:570
 AliHFPtSpectrum.cxx:571
 AliHFPtSpectrum.cxx:572
 AliHFPtSpectrum.cxx:573
 AliHFPtSpectrum.cxx:574
 AliHFPtSpectrum.cxx:575
 AliHFPtSpectrum.cxx:576
 AliHFPtSpectrum.cxx:577
 AliHFPtSpectrum.cxx:578
 AliHFPtSpectrum.cxx:579
 AliHFPtSpectrum.cxx:580
 AliHFPtSpectrum.cxx:581
 AliHFPtSpectrum.cxx:582
 AliHFPtSpectrum.cxx:583
 AliHFPtSpectrum.cxx:584
 AliHFPtSpectrum.cxx:585
 AliHFPtSpectrum.cxx:586
 AliHFPtSpectrum.cxx:587
 AliHFPtSpectrum.cxx:588
 AliHFPtSpectrum.cxx:589
 AliHFPtSpectrum.cxx:590
 AliHFPtSpectrum.cxx:591
 AliHFPtSpectrum.cxx:592
 AliHFPtSpectrum.cxx:593
 AliHFPtSpectrum.cxx:594
 AliHFPtSpectrum.cxx:595
 AliHFPtSpectrum.cxx:596
 AliHFPtSpectrum.cxx:597
 AliHFPtSpectrum.cxx:598
 AliHFPtSpectrum.cxx:599
 AliHFPtSpectrum.cxx:600
 AliHFPtSpectrum.cxx:601
 AliHFPtSpectrum.cxx:602
 AliHFPtSpectrum.cxx:603
 AliHFPtSpectrum.cxx:604
 AliHFPtSpectrum.cxx:605
 AliHFPtSpectrum.cxx:606
 AliHFPtSpectrum.cxx:607
 AliHFPtSpectrum.cxx:608
 AliHFPtSpectrum.cxx:609
 AliHFPtSpectrum.cxx:610
 AliHFPtSpectrum.cxx:611
 AliHFPtSpectrum.cxx:612
 AliHFPtSpectrum.cxx:613
 AliHFPtSpectrum.cxx:614
 AliHFPtSpectrum.cxx:615
 AliHFPtSpectrum.cxx:616
 AliHFPtSpectrum.cxx:617
 AliHFPtSpectrum.cxx:618
 AliHFPtSpectrum.cxx:619
 AliHFPtSpectrum.cxx:620
 AliHFPtSpectrum.cxx:621
 AliHFPtSpectrum.cxx:622
 AliHFPtSpectrum.cxx:623
 AliHFPtSpectrum.cxx:624
 AliHFPtSpectrum.cxx:625
 AliHFPtSpectrum.cxx:626
 AliHFPtSpectrum.cxx:627
 AliHFPtSpectrum.cxx:628
 AliHFPtSpectrum.cxx:629
 AliHFPtSpectrum.cxx:630
 AliHFPtSpectrum.cxx:631
 AliHFPtSpectrum.cxx:632
 AliHFPtSpectrum.cxx:633
 AliHFPtSpectrum.cxx:634
 AliHFPtSpectrum.cxx:635
 AliHFPtSpectrum.cxx:636
 AliHFPtSpectrum.cxx:637
 AliHFPtSpectrum.cxx:638
 AliHFPtSpectrum.cxx:639
 AliHFPtSpectrum.cxx:640
 AliHFPtSpectrum.cxx:641
 AliHFPtSpectrum.cxx:642
 AliHFPtSpectrum.cxx:643
 AliHFPtSpectrum.cxx:644
 AliHFPtSpectrum.cxx:645
 AliHFPtSpectrum.cxx:646
 AliHFPtSpectrum.cxx:647
 AliHFPtSpectrum.cxx:648
 AliHFPtSpectrum.cxx:649
 AliHFPtSpectrum.cxx:650
 AliHFPtSpectrum.cxx:651
 AliHFPtSpectrum.cxx:652
 AliHFPtSpectrum.cxx:653
 AliHFPtSpectrum.cxx:654
 AliHFPtSpectrum.cxx:655
 AliHFPtSpectrum.cxx:656
 AliHFPtSpectrum.cxx:657
 AliHFPtSpectrum.cxx:658
 AliHFPtSpectrum.cxx:659
 AliHFPtSpectrum.cxx:660
 AliHFPtSpectrum.cxx:661
 AliHFPtSpectrum.cxx:662
 AliHFPtSpectrum.cxx:663
 AliHFPtSpectrum.cxx:664
 AliHFPtSpectrum.cxx:665
 AliHFPtSpectrum.cxx:666
 AliHFPtSpectrum.cxx:667
 AliHFPtSpectrum.cxx:668
 AliHFPtSpectrum.cxx:669
 AliHFPtSpectrum.cxx:670
 AliHFPtSpectrum.cxx:671
 AliHFPtSpectrum.cxx:672
 AliHFPtSpectrum.cxx:673
 AliHFPtSpectrum.cxx:674
 AliHFPtSpectrum.cxx:675
 AliHFPtSpectrum.cxx:676
 AliHFPtSpectrum.cxx:677
 AliHFPtSpectrum.cxx:678
 AliHFPtSpectrum.cxx:679
 AliHFPtSpectrum.cxx:680
 AliHFPtSpectrum.cxx:681
 AliHFPtSpectrum.cxx:682
 AliHFPtSpectrum.cxx:683
 AliHFPtSpectrum.cxx:684
 AliHFPtSpectrum.cxx:685
 AliHFPtSpectrum.cxx:686
 AliHFPtSpectrum.cxx:687
 AliHFPtSpectrum.cxx:688
 AliHFPtSpectrum.cxx:689
 AliHFPtSpectrum.cxx:690
 AliHFPtSpectrum.cxx:691
 AliHFPtSpectrum.cxx:692
 AliHFPtSpectrum.cxx:693
 AliHFPtSpectrum.cxx:694
 AliHFPtSpectrum.cxx:695
 AliHFPtSpectrum.cxx:696
 AliHFPtSpectrum.cxx:697
 AliHFPtSpectrum.cxx:698
 AliHFPtSpectrum.cxx:699
 AliHFPtSpectrum.cxx:700
 AliHFPtSpectrum.cxx:701
 AliHFPtSpectrum.cxx:702
 AliHFPtSpectrum.cxx:703
 AliHFPtSpectrum.cxx:704
 AliHFPtSpectrum.cxx:705
 AliHFPtSpectrum.cxx:706
 AliHFPtSpectrum.cxx:707
 AliHFPtSpectrum.cxx:708
 AliHFPtSpectrum.cxx:709
 AliHFPtSpectrum.cxx:710
 AliHFPtSpectrum.cxx:711
 AliHFPtSpectrum.cxx:712
 AliHFPtSpectrum.cxx:713
 AliHFPtSpectrum.cxx:714
 AliHFPtSpectrum.cxx:715
 AliHFPtSpectrum.cxx:716
 AliHFPtSpectrum.cxx:717
 AliHFPtSpectrum.cxx:718
 AliHFPtSpectrum.cxx:719
 AliHFPtSpectrum.cxx:720
 AliHFPtSpectrum.cxx:721
 AliHFPtSpectrum.cxx:722
 AliHFPtSpectrum.cxx:723
 AliHFPtSpectrum.cxx:724
 AliHFPtSpectrum.cxx:725
 AliHFPtSpectrum.cxx:726
 AliHFPtSpectrum.cxx:727
 AliHFPtSpectrum.cxx:728
 AliHFPtSpectrum.cxx:729
 AliHFPtSpectrum.cxx:730
 AliHFPtSpectrum.cxx:731
 AliHFPtSpectrum.cxx:732
 AliHFPtSpectrum.cxx:733
 AliHFPtSpectrum.cxx:734
 AliHFPtSpectrum.cxx:735
 AliHFPtSpectrum.cxx:736
 AliHFPtSpectrum.cxx:737
 AliHFPtSpectrum.cxx:738
 AliHFPtSpectrum.cxx:739
 AliHFPtSpectrum.cxx:740
 AliHFPtSpectrum.cxx:741
 AliHFPtSpectrum.cxx:742
 AliHFPtSpectrum.cxx:743
 AliHFPtSpectrum.cxx:744
 AliHFPtSpectrum.cxx:745
 AliHFPtSpectrum.cxx:746
 AliHFPtSpectrum.cxx:747
 AliHFPtSpectrum.cxx:748
 AliHFPtSpectrum.cxx:749
 AliHFPtSpectrum.cxx:750
 AliHFPtSpectrum.cxx:751
 AliHFPtSpectrum.cxx:752
 AliHFPtSpectrum.cxx:753
 AliHFPtSpectrum.cxx:754
 AliHFPtSpectrum.cxx:755
 AliHFPtSpectrum.cxx:756
 AliHFPtSpectrum.cxx:757
 AliHFPtSpectrum.cxx:758
 AliHFPtSpectrum.cxx:759
 AliHFPtSpectrum.cxx:760
 AliHFPtSpectrum.cxx:761
 AliHFPtSpectrum.cxx:762
 AliHFPtSpectrum.cxx:763
 AliHFPtSpectrum.cxx:764
 AliHFPtSpectrum.cxx:765
 AliHFPtSpectrum.cxx:766
 AliHFPtSpectrum.cxx:767
 AliHFPtSpectrum.cxx:768
 AliHFPtSpectrum.cxx:769
 AliHFPtSpectrum.cxx:770
 AliHFPtSpectrum.cxx:771
 AliHFPtSpectrum.cxx:772
 AliHFPtSpectrum.cxx:773
 AliHFPtSpectrum.cxx:774
 AliHFPtSpectrum.cxx:775
 AliHFPtSpectrum.cxx:776
 AliHFPtSpectrum.cxx:777
 AliHFPtSpectrum.cxx:778
 AliHFPtSpectrum.cxx:779
 AliHFPtSpectrum.cxx:780
 AliHFPtSpectrum.cxx:781
 AliHFPtSpectrum.cxx:782
 AliHFPtSpectrum.cxx:783
 AliHFPtSpectrum.cxx:784
 AliHFPtSpectrum.cxx:785
 AliHFPtSpectrum.cxx:786
 AliHFPtSpectrum.cxx:787
 AliHFPtSpectrum.cxx:788
 AliHFPtSpectrum.cxx:789
 AliHFPtSpectrum.cxx:790
 AliHFPtSpectrum.cxx:791
 AliHFPtSpectrum.cxx:792
 AliHFPtSpectrum.cxx:793
 AliHFPtSpectrum.cxx:794
 AliHFPtSpectrum.cxx:795
 AliHFPtSpectrum.cxx:796
 AliHFPtSpectrum.cxx:797
 AliHFPtSpectrum.cxx:798
 AliHFPtSpectrum.cxx:799
 AliHFPtSpectrum.cxx:800
 AliHFPtSpectrum.cxx:801
 AliHFPtSpectrum.cxx:802
 AliHFPtSpectrum.cxx:803
 AliHFPtSpectrum.cxx:804
 AliHFPtSpectrum.cxx:805
 AliHFPtSpectrum.cxx:806
 AliHFPtSpectrum.cxx:807
 AliHFPtSpectrum.cxx:808
 AliHFPtSpectrum.cxx:809
 AliHFPtSpectrum.cxx:810
 AliHFPtSpectrum.cxx:811
 AliHFPtSpectrum.cxx:812
 AliHFPtSpectrum.cxx:813
 AliHFPtSpectrum.cxx:814
 AliHFPtSpectrum.cxx:815
 AliHFPtSpectrum.cxx:816
 AliHFPtSpectrum.cxx:817
 AliHFPtSpectrum.cxx:818
 AliHFPtSpectrum.cxx:819
 AliHFPtSpectrum.cxx:820
 AliHFPtSpectrum.cxx:821
 AliHFPtSpectrum.cxx:822
 AliHFPtSpectrum.cxx:823
 AliHFPtSpectrum.cxx:824
 AliHFPtSpectrum.cxx:825
 AliHFPtSpectrum.cxx:826
 AliHFPtSpectrum.cxx:827
 AliHFPtSpectrum.cxx:828
 AliHFPtSpectrum.cxx:829
 AliHFPtSpectrum.cxx:830
 AliHFPtSpectrum.cxx:831
 AliHFPtSpectrum.cxx:832
 AliHFPtSpectrum.cxx:833
 AliHFPtSpectrum.cxx:834
 AliHFPtSpectrum.cxx:835
 AliHFPtSpectrum.cxx:836
 AliHFPtSpectrum.cxx:837
 AliHFPtSpectrum.cxx:838
 AliHFPtSpectrum.cxx:839
 AliHFPtSpectrum.cxx:840
 AliHFPtSpectrum.cxx:841
 AliHFPtSpectrum.cxx:842
 AliHFPtSpectrum.cxx:843
 AliHFPtSpectrum.cxx:844
 AliHFPtSpectrum.cxx:845
 AliHFPtSpectrum.cxx:846
 AliHFPtSpectrum.cxx:847
 AliHFPtSpectrum.cxx:848
 AliHFPtSpectrum.cxx:849
 AliHFPtSpectrum.cxx:850
 AliHFPtSpectrum.cxx:851
 AliHFPtSpectrum.cxx:852
 AliHFPtSpectrum.cxx:853
 AliHFPtSpectrum.cxx:854
 AliHFPtSpectrum.cxx:855
 AliHFPtSpectrum.cxx:856
 AliHFPtSpectrum.cxx:857
 AliHFPtSpectrum.cxx:858
 AliHFPtSpectrum.cxx:859
 AliHFPtSpectrum.cxx:860
 AliHFPtSpectrum.cxx:861
 AliHFPtSpectrum.cxx:862
 AliHFPtSpectrum.cxx:863
 AliHFPtSpectrum.cxx:864
 AliHFPtSpectrum.cxx:865
 AliHFPtSpectrum.cxx:866
 AliHFPtSpectrum.cxx:867
 AliHFPtSpectrum.cxx:868
 AliHFPtSpectrum.cxx:869
 AliHFPtSpectrum.cxx:870
 AliHFPtSpectrum.cxx:871
 AliHFPtSpectrum.cxx:872
 AliHFPtSpectrum.cxx:873
 AliHFPtSpectrum.cxx:874
 AliHFPtSpectrum.cxx:875
 AliHFPtSpectrum.cxx:876
 AliHFPtSpectrum.cxx:877
 AliHFPtSpectrum.cxx:878
 AliHFPtSpectrum.cxx:879
 AliHFPtSpectrum.cxx:880
 AliHFPtSpectrum.cxx:881
 AliHFPtSpectrum.cxx:882
 AliHFPtSpectrum.cxx:883
 AliHFPtSpectrum.cxx:884
 AliHFPtSpectrum.cxx:885
 AliHFPtSpectrum.cxx:886
 AliHFPtSpectrum.cxx:887
 AliHFPtSpectrum.cxx:888
 AliHFPtSpectrum.cxx:889
 AliHFPtSpectrum.cxx:890
 AliHFPtSpectrum.cxx:891
 AliHFPtSpectrum.cxx:892
 AliHFPtSpectrum.cxx:893
 AliHFPtSpectrum.cxx:894
 AliHFPtSpectrum.cxx:895
 AliHFPtSpectrum.cxx:896
 AliHFPtSpectrum.cxx:897
 AliHFPtSpectrum.cxx:898
 AliHFPtSpectrum.cxx:899
 AliHFPtSpectrum.cxx:900
 AliHFPtSpectrum.cxx:901
 AliHFPtSpectrum.cxx:902
 AliHFPtSpectrum.cxx:903
 AliHFPtSpectrum.cxx:904
 AliHFPtSpectrum.cxx:905
 AliHFPtSpectrum.cxx:906
 AliHFPtSpectrum.cxx:907
 AliHFPtSpectrum.cxx:908
 AliHFPtSpectrum.cxx:909
 AliHFPtSpectrum.cxx:910
 AliHFPtSpectrum.cxx:911
 AliHFPtSpectrum.cxx:912
 AliHFPtSpectrum.cxx:913
 AliHFPtSpectrum.cxx:914
 AliHFPtSpectrum.cxx:915
 AliHFPtSpectrum.cxx:916
 AliHFPtSpectrum.cxx:917
 AliHFPtSpectrum.cxx:918
 AliHFPtSpectrum.cxx:919
 AliHFPtSpectrum.cxx:920
 AliHFPtSpectrum.cxx:921
 AliHFPtSpectrum.cxx:922
 AliHFPtSpectrum.cxx:923
 AliHFPtSpectrum.cxx:924
 AliHFPtSpectrum.cxx:925
 AliHFPtSpectrum.cxx:926
 AliHFPtSpectrum.cxx:927
 AliHFPtSpectrum.cxx:928
 AliHFPtSpectrum.cxx:929
 AliHFPtSpectrum.cxx:930
 AliHFPtSpectrum.cxx:931
 AliHFPtSpectrum.cxx:932
 AliHFPtSpectrum.cxx:933
 AliHFPtSpectrum.cxx:934
 AliHFPtSpectrum.cxx:935
 AliHFPtSpectrum.cxx:936
 AliHFPtSpectrum.cxx:937
 AliHFPtSpectrum.cxx:938
 AliHFPtSpectrum.cxx:939
 AliHFPtSpectrum.cxx:940
 AliHFPtSpectrum.cxx:941
 AliHFPtSpectrum.cxx:942
 AliHFPtSpectrum.cxx:943
 AliHFPtSpectrum.cxx:944
 AliHFPtSpectrum.cxx:945
 AliHFPtSpectrum.cxx:946
 AliHFPtSpectrum.cxx:947
 AliHFPtSpectrum.cxx:948
 AliHFPtSpectrum.cxx:949
 AliHFPtSpectrum.cxx:950
 AliHFPtSpectrum.cxx:951
 AliHFPtSpectrum.cxx:952
 AliHFPtSpectrum.cxx:953
 AliHFPtSpectrum.cxx:954
 AliHFPtSpectrum.cxx:955
 AliHFPtSpectrum.cxx:956
 AliHFPtSpectrum.cxx:957
 AliHFPtSpectrum.cxx:958
 AliHFPtSpectrum.cxx:959
 AliHFPtSpectrum.cxx:960
 AliHFPtSpectrum.cxx:961
 AliHFPtSpectrum.cxx:962
 AliHFPtSpectrum.cxx:963
 AliHFPtSpectrum.cxx:964
 AliHFPtSpectrum.cxx:965
 AliHFPtSpectrum.cxx:966
 AliHFPtSpectrum.cxx:967
 AliHFPtSpectrum.cxx:968
 AliHFPtSpectrum.cxx:969
 AliHFPtSpectrum.cxx:970
 AliHFPtSpectrum.cxx:971
 AliHFPtSpectrum.cxx:972
 AliHFPtSpectrum.cxx:973
 AliHFPtSpectrum.cxx:974
 AliHFPtSpectrum.cxx:975
 AliHFPtSpectrum.cxx:976
 AliHFPtSpectrum.cxx:977
 AliHFPtSpectrum.cxx:978
 AliHFPtSpectrum.cxx:979
 AliHFPtSpectrum.cxx:980
 AliHFPtSpectrum.cxx:981
 AliHFPtSpectrum.cxx:982
 AliHFPtSpectrum.cxx:983
 AliHFPtSpectrum.cxx:984
 AliHFPtSpectrum.cxx:985
 AliHFPtSpectrum.cxx:986
 AliHFPtSpectrum.cxx:987
 AliHFPtSpectrum.cxx:988
 AliHFPtSpectrum.cxx:989
 AliHFPtSpectrum.cxx:990
 AliHFPtSpectrum.cxx:991
 AliHFPtSpectrum.cxx:992
 AliHFPtSpectrum.cxx:993
 AliHFPtSpectrum.cxx:994
 AliHFPtSpectrum.cxx:995
 AliHFPtSpectrum.cxx:996
 AliHFPtSpectrum.cxx:997
 AliHFPtSpectrum.cxx:998
 AliHFPtSpectrum.cxx:999
 AliHFPtSpectrum.cxx:1000
 AliHFPtSpectrum.cxx:1001
 AliHFPtSpectrum.cxx:1002
 AliHFPtSpectrum.cxx:1003
 AliHFPtSpectrum.cxx:1004
 AliHFPtSpectrum.cxx:1005
 AliHFPtSpectrum.cxx:1006
 AliHFPtSpectrum.cxx:1007
 AliHFPtSpectrum.cxx:1008
 AliHFPtSpectrum.cxx:1009
 AliHFPtSpectrum.cxx:1010
 AliHFPtSpectrum.cxx:1011
 AliHFPtSpectrum.cxx:1012
 AliHFPtSpectrum.cxx:1013
 AliHFPtSpectrum.cxx:1014
 AliHFPtSpectrum.cxx:1015
 AliHFPtSpectrum.cxx:1016
 AliHFPtSpectrum.cxx:1017
 AliHFPtSpectrum.cxx:1018
 AliHFPtSpectrum.cxx:1019
 AliHFPtSpectrum.cxx:1020
 AliHFPtSpectrum.cxx:1021
 AliHFPtSpectrum.cxx:1022
 AliHFPtSpectrum.cxx:1023
 AliHFPtSpectrum.cxx:1024
 AliHFPtSpectrum.cxx:1025
 AliHFPtSpectrum.cxx:1026
 AliHFPtSpectrum.cxx:1027
 AliHFPtSpectrum.cxx:1028
 AliHFPtSpectrum.cxx:1029
 AliHFPtSpectrum.cxx:1030
 AliHFPtSpectrum.cxx:1031
 AliHFPtSpectrum.cxx:1032
 AliHFPtSpectrum.cxx:1033
 AliHFPtSpectrum.cxx:1034
 AliHFPtSpectrum.cxx:1035
 AliHFPtSpectrum.cxx:1036
 AliHFPtSpectrum.cxx:1037
 AliHFPtSpectrum.cxx:1038
 AliHFPtSpectrum.cxx:1039
 AliHFPtSpectrum.cxx:1040
 AliHFPtSpectrum.cxx:1041
 AliHFPtSpectrum.cxx:1042
 AliHFPtSpectrum.cxx:1043
 AliHFPtSpectrum.cxx:1044
 AliHFPtSpectrum.cxx:1045
 AliHFPtSpectrum.cxx:1046
 AliHFPtSpectrum.cxx:1047
 AliHFPtSpectrum.cxx:1048
 AliHFPtSpectrum.cxx:1049
 AliHFPtSpectrum.cxx:1050
 AliHFPtSpectrum.cxx:1051
 AliHFPtSpectrum.cxx:1052
 AliHFPtSpectrum.cxx:1053
 AliHFPtSpectrum.cxx:1054
 AliHFPtSpectrum.cxx:1055
 AliHFPtSpectrum.cxx:1056
 AliHFPtSpectrum.cxx:1057
 AliHFPtSpectrum.cxx:1058
 AliHFPtSpectrum.cxx:1059
 AliHFPtSpectrum.cxx:1060
 AliHFPtSpectrum.cxx:1061
 AliHFPtSpectrum.cxx:1062
 AliHFPtSpectrum.cxx:1063
 AliHFPtSpectrum.cxx:1064
 AliHFPtSpectrum.cxx:1065
 AliHFPtSpectrum.cxx:1066
 AliHFPtSpectrum.cxx:1067
 AliHFPtSpectrum.cxx:1068
 AliHFPtSpectrum.cxx:1069
 AliHFPtSpectrum.cxx:1070
 AliHFPtSpectrum.cxx:1071
 AliHFPtSpectrum.cxx:1072
 AliHFPtSpectrum.cxx:1073
 AliHFPtSpectrum.cxx:1074
 AliHFPtSpectrum.cxx:1075
 AliHFPtSpectrum.cxx:1076
 AliHFPtSpectrum.cxx:1077
 AliHFPtSpectrum.cxx:1078
 AliHFPtSpectrum.cxx:1079
 AliHFPtSpectrum.cxx:1080
 AliHFPtSpectrum.cxx:1081
 AliHFPtSpectrum.cxx:1082
 AliHFPtSpectrum.cxx:1083
 AliHFPtSpectrum.cxx:1084
 AliHFPtSpectrum.cxx:1085
 AliHFPtSpectrum.cxx:1086
 AliHFPtSpectrum.cxx:1087
 AliHFPtSpectrum.cxx:1088
 AliHFPtSpectrum.cxx:1089
 AliHFPtSpectrum.cxx:1090
 AliHFPtSpectrum.cxx:1091
 AliHFPtSpectrum.cxx:1092
 AliHFPtSpectrum.cxx:1093
 AliHFPtSpectrum.cxx:1094
 AliHFPtSpectrum.cxx:1095
 AliHFPtSpectrum.cxx:1096
 AliHFPtSpectrum.cxx:1097
 AliHFPtSpectrum.cxx:1098
 AliHFPtSpectrum.cxx:1099
 AliHFPtSpectrum.cxx:1100
 AliHFPtSpectrum.cxx:1101
 AliHFPtSpectrum.cxx:1102
 AliHFPtSpectrum.cxx:1103
 AliHFPtSpectrum.cxx:1104
 AliHFPtSpectrum.cxx:1105
 AliHFPtSpectrum.cxx:1106
 AliHFPtSpectrum.cxx:1107
 AliHFPtSpectrum.cxx:1108
 AliHFPtSpectrum.cxx:1109
 AliHFPtSpectrum.cxx:1110
 AliHFPtSpectrum.cxx:1111
 AliHFPtSpectrum.cxx:1112
 AliHFPtSpectrum.cxx:1113
 AliHFPtSpectrum.cxx:1114
 AliHFPtSpectrum.cxx:1115
 AliHFPtSpectrum.cxx:1116
 AliHFPtSpectrum.cxx:1117
 AliHFPtSpectrum.cxx:1118
 AliHFPtSpectrum.cxx:1119
 AliHFPtSpectrum.cxx:1120
 AliHFPtSpectrum.cxx:1121
 AliHFPtSpectrum.cxx:1122
 AliHFPtSpectrum.cxx:1123
 AliHFPtSpectrum.cxx:1124
 AliHFPtSpectrum.cxx:1125
 AliHFPtSpectrum.cxx:1126
 AliHFPtSpectrum.cxx:1127
 AliHFPtSpectrum.cxx:1128
 AliHFPtSpectrum.cxx:1129
 AliHFPtSpectrum.cxx:1130
 AliHFPtSpectrum.cxx:1131
 AliHFPtSpectrum.cxx:1132
 AliHFPtSpectrum.cxx:1133
 AliHFPtSpectrum.cxx:1134
 AliHFPtSpectrum.cxx:1135
 AliHFPtSpectrum.cxx:1136
 AliHFPtSpectrum.cxx:1137
 AliHFPtSpectrum.cxx:1138
 AliHFPtSpectrum.cxx:1139
 AliHFPtSpectrum.cxx:1140
 AliHFPtSpectrum.cxx:1141
 AliHFPtSpectrum.cxx:1142
 AliHFPtSpectrum.cxx:1143
 AliHFPtSpectrum.cxx:1144
 AliHFPtSpectrum.cxx:1145
 AliHFPtSpectrum.cxx:1146
 AliHFPtSpectrum.cxx:1147
 AliHFPtSpectrum.cxx:1148
 AliHFPtSpectrum.cxx:1149
 AliHFPtSpectrum.cxx:1150
 AliHFPtSpectrum.cxx:1151
 AliHFPtSpectrum.cxx:1152
 AliHFPtSpectrum.cxx:1153
 AliHFPtSpectrum.cxx:1154
 AliHFPtSpectrum.cxx:1155
 AliHFPtSpectrum.cxx:1156
 AliHFPtSpectrum.cxx:1157
 AliHFPtSpectrum.cxx:1158
 AliHFPtSpectrum.cxx:1159
 AliHFPtSpectrum.cxx:1160
 AliHFPtSpectrum.cxx:1161
 AliHFPtSpectrum.cxx:1162
 AliHFPtSpectrum.cxx:1163
 AliHFPtSpectrum.cxx:1164
 AliHFPtSpectrum.cxx:1165
 AliHFPtSpectrum.cxx:1166
 AliHFPtSpectrum.cxx:1167
 AliHFPtSpectrum.cxx:1168
 AliHFPtSpectrum.cxx:1169
 AliHFPtSpectrum.cxx:1170
 AliHFPtSpectrum.cxx:1171
 AliHFPtSpectrum.cxx:1172
 AliHFPtSpectrum.cxx:1173
 AliHFPtSpectrum.cxx:1174
 AliHFPtSpectrum.cxx:1175
 AliHFPtSpectrum.cxx:1176
 AliHFPtSpectrum.cxx:1177
 AliHFPtSpectrum.cxx:1178
 AliHFPtSpectrum.cxx:1179
 AliHFPtSpectrum.cxx:1180
 AliHFPtSpectrum.cxx:1181
 AliHFPtSpectrum.cxx:1182
 AliHFPtSpectrum.cxx:1183
 AliHFPtSpectrum.cxx:1184
 AliHFPtSpectrum.cxx:1185
 AliHFPtSpectrum.cxx:1186
 AliHFPtSpectrum.cxx:1187
 AliHFPtSpectrum.cxx:1188
 AliHFPtSpectrum.cxx:1189
 AliHFPtSpectrum.cxx:1190
 AliHFPtSpectrum.cxx:1191
 AliHFPtSpectrum.cxx:1192
 AliHFPtSpectrum.cxx:1193
 AliHFPtSpectrum.cxx:1194
 AliHFPtSpectrum.cxx:1195
 AliHFPtSpectrum.cxx:1196
 AliHFPtSpectrum.cxx:1197
 AliHFPtSpectrum.cxx:1198
 AliHFPtSpectrum.cxx:1199
 AliHFPtSpectrum.cxx:1200
 AliHFPtSpectrum.cxx:1201
 AliHFPtSpectrum.cxx:1202
 AliHFPtSpectrum.cxx:1203
 AliHFPtSpectrum.cxx:1204
 AliHFPtSpectrum.cxx:1205
 AliHFPtSpectrum.cxx:1206
 AliHFPtSpectrum.cxx:1207
 AliHFPtSpectrum.cxx:1208
 AliHFPtSpectrum.cxx:1209
 AliHFPtSpectrum.cxx:1210
 AliHFPtSpectrum.cxx:1211
 AliHFPtSpectrum.cxx:1212
 AliHFPtSpectrum.cxx:1213
 AliHFPtSpectrum.cxx:1214
 AliHFPtSpectrum.cxx:1215
 AliHFPtSpectrum.cxx:1216
 AliHFPtSpectrum.cxx:1217
 AliHFPtSpectrum.cxx:1218
 AliHFPtSpectrum.cxx:1219
 AliHFPtSpectrum.cxx:1220
 AliHFPtSpectrum.cxx:1221
 AliHFPtSpectrum.cxx:1222
 AliHFPtSpectrum.cxx:1223
 AliHFPtSpectrum.cxx:1224
 AliHFPtSpectrum.cxx:1225
 AliHFPtSpectrum.cxx:1226
 AliHFPtSpectrum.cxx:1227
 AliHFPtSpectrum.cxx:1228
 AliHFPtSpectrum.cxx:1229
 AliHFPtSpectrum.cxx:1230
 AliHFPtSpectrum.cxx:1231
 AliHFPtSpectrum.cxx:1232
 AliHFPtSpectrum.cxx:1233
 AliHFPtSpectrum.cxx:1234
 AliHFPtSpectrum.cxx:1235
 AliHFPtSpectrum.cxx:1236
 AliHFPtSpectrum.cxx:1237
 AliHFPtSpectrum.cxx:1238
 AliHFPtSpectrum.cxx:1239
 AliHFPtSpectrum.cxx:1240
 AliHFPtSpectrum.cxx:1241
 AliHFPtSpectrum.cxx:1242
 AliHFPtSpectrum.cxx:1243
 AliHFPtSpectrum.cxx:1244
 AliHFPtSpectrum.cxx:1245
 AliHFPtSpectrum.cxx:1246
 AliHFPtSpectrum.cxx:1247
 AliHFPtSpectrum.cxx:1248
 AliHFPtSpectrum.cxx:1249
 AliHFPtSpectrum.cxx:1250
 AliHFPtSpectrum.cxx:1251
 AliHFPtSpectrum.cxx:1252
 AliHFPtSpectrum.cxx:1253
 AliHFPtSpectrum.cxx:1254
 AliHFPtSpectrum.cxx:1255
 AliHFPtSpectrum.cxx:1256
 AliHFPtSpectrum.cxx:1257
 AliHFPtSpectrum.cxx:1258
 AliHFPtSpectrum.cxx:1259
 AliHFPtSpectrum.cxx:1260
 AliHFPtSpectrum.cxx:1261
 AliHFPtSpectrum.cxx:1262
 AliHFPtSpectrum.cxx:1263
 AliHFPtSpectrum.cxx:1264
 AliHFPtSpectrum.cxx:1265
 AliHFPtSpectrum.cxx:1266
 AliHFPtSpectrum.cxx:1267
 AliHFPtSpectrum.cxx:1268
 AliHFPtSpectrum.cxx:1269
 AliHFPtSpectrum.cxx:1270
 AliHFPtSpectrum.cxx:1271
 AliHFPtSpectrum.cxx:1272
 AliHFPtSpectrum.cxx:1273
 AliHFPtSpectrum.cxx:1274
 AliHFPtSpectrum.cxx:1275
 AliHFPtSpectrum.cxx:1276
 AliHFPtSpectrum.cxx:1277
 AliHFPtSpectrum.cxx:1278
 AliHFPtSpectrum.cxx:1279
 AliHFPtSpectrum.cxx:1280
 AliHFPtSpectrum.cxx:1281
 AliHFPtSpectrum.cxx:1282
 AliHFPtSpectrum.cxx:1283
 AliHFPtSpectrum.cxx:1284
 AliHFPtSpectrum.cxx:1285
 AliHFPtSpectrum.cxx:1286
 AliHFPtSpectrum.cxx:1287
 AliHFPtSpectrum.cxx:1288
 AliHFPtSpectrum.cxx:1289
 AliHFPtSpectrum.cxx:1290
 AliHFPtSpectrum.cxx:1291
 AliHFPtSpectrum.cxx:1292
 AliHFPtSpectrum.cxx:1293
 AliHFPtSpectrum.cxx:1294
 AliHFPtSpectrum.cxx:1295
 AliHFPtSpectrum.cxx:1296
 AliHFPtSpectrum.cxx:1297
 AliHFPtSpectrum.cxx:1298
 AliHFPtSpectrum.cxx:1299
 AliHFPtSpectrum.cxx:1300
 AliHFPtSpectrum.cxx:1301
 AliHFPtSpectrum.cxx:1302
 AliHFPtSpectrum.cxx:1303
 AliHFPtSpectrum.cxx:1304
 AliHFPtSpectrum.cxx:1305
 AliHFPtSpectrum.cxx:1306
 AliHFPtSpectrum.cxx:1307
 AliHFPtSpectrum.cxx:1308
 AliHFPtSpectrum.cxx:1309
 AliHFPtSpectrum.cxx:1310
 AliHFPtSpectrum.cxx:1311
 AliHFPtSpectrum.cxx:1312
 AliHFPtSpectrum.cxx:1313
 AliHFPtSpectrum.cxx:1314
 AliHFPtSpectrum.cxx:1315
 AliHFPtSpectrum.cxx:1316
 AliHFPtSpectrum.cxx:1317
 AliHFPtSpectrum.cxx:1318
 AliHFPtSpectrum.cxx:1319
 AliHFPtSpectrum.cxx:1320
 AliHFPtSpectrum.cxx:1321
 AliHFPtSpectrum.cxx:1322
 AliHFPtSpectrum.cxx:1323
 AliHFPtSpectrum.cxx:1324
 AliHFPtSpectrum.cxx:1325
 AliHFPtSpectrum.cxx:1326
 AliHFPtSpectrum.cxx:1327
 AliHFPtSpectrum.cxx:1328
 AliHFPtSpectrum.cxx:1329
 AliHFPtSpectrum.cxx:1330
 AliHFPtSpectrum.cxx:1331
 AliHFPtSpectrum.cxx:1332
 AliHFPtSpectrum.cxx:1333
 AliHFPtSpectrum.cxx:1334
 AliHFPtSpectrum.cxx:1335
 AliHFPtSpectrum.cxx:1336
 AliHFPtSpectrum.cxx:1337
 AliHFPtSpectrum.cxx:1338
 AliHFPtSpectrum.cxx:1339
 AliHFPtSpectrum.cxx:1340
 AliHFPtSpectrum.cxx:1341
 AliHFPtSpectrum.cxx:1342
 AliHFPtSpectrum.cxx:1343
 AliHFPtSpectrum.cxx:1344
 AliHFPtSpectrum.cxx:1345
 AliHFPtSpectrum.cxx:1346
 AliHFPtSpectrum.cxx:1347
 AliHFPtSpectrum.cxx:1348
 AliHFPtSpectrum.cxx:1349
 AliHFPtSpectrum.cxx:1350
 AliHFPtSpectrum.cxx:1351
 AliHFPtSpectrum.cxx:1352
 AliHFPtSpectrum.cxx:1353
 AliHFPtSpectrum.cxx:1354
 AliHFPtSpectrum.cxx:1355
 AliHFPtSpectrum.cxx:1356
 AliHFPtSpectrum.cxx:1357
 AliHFPtSpectrum.cxx:1358
 AliHFPtSpectrum.cxx:1359
 AliHFPtSpectrum.cxx:1360
 AliHFPtSpectrum.cxx:1361
 AliHFPtSpectrum.cxx:1362
 AliHFPtSpectrum.cxx:1363
 AliHFPtSpectrum.cxx:1364
 AliHFPtSpectrum.cxx:1365
 AliHFPtSpectrum.cxx:1366
 AliHFPtSpectrum.cxx:1367
 AliHFPtSpectrum.cxx:1368
 AliHFPtSpectrum.cxx:1369
 AliHFPtSpectrum.cxx:1370
 AliHFPtSpectrum.cxx:1371
 AliHFPtSpectrum.cxx:1372
 AliHFPtSpectrum.cxx:1373
 AliHFPtSpectrum.cxx:1374
 AliHFPtSpectrum.cxx:1375
 AliHFPtSpectrum.cxx:1376
 AliHFPtSpectrum.cxx:1377
 AliHFPtSpectrum.cxx:1378
 AliHFPtSpectrum.cxx:1379
 AliHFPtSpectrum.cxx:1380
 AliHFPtSpectrum.cxx:1381
 AliHFPtSpectrum.cxx:1382
 AliHFPtSpectrum.cxx:1383
 AliHFPtSpectrum.cxx:1384
 AliHFPtSpectrum.cxx:1385
 AliHFPtSpectrum.cxx:1386
 AliHFPtSpectrum.cxx:1387
 AliHFPtSpectrum.cxx:1388
 AliHFPtSpectrum.cxx:1389
 AliHFPtSpectrum.cxx:1390
 AliHFPtSpectrum.cxx:1391
 AliHFPtSpectrum.cxx:1392
 AliHFPtSpectrum.cxx:1393
 AliHFPtSpectrum.cxx:1394
 AliHFPtSpectrum.cxx:1395
 AliHFPtSpectrum.cxx:1396
 AliHFPtSpectrum.cxx:1397
 AliHFPtSpectrum.cxx:1398
 AliHFPtSpectrum.cxx:1399
 AliHFPtSpectrum.cxx:1400
 AliHFPtSpectrum.cxx:1401
 AliHFPtSpectrum.cxx:1402
 AliHFPtSpectrum.cxx:1403
 AliHFPtSpectrum.cxx:1404
 AliHFPtSpectrum.cxx:1405
 AliHFPtSpectrum.cxx:1406
 AliHFPtSpectrum.cxx:1407
 AliHFPtSpectrum.cxx:1408
 AliHFPtSpectrum.cxx:1409
 AliHFPtSpectrum.cxx:1410
 AliHFPtSpectrum.cxx:1411
 AliHFPtSpectrum.cxx:1412
 AliHFPtSpectrum.cxx:1413
 AliHFPtSpectrum.cxx:1414
 AliHFPtSpectrum.cxx:1415
 AliHFPtSpectrum.cxx:1416
 AliHFPtSpectrum.cxx:1417
 AliHFPtSpectrum.cxx:1418
 AliHFPtSpectrum.cxx:1419
 AliHFPtSpectrum.cxx:1420
 AliHFPtSpectrum.cxx:1421
 AliHFPtSpectrum.cxx:1422
 AliHFPtSpectrum.cxx:1423
 AliHFPtSpectrum.cxx:1424
 AliHFPtSpectrum.cxx:1425
 AliHFPtSpectrum.cxx:1426
 AliHFPtSpectrum.cxx:1427
 AliHFPtSpectrum.cxx:1428
 AliHFPtSpectrum.cxx:1429
 AliHFPtSpectrum.cxx:1430
 AliHFPtSpectrum.cxx:1431
 AliHFPtSpectrum.cxx:1432
 AliHFPtSpectrum.cxx:1433
 AliHFPtSpectrum.cxx:1434
 AliHFPtSpectrum.cxx:1435
 AliHFPtSpectrum.cxx:1436
 AliHFPtSpectrum.cxx:1437
 AliHFPtSpectrum.cxx:1438
 AliHFPtSpectrum.cxx:1439
 AliHFPtSpectrum.cxx:1440
 AliHFPtSpectrum.cxx:1441
 AliHFPtSpectrum.cxx:1442
 AliHFPtSpectrum.cxx:1443
 AliHFPtSpectrum.cxx:1444
 AliHFPtSpectrum.cxx:1445
 AliHFPtSpectrum.cxx:1446
 AliHFPtSpectrum.cxx:1447
 AliHFPtSpectrum.cxx:1448
 AliHFPtSpectrum.cxx:1449
 AliHFPtSpectrum.cxx:1450
 AliHFPtSpectrum.cxx:1451
 AliHFPtSpectrum.cxx:1452
 AliHFPtSpectrum.cxx:1453
 AliHFPtSpectrum.cxx:1454
 AliHFPtSpectrum.cxx:1455
 AliHFPtSpectrum.cxx:1456
 AliHFPtSpectrum.cxx:1457
 AliHFPtSpectrum.cxx:1458
 AliHFPtSpectrum.cxx:1459
 AliHFPtSpectrum.cxx:1460
 AliHFPtSpectrum.cxx:1461
 AliHFPtSpectrum.cxx:1462
 AliHFPtSpectrum.cxx:1463
 AliHFPtSpectrum.cxx:1464
 AliHFPtSpectrum.cxx:1465
 AliHFPtSpectrum.cxx:1466
 AliHFPtSpectrum.cxx:1467
 AliHFPtSpectrum.cxx:1468
 AliHFPtSpectrum.cxx:1469
 AliHFPtSpectrum.cxx:1470
 AliHFPtSpectrum.cxx:1471
 AliHFPtSpectrum.cxx:1472
 AliHFPtSpectrum.cxx:1473
 AliHFPtSpectrum.cxx:1474
 AliHFPtSpectrum.cxx:1475
 AliHFPtSpectrum.cxx:1476
 AliHFPtSpectrum.cxx:1477
 AliHFPtSpectrum.cxx:1478
 AliHFPtSpectrum.cxx:1479
 AliHFPtSpectrum.cxx:1480
 AliHFPtSpectrum.cxx:1481
 AliHFPtSpectrum.cxx:1482
 AliHFPtSpectrum.cxx:1483
 AliHFPtSpectrum.cxx:1484
 AliHFPtSpectrum.cxx:1485
 AliHFPtSpectrum.cxx:1486
 AliHFPtSpectrum.cxx:1487
 AliHFPtSpectrum.cxx:1488
 AliHFPtSpectrum.cxx:1489
 AliHFPtSpectrum.cxx:1490
 AliHFPtSpectrum.cxx:1491
 AliHFPtSpectrum.cxx:1492
 AliHFPtSpectrum.cxx:1493
 AliHFPtSpectrum.cxx:1494
 AliHFPtSpectrum.cxx:1495
 AliHFPtSpectrum.cxx:1496
 AliHFPtSpectrum.cxx:1497
 AliHFPtSpectrum.cxx:1498
 AliHFPtSpectrum.cxx:1499
 AliHFPtSpectrum.cxx:1500
 AliHFPtSpectrum.cxx:1501
 AliHFPtSpectrum.cxx:1502
 AliHFPtSpectrum.cxx:1503
 AliHFPtSpectrum.cxx:1504
 AliHFPtSpectrum.cxx:1505
 AliHFPtSpectrum.cxx:1506
 AliHFPtSpectrum.cxx:1507
 AliHFPtSpectrum.cxx:1508
 AliHFPtSpectrum.cxx:1509
 AliHFPtSpectrum.cxx:1510
 AliHFPtSpectrum.cxx:1511
 AliHFPtSpectrum.cxx:1512
 AliHFPtSpectrum.cxx:1513
 AliHFPtSpectrum.cxx:1514
 AliHFPtSpectrum.cxx:1515
 AliHFPtSpectrum.cxx:1516
 AliHFPtSpectrum.cxx:1517
 AliHFPtSpectrum.cxx:1518
 AliHFPtSpectrum.cxx:1519
 AliHFPtSpectrum.cxx:1520
 AliHFPtSpectrum.cxx:1521
 AliHFPtSpectrum.cxx:1522
 AliHFPtSpectrum.cxx:1523
 AliHFPtSpectrum.cxx:1524
 AliHFPtSpectrum.cxx:1525
 AliHFPtSpectrum.cxx:1526
 AliHFPtSpectrum.cxx:1527
 AliHFPtSpectrum.cxx:1528
 AliHFPtSpectrum.cxx:1529
 AliHFPtSpectrum.cxx:1530
 AliHFPtSpectrum.cxx:1531
 AliHFPtSpectrum.cxx:1532
 AliHFPtSpectrum.cxx:1533
 AliHFPtSpectrum.cxx:1534
 AliHFPtSpectrum.cxx:1535
 AliHFPtSpectrum.cxx:1536
 AliHFPtSpectrum.cxx:1537
 AliHFPtSpectrum.cxx:1538
 AliHFPtSpectrum.cxx:1539
 AliHFPtSpectrum.cxx:1540
 AliHFPtSpectrum.cxx:1541
 AliHFPtSpectrum.cxx:1542
 AliHFPtSpectrum.cxx:1543
 AliHFPtSpectrum.cxx:1544
 AliHFPtSpectrum.cxx:1545
 AliHFPtSpectrum.cxx:1546
 AliHFPtSpectrum.cxx:1547
 AliHFPtSpectrum.cxx:1548
 AliHFPtSpectrum.cxx:1549
 AliHFPtSpectrum.cxx:1550
 AliHFPtSpectrum.cxx:1551
 AliHFPtSpectrum.cxx:1552
 AliHFPtSpectrum.cxx:1553
 AliHFPtSpectrum.cxx:1554
 AliHFPtSpectrum.cxx:1555
 AliHFPtSpectrum.cxx:1556
 AliHFPtSpectrum.cxx:1557
 AliHFPtSpectrum.cxx:1558
 AliHFPtSpectrum.cxx:1559
 AliHFPtSpectrum.cxx:1560
 AliHFPtSpectrum.cxx:1561
 AliHFPtSpectrum.cxx:1562
 AliHFPtSpectrum.cxx:1563
 AliHFPtSpectrum.cxx:1564
 AliHFPtSpectrum.cxx:1565
 AliHFPtSpectrum.cxx:1566
 AliHFPtSpectrum.cxx:1567
 AliHFPtSpectrum.cxx:1568
 AliHFPtSpectrum.cxx:1569
 AliHFPtSpectrum.cxx:1570
 AliHFPtSpectrum.cxx:1571
 AliHFPtSpectrum.cxx:1572
 AliHFPtSpectrum.cxx:1573
 AliHFPtSpectrum.cxx:1574
 AliHFPtSpectrum.cxx:1575
 AliHFPtSpectrum.cxx:1576
 AliHFPtSpectrum.cxx:1577
 AliHFPtSpectrum.cxx:1578
 AliHFPtSpectrum.cxx:1579
 AliHFPtSpectrum.cxx:1580
 AliHFPtSpectrum.cxx:1581
 AliHFPtSpectrum.cxx:1582
 AliHFPtSpectrum.cxx:1583
 AliHFPtSpectrum.cxx:1584
 AliHFPtSpectrum.cxx:1585
 AliHFPtSpectrum.cxx:1586
 AliHFPtSpectrum.cxx:1587
 AliHFPtSpectrum.cxx:1588
 AliHFPtSpectrum.cxx:1589
 AliHFPtSpectrum.cxx:1590
 AliHFPtSpectrum.cxx:1591
 AliHFPtSpectrum.cxx:1592
 AliHFPtSpectrum.cxx:1593
 AliHFPtSpectrum.cxx:1594
 AliHFPtSpectrum.cxx:1595
 AliHFPtSpectrum.cxx:1596
 AliHFPtSpectrum.cxx:1597
 AliHFPtSpectrum.cxx:1598
 AliHFPtSpectrum.cxx:1599
 AliHFPtSpectrum.cxx:1600
 AliHFPtSpectrum.cxx:1601
 AliHFPtSpectrum.cxx:1602
 AliHFPtSpectrum.cxx:1603
 AliHFPtSpectrum.cxx:1604
 AliHFPtSpectrum.cxx:1605
 AliHFPtSpectrum.cxx:1606
 AliHFPtSpectrum.cxx:1607
 AliHFPtSpectrum.cxx:1608
 AliHFPtSpectrum.cxx:1609
 AliHFPtSpectrum.cxx:1610
 AliHFPtSpectrum.cxx:1611
 AliHFPtSpectrum.cxx:1612
 AliHFPtSpectrum.cxx:1613
 AliHFPtSpectrum.cxx:1614
 AliHFPtSpectrum.cxx:1615
 AliHFPtSpectrum.cxx:1616
 AliHFPtSpectrum.cxx:1617
 AliHFPtSpectrum.cxx:1618
 AliHFPtSpectrum.cxx:1619
 AliHFPtSpectrum.cxx:1620
 AliHFPtSpectrum.cxx:1621
 AliHFPtSpectrum.cxx:1622
 AliHFPtSpectrum.cxx:1623
 AliHFPtSpectrum.cxx:1624
 AliHFPtSpectrum.cxx:1625
 AliHFPtSpectrum.cxx:1626
 AliHFPtSpectrum.cxx:1627
 AliHFPtSpectrum.cxx:1628
 AliHFPtSpectrum.cxx:1629
 AliHFPtSpectrum.cxx:1630
 AliHFPtSpectrum.cxx:1631
 AliHFPtSpectrum.cxx:1632
 AliHFPtSpectrum.cxx:1633
 AliHFPtSpectrum.cxx:1634
 AliHFPtSpectrum.cxx:1635
 AliHFPtSpectrum.cxx:1636
 AliHFPtSpectrum.cxx:1637
 AliHFPtSpectrum.cxx:1638
 AliHFPtSpectrum.cxx:1639
 AliHFPtSpectrum.cxx:1640
 AliHFPtSpectrum.cxx:1641
 AliHFPtSpectrum.cxx:1642
 AliHFPtSpectrum.cxx:1643
 AliHFPtSpectrum.cxx:1644
 AliHFPtSpectrum.cxx:1645
 AliHFPtSpectrum.cxx:1646
 AliHFPtSpectrum.cxx:1647
 AliHFPtSpectrum.cxx:1648
 AliHFPtSpectrum.cxx:1649
 AliHFPtSpectrum.cxx:1650
 AliHFPtSpectrum.cxx:1651
 AliHFPtSpectrum.cxx:1652
 AliHFPtSpectrum.cxx:1653
 AliHFPtSpectrum.cxx:1654
 AliHFPtSpectrum.cxx:1655
 AliHFPtSpectrum.cxx:1656
 AliHFPtSpectrum.cxx:1657
 AliHFPtSpectrum.cxx:1658
 AliHFPtSpectrum.cxx:1659
 AliHFPtSpectrum.cxx:1660
 AliHFPtSpectrum.cxx:1661
 AliHFPtSpectrum.cxx:1662
 AliHFPtSpectrum.cxx:1663
 AliHFPtSpectrum.cxx:1664
 AliHFPtSpectrum.cxx:1665
 AliHFPtSpectrum.cxx:1666
 AliHFPtSpectrum.cxx:1667
 AliHFPtSpectrum.cxx:1668
 AliHFPtSpectrum.cxx:1669
 AliHFPtSpectrum.cxx:1670
 AliHFPtSpectrum.cxx:1671
 AliHFPtSpectrum.cxx:1672
 AliHFPtSpectrum.cxx:1673
 AliHFPtSpectrum.cxx:1674
 AliHFPtSpectrum.cxx:1675
 AliHFPtSpectrum.cxx:1676
 AliHFPtSpectrum.cxx:1677
 AliHFPtSpectrum.cxx:1678
 AliHFPtSpectrum.cxx:1679
 AliHFPtSpectrum.cxx:1680
 AliHFPtSpectrum.cxx:1681
 AliHFPtSpectrum.cxx:1682
 AliHFPtSpectrum.cxx:1683
 AliHFPtSpectrum.cxx:1684
 AliHFPtSpectrum.cxx:1685
 AliHFPtSpectrum.cxx:1686
 AliHFPtSpectrum.cxx:1687
 AliHFPtSpectrum.cxx:1688
 AliHFPtSpectrum.cxx:1689
 AliHFPtSpectrum.cxx:1690
 AliHFPtSpectrum.cxx:1691
 AliHFPtSpectrum.cxx:1692
 AliHFPtSpectrum.cxx:1693
 AliHFPtSpectrum.cxx:1694
 AliHFPtSpectrum.cxx:1695
 AliHFPtSpectrum.cxx:1696
 AliHFPtSpectrum.cxx:1697
 AliHFPtSpectrum.cxx:1698
 AliHFPtSpectrum.cxx:1699
 AliHFPtSpectrum.cxx:1700
 AliHFPtSpectrum.cxx:1701
 AliHFPtSpectrum.cxx:1702
 AliHFPtSpectrum.cxx:1703
 AliHFPtSpectrum.cxx:1704
 AliHFPtSpectrum.cxx:1705
 AliHFPtSpectrum.cxx:1706
 AliHFPtSpectrum.cxx:1707
 AliHFPtSpectrum.cxx:1708
 AliHFPtSpectrum.cxx:1709
 AliHFPtSpectrum.cxx:1710
 AliHFPtSpectrum.cxx:1711
 AliHFPtSpectrum.cxx:1712
 AliHFPtSpectrum.cxx:1713
 AliHFPtSpectrum.cxx:1714
 AliHFPtSpectrum.cxx:1715
 AliHFPtSpectrum.cxx:1716
 AliHFPtSpectrum.cxx:1717
 AliHFPtSpectrum.cxx:1718
 AliHFPtSpectrum.cxx:1719
 AliHFPtSpectrum.cxx:1720
 AliHFPtSpectrum.cxx:1721
 AliHFPtSpectrum.cxx:1722
 AliHFPtSpectrum.cxx:1723
 AliHFPtSpectrum.cxx:1724
 AliHFPtSpectrum.cxx:1725
 AliHFPtSpectrum.cxx:1726
 AliHFPtSpectrum.cxx:1727
 AliHFPtSpectrum.cxx:1728
 AliHFPtSpectrum.cxx:1729
 AliHFPtSpectrum.cxx:1730
 AliHFPtSpectrum.cxx:1731
 AliHFPtSpectrum.cxx:1732
 AliHFPtSpectrum.cxx:1733
 AliHFPtSpectrum.cxx:1734
 AliHFPtSpectrum.cxx:1735
 AliHFPtSpectrum.cxx:1736
 AliHFPtSpectrum.cxx:1737
 AliHFPtSpectrum.cxx:1738
 AliHFPtSpectrum.cxx:1739
 AliHFPtSpectrum.cxx:1740
 AliHFPtSpectrum.cxx:1741
 AliHFPtSpectrum.cxx:1742
 AliHFPtSpectrum.cxx:1743
 AliHFPtSpectrum.cxx:1744
 AliHFPtSpectrum.cxx:1745
 AliHFPtSpectrum.cxx:1746
 AliHFPtSpectrum.cxx:1747
 AliHFPtSpectrum.cxx:1748
 AliHFPtSpectrum.cxx:1749
 AliHFPtSpectrum.cxx:1750
 AliHFPtSpectrum.cxx:1751
 AliHFPtSpectrum.cxx:1752
 AliHFPtSpectrum.cxx:1753
 AliHFPtSpectrum.cxx:1754
 AliHFPtSpectrum.cxx:1755
 AliHFPtSpectrum.cxx:1756
 AliHFPtSpectrum.cxx:1757
 AliHFPtSpectrum.cxx:1758
 AliHFPtSpectrum.cxx:1759
 AliHFPtSpectrum.cxx:1760
 AliHFPtSpectrum.cxx:1761
 AliHFPtSpectrum.cxx:1762
 AliHFPtSpectrum.cxx:1763
 AliHFPtSpectrum.cxx:1764
 AliHFPtSpectrum.cxx:1765
 AliHFPtSpectrum.cxx:1766
 AliHFPtSpectrum.cxx:1767
 AliHFPtSpectrum.cxx:1768
 AliHFPtSpectrum.cxx:1769
 AliHFPtSpectrum.cxx:1770
 AliHFPtSpectrum.cxx:1771
 AliHFPtSpectrum.cxx:1772
 AliHFPtSpectrum.cxx:1773
 AliHFPtSpectrum.cxx:1774
 AliHFPtSpectrum.cxx:1775
 AliHFPtSpectrum.cxx:1776
 AliHFPtSpectrum.cxx:1777
 AliHFPtSpectrum.cxx:1778
 AliHFPtSpectrum.cxx:1779
 AliHFPtSpectrum.cxx:1780
 AliHFPtSpectrum.cxx:1781
 AliHFPtSpectrum.cxx:1782
 AliHFPtSpectrum.cxx:1783
 AliHFPtSpectrum.cxx:1784
 AliHFPtSpectrum.cxx:1785
 AliHFPtSpectrum.cxx:1786
 AliHFPtSpectrum.cxx:1787
 AliHFPtSpectrum.cxx:1788
 AliHFPtSpectrum.cxx:1789
 AliHFPtSpectrum.cxx:1790
 AliHFPtSpectrum.cxx:1791
 AliHFPtSpectrum.cxx:1792
 AliHFPtSpectrum.cxx:1793
 AliHFPtSpectrum.cxx:1794
 AliHFPtSpectrum.cxx:1795
 AliHFPtSpectrum.cxx:1796