ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, 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 to generate correlated Heavy Flavor hadron pairs (one or several pairs
// per event) using paramtrized kinematics of quark pairs from some generator
// and quark fragmentation functions.
// Is a generalisation of AliGenParam class for correlated pairs of hadrons.
// In this version quark pairs and fragmentation functions are obtained from
// ~2.10^6 Pythia6.214 events generated with kCharmppMNRwmi & kBeautyppMNRwmi, 
// CTEQ5L PDF and Pt_hard = 2.76 GeV/c for p-p collisions at 2.76, 7, 8, 10 and 14 TeV,
// and with kCharmppMNR (Pt_hard = 2.10 GeV/c) & kBeautyppMNR (Pt_hard = 2.75 GeV/c), 
// CTEQ4L PDF for Pb-Pb at 2.76 and 3.94 TeV, for p-Pb & Pb-p at 5 and 8.8 TeV. 
// Decays are performed by Pythia.
// Author: S. Grigoryan, LPC Clermont-Fd & YerPhI, Smbat.Grigoryan@cern.ch
// July 07: added quarks in the stack (B. Vulpescu)
// April 09: added energy choice between 10 and 14 TeV (S. Grigoryan)
// Sept 09: added hadron pair composition probabilities via 2D histo (X.M. Zhang)
// Oct 09: added energy choice between 7, 10, 14 TeV (for p-p), 4 TeV (for Pb-Pb),
// 9 TeV (for p-Pb) and -9 TeV (for Pb-p) (S. Grigoryan)
// April 10: removed "static" from definition of some variables (B. Vulpescu)
// May 11: added Flag for transportation of background particles while using 
// SetForceDecay() function (L. Manceau)
// June 11: added modifications allowing the setting of cuts on HF-hadron children.
// Quarks, hadrons and decay particles are loaded in the stack outside the loop
// of HF-hadrons, when the cuts on their children are satisfied (L. Manceau)
// Oct 11: added Pb-Pb at 2.76 TeV (S. Grigoryan)
// June 12: added p-Pb & Pb-p at 5 TeV (S. Grigoryan)
// April 13: added p-p at 2.76 and 8 TeV (S. Grigoryan)
// 
//-------------------------------------------------------------------------
// How it works (for the given flavor and p-p energy):
//
// 1) Reads QQbar kinematical grid (TTree) from the Input file and generates
// quark pairs according to the weights of the cells.
// It is a 5D grid in y1,y2,pt1,pt2 and deltaphi, with occupancy weights
// of the cells obtained from Pythia (see details in GetQuarkPair).
// 2) Reads "soft" and "hard" fragmentation functions (12 2D-histograms each,
// for 12 pt bins) from the Input file, applies to quarks and produces hadrons
// (only lower states, with proportions of species obtained from Pythia).
// Fragmentation functions are the same for all hadron species and depend
// on 2 variables - light cone energy-momentum fractions:
//     z1=(E_H + Pz_H)/(E_Q + Pz_Q),  z2=(E_H - Pz_H)/(E_Q - Pz_Q).
// "soft" & "hard" FFs correspond to "slower" & "faster" quark of a pair 
// (see details in GetHadronPair). Fragmentation does not depend on p-p energy.
// 3) Decays the hadrons and saves all the particles in the event stack in the
// following order: HF hadron from Q, then its decay products, then HF hadron
// from Qbar, then its decay productes, then next HF hadon pair (if any) 
// in the same way, etc... 
// 4) It is fast, e.g., generates the same number of events with a beauty pair 
//  ~15 times faster than AliGenPythia with kBeautyppMNRwmi (w/o tracking)
//
// An Input file for each quark flavor and p-p energy is in EVGEN/dataCorrHF/
// One can use also user-defined Input files.
//
// More details could be found in my presentation at DiMuonNet Workshop, Dec 2006: 
// http://www-dapnia.cea.fr/Sphn/Alice/DiMuonNet.
//
//-------------------------------------------------------------------------
// How to use it:
//
// add the following typical lines in Config.C
/*
    // An example for correlated charm or beauty hadron pair production at 14 TeV

    // AliGenCorrHF *gener = new AliGenCorrHF(1, 4, 14);  // for charm, 1 pair per event
    AliGenCorrHF *gener = new AliGenCorrHF(1, 5, 14);  // for beauty, 1 pair per event

    gener->SetMomentumRange(0,9999);
    gener->SetCutOnChild(0);          // 1/0 means cuts on children enable/disable
    gener->SetChildThetaRange(171.0,178.0);
    gener->SetOrigin(0,0,0);          //vertex position    
    gener->SetSigma(0,0,0);           //Sigma in (X,Y,Z) (cm) on IP position
    gener->SetForceDecay(kSemiMuonic);
    gener->SetSelectAll(kTRUE);      //Force the transport of all particles. 
                                     //Necessary while using a different 
                                     //option than kAll for SetForceDecay
    gener->SetTrackingFlag(1);       //1: Decay during transport, 
                                     //0: No Decay during transport
    gener->Init();
*/
// One can include AliGenCorrHF in an AliGenCocktail generator.
//--------------------------------------------------------------------------

#include <Riostream.h>
#include <TCanvas.h>
#include <TClonesArray.h>
#include <TDatabasePDG.h>
#include <TFile.h>
#include <TH2F.h>
#include <TLorentzVector.h>
#include <TMath.h>
#include <TParticle.h>
#include <TParticlePDG.h>
#include <TROOT.h>
#include <TRandom.h>
#include <TTree.h>
#include <TVirtualMC.h>
#include <TVector3.h>

#include "AliGenCorrHF.h"
#include "AliLog.h"
#include "AliConst.h"
#include "AliDecayer.h"
#include "AliMC.h"
#include "AliRun.h"
#include "AliGenEventHeader.h"

ClassImp(AliGenCorrHF)

  //Begin_Html
  /*
    <img src="picts/AliGenCorrHF.gif">
  */
  //End_Html

Double_t AliGenCorrHF::fgdph[19] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100, 110, 120, 130, 140, 150, 160, 170, 180};
Double_t AliGenCorrHF::fgy[31] = {-10,-7, -6.5, -6, -5.5, -5, -4.5, -4, -3.5, -3, -2.5, -2,- 1.5, -1, -0.5, 0, 0.5, 1, 1.5, 2, 2.5, 3, 3.5, 4, 4.5, 5, 5.5, 6, 6.5, 7, 10};
Double_t AliGenCorrHF::fgpt[51] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 3.5, 4, 4.5, 5, 5.5, 6, 6.6, 7.2, 7.8, 8.4, 9, 9.6, 10.3, 11.1, 12, 13, 14, 15, 16, 17, 18, 19, 20.1, 21.5, 23, 24.5, 26, 27.5, 29.1, 31, 33, 35, 37, 39.2, 42, 45, 48, 51, 55.2, 60, 65, 71, 81, 100};
Int_t AliGenCorrHF::fgnptbins = 12;
Double_t AliGenCorrHF::fgptbmin[12] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9};
Double_t AliGenCorrHF::fgptbmax[12] = {0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9, 100};

//____________________________________________________________
    AliGenCorrHF::AliGenCorrHF():
	fFileName(0),
	fFile(0),
	fQuark(0),
	fEnergy(0),
	fBias(0.),
	fTrials(0),
	fSelectAll(kFALSE),
	fDecayer(0),
	fgIntegral(0)
{
// Default constructor
}

//____________________________________________________________
AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t idquark, Int_t energy):
    AliGenMC(npart),
    fFileName(0),
    fFile(0),
    fQuark(idquark),
    fEnergy(energy),
    fBias(0.),
    fTrials(0),
    fSelectAll(kFALSE),
    fDecayer(0),
    fgIntegral(0)
{
// Constructor using particle number, quark type, energy & default InputFile
//
    if (fQuark == 5) {
      if (fEnergy == 7)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP7PythiaMNRwmi.root";
      else if (fEnergy == 8)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP8PythiaMNRwmi.root";
      else if (fEnergy == 10)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP10PythiaMNRwmi.root";
      else if (fEnergy == 14)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP14PythiaMNRwmi.root";
      else if (fEnergy == 2)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP276PythiaMNRwmi.root";
      else if (fEnergy == 3)
	   fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb276PythiaMNR.root";
      else if (fEnergy == 4)
	   fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
      else if (fEnergy == 5 || fEnergy == -5)
	   fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb5PythiaMNR.root";
      else if (fEnergy == 9 || fEnergy == -9)
	   fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb88PythiaMNR.root";
      else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
    }
    else {
      fQuark = 4;
      if (fEnergy == 7)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP7PythiaMNRwmi.root";
      else if (fEnergy == 8)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP8PythiaMNRwmi.root";
      else if (fEnergy == 10)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP10PythiaMNRwmi.root";
      else if (fEnergy == 14)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP14PythiaMNRwmi.root";
      else if (fEnergy == 2)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP276PythiaMNRwmi.root";
      else if (fEnergy == 3)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb276PythiaMNR.root";
      else if (fEnergy == 4)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
      else if (fEnergy == 5 || fEnergy == -5)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb5PythiaMNR.root";
      else if (fEnergy == 9 || fEnergy == -9)
           fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb88PythiaMNR.root";
      else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
    }
    fName = "Default";
    fTitle= "Generator for correlated pairs of HF hadrons";
      
    fChildSelect.Set(5);
    for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
    SetForceDecay();
    SetCutOnChild();
    SetChildMomentumRange();
    SetChildPtRange();
    SetChildPhiRange();
    SetChildThetaRange(); 
}

//___________________________________________________________________
AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t idquark, Int_t energy):
    AliGenMC(npart),
    fFileName(tname),
    fFile(0),
    fQuark(idquark),
    fEnergy(energy),
    fBias(0.),
    fTrials(0),
    fSelectAll(kFALSE),
    fDecayer(0),
    fgIntegral(0)
{
// Constructor using particle number, quark type, energy & user-defined InputFile
//
    if (fQuark != 5) fQuark = 4;
    fName = "UserDefined";
    fTitle= "Generator for correlated pairs of HF hadrons";
      
    fChildSelect.Set(5);
    for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
    SetForceDecay();
    SetCutOnChild();
    SetChildMomentumRange();
    SetChildPtRange();
    SetChildPhiRange();
    SetChildThetaRange(); 
}

//____________________________________________________________
AliGenCorrHF::~AliGenCorrHF()
{
// Destructor
  delete fFile;
}

//____________________________________________________________
void AliGenCorrHF::Init()
{
// Initialisation
  AliInfo(Form("Number of HF-hadron pairs = %d",fNpart)); 
  AliInfo(Form(" QQbar kinematics and fragm. functions from:  %s",fFileName.Data() )); 
    fFile = TFile::Open(fFileName.Data());
    if(!fFile->IsOpen()){
      AliError(Form("Could not open file %s",fFileName.Data() ));
    }

    ComputeIntegral(fFile);
    
    fParentWeight = 1./fNpart;   // fNpart is number of HF-hadron pairs

// particle decay related initialization

    if (TVirtualMC::GetMC()) fDecayer = TVirtualMC::GetMC()->GetDecayer();
    fDecayer->SetForceDecay(fForceDecay);
    fDecayer->Init();

//
    AliGenMC::Init();
}
//____________________________________________________________
void AliGenCorrHF::Generate()
{
//
// Generate fNpart of correlated HF hadron pairs per event
// in the the desired theta and momentum windows (phi = 0 - 2pi).
//

//  Reinitialize decayer

  fDecayer->SetForceDecay(fForceDecay);
  fDecayer->Init();

  Float_t polar[2][3];        // Polarisation of the parent particle (for GEANT tracking)
  Float_t origin0[2][3];      // Origin of the generated parent particle (for GEANT tracking)
  Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
  Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
  Float_t p[2][3];            // Momenta
  Int_t nt, i, j, ihad, ipa, ipa0, ipa1, ihadron[2], iquark[2];
  Float_t  wgtp[2], wgtch[2], random[6];
  Float_t pq[2][3], pc[3];    // Momenta of the two quarks
  Double_t tanhy2, qm = 0;
  Int_t np[2];
  Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2], phiq[2];
  Int_t ncsel[2];
  Int_t** pSelected = new Int_t* [2];
  Int_t** trackIt = new Int_t* [2];

  for (i=0; i<2; i++) { 
    ptq[i]     =0; 
    yq[i]      =0; 
    pth[i]     =0; 
    plh[i]     =0;
    phih[i]    =0; 
    phiq[i]    =0;
    ihadron[i] =0; 
    iquark[i]  =0;
    for (j=0; j<3; j++) polar[i][j]=0;
  }

  // same quarks mass as in the fragmentation functions
  if (fQuark == 4) qm = 1.20;
  else             qm = 4.75;
  
  TClonesArray *particleshad1 = new TClonesArray("TParticle",1000);
  TClonesArray *particleshad2 = new TClonesArray("TParticle",1000);
  
  TList *particleslist = new TList();
  particleslist->Add(particleshad1);
  particleslist->Add(particleshad2);
  
  TDatabasePDG *pDataBase = TDatabasePDG::Instance();

  // Calculating vertex position per event
  if (fVertexSmear==kPerEvent) {
    Vertex();
    for (i=0;i<2;i++){
      for (j=0;j<3;j++) origin0[i][j]=fVertex[j];
    }
  }
  else {
    for (i=0;i<2;i++){
      for (j=0;j<3;j++) origin0[i][j]=fOrigin[j];
    }
  }
  
  ipa  = 0;
  ipa1 = 0;
  ipa0 = 0;
  
  // Generating fNpart HF-hadron pairs
  fNprimaries = 0;
 
  while (ipa<2*fNpart) {

    GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
    
    GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
    
    // Boost particles from c.m.s. to ALICE lab frame for p-Pb & Pb-p collisions
    if (fEnergy == 5 || fEnergy == -5 || fEnergy == 9 || fEnergy == -9) {
      Double_t dyBoost = 0.47;
      Double_t beta  = TMath::TanH(dyBoost);
      Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
      Double_t gb    = gamma * beta;
      yq[0] += dyBoost;
      yq[1] += dyBoost;
      plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
      plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
      if (fEnergy == 5 || fEnergy == 9) {
	yq[0] *= -1;
	yq[1] *= -1;
	plh[0] *= -1;
	plh[1] *= -1;
      }
    }      
    
    // Cuts from AliGenerator
    
    // Cut on theta
    theta=TMath::ATan2(pth[0],plh[0]);
    if (theta<fThetaMin || theta>fThetaMax) continue;
    theta=TMath::ATan2(pth[1],plh[1]);
    if (theta<fThetaMin || theta>fThetaMax) continue;
    
    // Cut on momentum
    ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
    if (ph[0]<fPMin || ph[0]>fPMax) continue;
    ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
    if (ph[1]<fPMin || ph[1]>fPMax) continue;
    
    // Add the quarks in the stack
    
    phiq[0] = Rndm()*k2PI;
    if (Rndm() < 0.5) {
      phiq[1] = phiq[0] + dphi*kDegrad; 
    } else {
      phiq[1] = phiq[0] - dphi*kDegrad; 
    }    
    if (phiq[1] > k2PI) phiq[1] -= k2PI;
    if (phiq[1] < 0   ) phiq[1] += k2PI;
    
    // quarks pdg
    iquark[0] = +fQuark;
    iquark[1] = -fQuark;
    
    // px and py
    TVector2 qvect1 = TVector2();
    TVector2 qvect2 = TVector2();
    qvect1.SetMagPhi(ptq[0],phiq[0]);
    qvect2.SetMagPhi(ptq[1],phiq[1]);
    pq[0][0] = qvect1.Px();
    pq[0][1] = qvect1.Py();
    pq[1][0] = qvect2.Px();
    pq[1][1] = qvect2.Py();
    
    // pz
    tanhy2 = TMath::TanH(yq[0]);
    tanhy2 *= tanhy2;
    pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
    pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
    tanhy2 = TMath::TanH(yq[1]);
    tanhy2 *= tanhy2;
    pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
    pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
    
    // Here we assume that  |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
    // which is a good approximation for heavy flavors in Pythia
    // ... moreover, same phi angles as for the quarks ...
    
    phih[0] = phiq[0];    
    phih[1] = phiq[1];    
    
    ipa1 = 0;

    for (ihad = 0; ihad < 2; ihad++) {
      while(1) {
	
	ipa0=ipa1;
	
	// particle type 
	fChildWeight=(fDecayer->GetPartialBranchingRatio(ihadron[ihad]))*fParentWeight;
	wgtp[ihad]=fParentWeight;
	wgtch[ihad]=fChildWeight;
	TParticlePDG *particle = pDataBase->GetParticle(ihadron[ihad]);
	Float_t am = particle->Mass();
	phi = phih[ihad];
	pt  = pth[ihad];
	pl  = plh[ihad];
	ptot=TMath::Sqrt(pt*pt+pl*pl);
	
	p[ihad][0]=pt*TMath::Cos(phi);
	p[ihad][1]=pt*TMath::Sin(phi);
	p[ihad][2]=pl;
	
	if(fVertexSmear==kPerTrack) {
	  Rndm(random,6);
	  for (j=0;j<3;j++) {
	    origin0[ihad][j]=
	      fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
	      TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
	  }
	}
	
	// Looking at fForceDecay : 
	// if fForceDecay != none Primary particle decays using 
	// AliPythia and children are tracked by GEANT
	//
	// if fForceDecay == none Primary particle is tracked by GEANT 
	// (In the latest, make sure that GEANT actually does all the decays you want)	  

	if (fForceDecay != kNoDecay) {
	  // Using lujet to decay particle
	  Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
	  TLorentzVector pmom(p[ihad][0], p[ihad][1], p[ihad][2], energy);
	  fDecayer->Decay(ihadron[ihad],&pmom);
	  
	  // select decay particles
	 
	  np[ihad]=fDecayer->ImportParticles((TClonesArray *)particleslist->At(ihad));

	  //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
	  
	  if (fGeometryAcceptance) 
	    if (!CheckAcceptanceGeometry(np[ihad],(TClonesArray*)particleslist->At(ihad))) continue;
	  
	  trackIt[ihad]     = new Int_t [np[ihad]];
	  pSelected[ihad]   = new Int_t [np[ihad]];
	  Int_t* pFlag      = new Int_t [np[ihad]];
	  
	  for (i=0; i<np[ihad]; i++) {
	    pFlag[i]     =  0;
	    pSelected[ihad][i] =  0;
	  }

	  if (np[ihad] >1) {
	    TParticle* iparticle = 0;
	    Int_t ipF, ipL;
	    
	    for (i = 1; i<np[ihad] ; i++) {
	      trackIt[ihad][i] = 1;
	      iparticle = 
		(TParticle *) ((TClonesArray *) particleslist->At(ihad))->At(i);
	      Int_t kf = iparticle->GetPdgCode();
	      Int_t ks = iparticle->GetStatusCode();
	      // flagged particle
	      if (pFlag[i] == 1) {
		ipF = iparticle->GetFirstDaughter();
		ipL = iparticle->GetLastDaughter();	
		if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
		continue;
	      }
	      
	      // flag decay products of particles with long life-time (c tau > .3 mum)
	      if (ks != 1) { 
		Double_t lifeTime = fDecayer->GetLifetime(kf);
		if (lifeTime > (Double_t) fMaxLifeTime) {
		  ipF = iparticle->GetFirstDaughter();
		  ipL = iparticle->GetLastDaughter();	
		  if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
		} else {
		  trackIt[ihad][i]     = 0;
		  pSelected[ihad][i]   = 1;
		}
	      } // ks==1 ?
	      //
	      // children
	      if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[ihad][i])
		{      
		  if (fCutOnChild) {
		    pc[0]=iparticle->Px();
		    pc[1]=iparticle->Py();
		    pc[2]=iparticle->Pz();
		    //printf("px %f py %f pz %f\n",pc[0],pc[1],pc[2]);
		    Bool_t  childok = KinematicSelection(iparticle, 1);
		    if(childok) {
		      pSelected[ihad][i]  = 1;
		      ncsel[ihad]++;
		    } else {
		      ncsel[ihad]=-1;
		      break;
		    } // child kine cuts
		  } else {
		    pSelected[ihad][i]  = 1;
		    ncsel[ihad]++;
		  } // if child selection
		} // select muon
	    } // decay particle loop
	  } // if decay products

	  if ((fCutOnChild && ncsel[ihad] >0) || !fCutOnChild) ipa1++;

	  if (pFlag) delete[] pFlag;

	} // kinematic selection
	else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
	  {
	    gAlice->GetMCApp()->
	      PushTrack(fTrackIt,-1,ihadron[ihad],p[ihad],origin0[ihad],polar[ihad],0,kPPrimary,nt,wgtp[ihad]);
	    ipa1++; 
	    fNprimaries++;
	    
	  }
	break;
      } // while(1) loop
      if (ipa1<ipa0+1){
	ipa1=0; 
	if (pSelected[ihad]) delete pSelected[ihad];
	if (trackIt[ihad])   delete trackIt[ihad];
	particleshad1->Clear();
	particleshad2->Clear();
	break;
      }//go out of loop and generate new pair if at least one hadron is rejected
    } // hadron pair loop
    if(ipa1==2){ 
 
      ipa=ipa+ipa1;
 
      if(fForceDecay != kNoDecay){
	for(ihad=0;ihad<2;ihad++){

	  //load tracks in the stack if both hadrons of the pair accepted
	  LoadTracks(iquark[ihad],pq[ihad],ihadron[ihad],p[ihad],np[ihad],
		     (TClonesArray *)particleslist->At(ihad),origin0[ihad],
		     polar[ihad],wgtp[ihad],wgtch[ihad],nt,ncsel[ihad],
		     pSelected[ihad],trackIt[ihad]);

	  if (pSelected[ihad]) delete pSelected[ihad];
	  if (trackIt[ihad])   delete trackIt[ihad];

	}
	  particleshad1->Clear();
	  particleshad2->Clear();
      }
    }
  }   // while (ipa<2*fNpart) loop
  
  SetHighWaterMark(nt);
  
  AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
  header->SetPrimaryVertex(fVertex);
  header->SetNProduced(fNprimaries);
  AddHeader(header);
  
  
  delete particleshad1;
  delete particleshad2;
  delete particleslist;
 
  delete[] pSelected;
  delete[] trackIt;
}
//____________________________________________________________________________________
void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
{  
// Composition of a lower state charm hadron pair from a ccbar quark pair
   Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};

   Double_t id3, id4;
   hProbHH->GetRandom2(id3, id4);
   pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
   pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];

   return;
}

void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
{  
// Composition of a lower state beauty hadron pair from a bbbar quark pair
   // B-Bbar mixing will be done by Pythia at their decay point
   Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};

   Double_t id3, id4;
   hProbHH->GetRandom2(id3, id4);
   pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
   pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];

   if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
   if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;

   return;
}

//____________________________________________________________________________________
Double_t AliGenCorrHF::ComputeIntegral(TFile* fG)       // needed by GetQuarkPair
{
   // Read QQbar kinematical 5D grid's cell occupancy weights
   Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
   TTree* tG = (TTree*) fG->Get("tGqq");
   tG->GetBranch("cell")->SetAddress(&cell);
   Int_t nbins = tG->GetEntries();

   //   delete previously computed integral (if any)
   if(fgIntegral) delete [] fgIntegral;

   fgIntegral = new Double_t[nbins+1];
   fgIntegral[0] = 0;
   Int_t bin;
   for(bin=0;bin<nbins;bin++) {
     tG->GetEvent(bin);
     fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
   }
   //   Normalize integral to 1
   if (fgIntegral[nbins] == 0 ) {
      return 0;
   }
   for (bin=1;bin<=nbins;bin++)  fgIntegral[bin] /= fgIntegral[nbins];

   return fgIntegral[nbins];
}


//____________________________________________________________________________________
void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)              
                                 // modification of ROOT's TH3::GetRandom3 for 5D
{
   // Read QQbar kinematical 5D grid's cell coordinates
   Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
   TTree* tG = (TTree*) fG->Get("tGqq");
   tG->GetBranch("cell")->SetAddress(&cell);
   Int_t nbins = tG->GetEntries();
   Double_t rand[6];
   gRandom->RndmArray(6,rand);
   Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
   tG->GetEvent(ibin);
   y1   = fgy[cell[1]]  + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
   y2   = fgy[cell[2]]  + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
   pt1  = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
   pt2  = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
   dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
}

//____________________________________________________________________________________
void AliGenCorrHF::GetHadronPair(TFile* fG, Int_t idq, Double_t y1, Double_t y2, Double_t pt1, Double_t pt2, Int_t &id3, Int_t &id4, Double_t &pz3, Double_t &pz4, Double_t &pt3, Double_t &pt4) 
{
    // Generate a hadron pair
    void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
    fIpParaFunc = IpCharm;
    Double_t mq = 1.2;              // c & b quark masses (used in AliPythia)
    if (idq == 5) {
      fIpParaFunc = IpBeauty;
      mq = 4.75;
    }
    Double_t z11 = 0.;
    Double_t z12 = 0.;
    Double_t z21 = 0.;
    Double_t z22 = 0.;
    Double_t pz1, pz2, e1, e2, mh, ptemp, rand[2];
    char tag[100]; 
    TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
    for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
      snprintf(tag,100, "h2h_pt%d",ipt); 
      h2h[ipt] = (TH2F*) fG->Get(tag); 
      snprintf(tag,100, "h2s_pt%d",ipt); 
      h2s[ipt] = (TH2F*) fG->Get(tag); 
    }

       if (y1*y2 < 0) {
	 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
	   if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
      	     h2h[ipt]->GetRandom2(z11, z21);
	   if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
      	     h2h[ipt]->GetRandom2(z12, z22); 
	 }
       }
       else {
	 if (TMath::Abs(y1) > TMath::Abs(y2)) {
	   for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
	     if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
	       h2h[ipt]->GetRandom2(z11, z21);
	     if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
	       h2s[ipt]->GetRandom2(z12, z22); 
	   }
	 }
	 else {
	   for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
	     if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
	       h2s[ipt]->GetRandom2(z11, z21);
	     if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
	       h2h[ipt]->GetRandom2(z12, z22); 
	   }
	 }
       }
      gRandom->RndmArray(2,rand);
      ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
      pz1   = ptemp*TMath::SinH(y1); 
      e1    = ptemp*TMath::CosH(y1); 
      ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
      pz2   = ptemp*TMath::SinH(y2); 
      e2    = ptemp*TMath::CosH(y2); 

      hProbHH = (TH2F*)fG->Get("hProbHH");
      fIpParaFunc(hProbHH, id3, id4);
      mh    = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
      ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
      if (idq==5) pt3   = pt1;                // an approximation at low pt, try better
      else        pt3   = rand[0];            // pt3=pt1 gives less D-hadrons at low pt 
      if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
      if (pz1 > 0)   pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
      else           pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
      e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);

      mh    = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
      ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
      if (idq==5) pt4   = pt2;                // an approximation at low pt, try better
      else        pt4   = rand[1];
      if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
      if (pz2 > 0)   pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
      else           pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
      e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);

      // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
      Float_t ycorr = 0.2, y3, y4;
      gRandom->RndmArray(2,rand);
      y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
      y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
      if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
	ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
	y3  = 4*(1 - 2*rand[1]);
	pz3 = ptemp*TMath::SinH(y3);
	pz4 = pz3;
      }
}

		
//____________________________________________________________________________________
void AliGenCorrHF::LoadTracks(Int_t iquark, Float_t *pq, 
			      Int_t iPart, Float_t *p, 
			      Int_t np, TClonesArray *particles,
			      Float_t *origin0, Float_t *polar, 
			      Float_t wgtp, Float_t wgtch,
			      Int_t &nt, Int_t ncsel, Int_t *pSelected, 
			      Int_t *trackIt){
  Int_t i; 
  Int_t ntq=-1;
  Int_t* pParent = new Int_t[np];
  Float_t pc[3], och[3];
  Int_t iparent;

  for(i=0;i<np;i++) pParent[i]=-1;

  if ((fCutOnChild && ncsel >0) || !fCutOnChild){  
    // Parents
    // quark
    PushTrack(0, -1, iquark, pq, origin0, polar, 0, kPPrimary, nt, wgtp);
    KeepTrack(nt);
    ntq = nt;
    // hadron
    PushTrack(0, ntq, iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
    pParent[0] = nt;
    KeepTrack(nt); 
    fNprimaries++;

    // Decay Products  
    for (i = 1; i < np; i++) {
      if (pSelected[i]) {

	TParticle* iparticle = (TParticle *) particles->At(i);
	Int_t kf  = iparticle->GetPdgCode();
	Int_t jpa = iparticle->GetFirstMother()-1;
        Int_t ksc  = iparticle->GetStatusCode();
	// RS: note, the conversion mm->cm is done now in the decayer. The time is ignored here!
	och[0] = origin0[0]+iparticle->Vx();
	och[1] = origin0[1]+iparticle->Vy();
	och[2] = origin0[2]+iparticle->Vz();
	pc[0]  = iparticle->Px();
	pc[1]  = iparticle->Py();
	pc[2]  = iparticle->Pz();
	
	if (jpa > -1) {
	  iparent = pParent[jpa];
	} else {
	  iparent = -1;
	}
	
	PushTrack(fTrackIt*trackIt[i], iparent, kf,
		  pc, och, polar,
		  0, kPDecay, nt, wgtch,ksc);
	pParent[i] = nt;
	KeepTrack(nt); 
	fNprimaries++;

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