]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenCorrHF.cxx
Possibility to introduce centrality dependence into the existing AliGenMUONCocktailpp.
[u/mrichter/AliRoot.git] / EVGEN / AliGenCorrHF.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // Class to generate correlated Heavy Flavor hadron pairs (one or several pairs
19 // per event) using paramtrized kinematics of quark pairs from some generator
20 // and quark fragmentation functions.
21 // Is a generalisation of AliGenParam class for correlated pairs of hadrons.
22 // In this version quark pairs and fragmentation functions are obtained from
23 // ~2.10^6 Pythia6.214 events generated with kCharmppMNRwmi & kBeautyppMNRwmi, 
24 // CTEQ5L PDF and Pt_hard = 2.76 GeV/c for p-p collisions at 7, 10 and 14 TeV,
25 // and with kCharmppMNR (Pt_hard = 2.10 GeV/c) & kBeautyppMNR (Pt_hard = 2.75 GeV/c), 
26 // CTEQ4L PDF for Pb-Pb at 3.94 TeV, for p-Pb & Pb-p at 8.8 TeV. 
27 // Decays are performed by Pythia.
28 // Author: S. Grigoryan, LPC Clermont-Fd & YerPhI, Smbat.Grigoryan@cern.ch
29 // July 07: added quarks in the stack (B. Vulpescu)
30 // April 09: added energy choice between 10 and 14 TeV (S. Grigoryan)
31 // Sept 09: added hadron pair composition probabilities via 2D histo (X.M. Zhang)
32 // Oct 09: added energy choice between 7, 10, 14 TeV (for p-p), 4 TeV (for Pb-Pb),
33 // 9 TeV (for p-Pb) and -9 TeV (for Pb-p) (S. Grigoryan)
34 // April 10: removed "static" from definition of some variables (B. Vulpescu)
35 //-------------------------------------------------------------------------
36 // How it works (for the given flavor and p-p energy):
37 //
38 // 1) Reads QQbar kinematical grid (TTree) from the Input file and generates
39 // quark pairs according to the weights of the cells.
40 // It is a 5D grid in y1,y2,pt1,pt2 and deltaphi, with occupancy weights
41 // of the cells obtained from Pythia (see details in GetQuarkPair).
42 // 2) Reads "soft" and "hard" fragmentation functions (12 2D-histograms each,
43 // for 12 pt bins) from the Input file, applies to quarks and produces hadrons
44 // (only lower states, with proportions of species obtained from Pythia).
45 // Fragmentation functions are the same for all hadron species and depend
46 // on 2 variables - light cone energy-momentum fractions:
47 //     z1=(E_H + Pz_H)/(E_Q + Pz_Q),  z2=(E_H - Pz_H)/(E_Q - Pz_Q).
48 // "soft" & "hard" FFs correspond to "slower" & "faster" quark of a pair 
49 // (see details in GetHadronPair). Fragmentation does not depend on p-p energy.
50 // 3) Decays the hadrons and saves all the particles in the event stack in the
51 // following order: HF hadron from Q, then its decay products, then HF hadron
52 // from Qbar, then its decay productes, then next HF hadon pair (if any) 
53 // in the same way, etc... 
54 // 4) It is fast, e.g., generates the same number of events with a beauty pair 
55 //  ~15 times faster than AliGenPythia with kBeautyppMNRwmi (w/o tracking)
56 //
57 // An Input file for each quark flavor and p-p energy is in EVGEN/dataCorrHF/
58 // One can use also user-defined Input files.
59 //
60 // More details could be found in my presentation at DiMuonNet Workshop, Dec 2006: 
61 // http://www-dapnia.cea.fr/Sphn/Alice/DiMuonNet.
62 //
63 //-------------------------------------------------------------------------
64 // How to use it:
65 //
66 // add the following typical lines in Config.C
67 /*
68   if (!strcmp(option,"corr")) {
69     // An example for correlated charm or beauty hadron pair production at 14 TeV
70
71     // AliGenCorrHF *gener = new AliGenCorrHF(1, 4, 14);  // for charm, 1 pair per event
72     AliGenCorrHF *gener = new AliGenCorrHF(1, 5, 14);  // for beauty, 1 pair per event
73
74     gener->SetMomentumRange(0,9999);
75     gener->SetCutOnChild(0);          // 1/0 means cuts on children enable/disable
76     gener->SetChildThetaRange(171.0,178.0);
77     gener->SetOrigin(0,0,0);          //vertex position    
78     gener->SetSigma(0,0,0);           //Sigma in (X,Y,Z) (cm) on IP position
79     gener->SetForceDecay(kSemiMuonic);
80     gener->SetTrackingFlag(0);
81     gener->Init();
82 }
83 */
84 // and in aliroot do e.g. gAlice->Run(10,"Config.C") to produce 10 events.
85 // One can include AliGenCorrHF in an AliGenCocktail generator.
86 //--------------------------------------------------------------------------
87
88 #include <Riostream.h>
89 #include <TCanvas.h>
90 #include <TClonesArray.h>
91 #include <TDatabasePDG.h>
92 #include <TFile.h>
93 #include <TH2F.h>
94 #include <TLorentzVector.h>
95 #include <TMath.h>
96 #include <TParticle.h>
97 #include <TParticlePDG.h>
98 #include <TROOT.h>
99 #include <TRandom.h>
100 #include <TTree.h>
101 #include <TVirtualMC.h>
102 #include <TVector3.h>
103
104 #include "AliGenCorrHF.h"
105 #include "AliLog.h"
106 #include "AliConst.h"
107 #include "AliDecayer.h"
108 #include "AliMC.h"
109 #include "AliRun.h"
110 #include "AliGenEventHeader.h"
111
112 ClassImp(AliGenCorrHF)
113
114   //Begin_Html
115   /*
116     <img src="picts/AliGenCorrHF.gif">
117   */
118   //End_Html
119
120 Double_t AliGenCorrHF::fgdph[19] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100, 110, 120, 130, 140, 150, 160, 170, 180};
121 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};
122 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};
123 Int_t AliGenCorrHF::fgnptbins = 12;
124 Double_t AliGenCorrHF::fgptbmin[12] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9};
125 Double_t AliGenCorrHF::fgptbmax[12] = {0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9, 100};
126
127 //____________________________________________________________
128     AliGenCorrHF::AliGenCorrHF():
129         fFileName(0),
130         fFile(0),
131         fQuark(0),
132         fEnergy(0),
133         fBias(0.),
134         fTrials(0),
135         fSelectAll(kFALSE),
136         fDecayer(0),
137         fgIntegral(0)
138 {
139 // Default constructor
140 }
141
142 //____________________________________________________________
143 AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t idquark, Int_t energy):
144     AliGenMC(npart),
145     fFileName(0),
146     fFile(0),
147     fQuark(idquark),
148     fEnergy(energy),
149     fBias(0.),
150     fTrials(0),
151     fSelectAll(kFALSE),
152     fDecayer(0),
153     fgIntegral(0)
154 {
155 // Constructor using particle number, quark type, energy & default InputFile
156 //
157     if (fQuark == 5) {
158       if (fEnergy == 7)
159            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP7PythiaMNRwmi.root";
160       else if (fEnergy == 10)
161            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP10PythiaMNRwmi.root";
162       else if (fEnergy == 14)
163            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP14PythiaMNRwmi.root";
164       else if (fEnergy == 4)
165            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
166       else if (fEnergy == 9 || fEnergy == -9)
167            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb88PythiaMNR.root";
168       else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
169     }
170     else {
171       fQuark = 4;
172       if (fEnergy == 7)
173            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP7PythiaMNRwmi.root";
174       else if (fEnergy == 10)
175            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP10PythiaMNRwmi.root";
176       else if (fEnergy == 14)
177            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP14PythiaMNRwmi.root";
178       else if (fEnergy == 4)
179            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
180       else if (fEnergy == 9 || fEnergy == -9)
181            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb88PythiaMNR.root";
182       else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
183     }
184     fName = "Default";
185     fTitle= "Generator for correlated pairs of HF hadrons";
186       
187     fChildSelect.Set(5);
188     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
189     SetForceDecay();
190     SetCutOnChild();
191     SetChildMomentumRange();
192     SetChildPtRange();
193     SetChildPhiRange();
194     SetChildThetaRange(); 
195 }
196
197 //___________________________________________________________________
198 AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t idquark, Int_t energy):
199     AliGenMC(npart),
200     fFileName(tname),
201     fFile(0),
202     fQuark(idquark),
203     fEnergy(energy),
204     fBias(0.),
205     fTrials(0),
206     fSelectAll(kFALSE),
207     fDecayer(0),
208     fgIntegral(0)
209 {
210 // Constructor using particle number, quark type, energy & user-defined InputFile
211 //
212     if (fQuark != 5) fQuark = 4;
213     fName = "UserDefined";
214     fTitle= "Generator for correlated pairs of HF hadrons";
215       
216     fChildSelect.Set(5);
217     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
218     SetForceDecay();
219     SetCutOnChild();
220     SetChildMomentumRange();
221     SetChildPtRange();
222     SetChildPhiRange();
223     SetChildThetaRange(); 
224 }
225
226 //____________________________________________________________
227 AliGenCorrHF::~AliGenCorrHF()
228 {
229 // Destructor
230   delete fFile;
231 }
232
233 //____________________________________________________________
234 void AliGenCorrHF::Init()
235 {
236 // Initialisation
237
238   AliInfo(Form(" QQbar kinematics and fragm. functions from:  %s",fFileName.Data() )); 
239     fFile = TFile::Open(fFileName.Data());
240     if(!fFile->IsOpen()){
241       AliError(Form("Could not open file %s",fFileName.Data() ));
242     }
243
244     ComputeIntegral(fFile);
245     
246     fParentWeight = 1./fNpart;   // fNpart is number of HF-hadron pairs
247
248 // particle decay related initialization
249
250     if (gMC) fDecayer = gMC->GetDecayer();
251     fDecayer->SetForceDecay(fForceDecay);
252     fDecayer->Init();
253
254 //
255     AliGenMC::Init();
256 }
257 //____________________________________________________________
258 void AliGenCorrHF::Generate()
259 {
260 //
261 // Generate fNpart of correlated HF hadron pairs per event
262 // in the the desired theta and momentum windows (phi = 0 - 2pi).
263 //
264
265 //  Reinitialize decayer
266
267   fDecayer->SetForceDecay(fForceDecay);
268   fDecayer->Init();
269   
270   //
271   Float_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
272   Float_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
273   Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
274   Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
275   Float_t p[3], pc[3], och[3];// Momentum, polarisation and origin of the children particles from lujet
276   Int_t nt, i, j, ipa, ihadron[2], iquark[2];
277   Float_t  wgtp, wgtch, random[6];
278   Float_t pq[2][3];           // Momenta of the two quarks
279   Int_t ntq[2] = {-1, -1};
280   Double_t tanhy2, qm = 0;
281
282   Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2], phiq[2];
283   for (i=0; i<2; i++) { 
284     ptq[i]     =0; 
285     yq[i]      =0; 
286     pth[i]     =0; 
287     plh[i]     =0;
288     phih[i]    =0; 
289     phiq[i]    =0;
290     ihadron[i] =0; 
291     iquark[i]  =0;
292   }
293
294   // same quarks mass as in the fragmentation functions
295   if (fQuark == 4) qm = 1.20;
296   else             qm = 4.75;
297
298   TClonesArray *particles = new TClonesArray("TParticle",1000);
299   
300   TDatabasePDG *pDataBase = TDatabasePDG::Instance();
301   //
302   
303   // Calculating vertex position per event
304   for (j=0;j<3;j++) origin0[j]=fOrigin[j];
305   if (fVertexSmear==kPerEvent) {
306     Vertex();
307     for (j=0;j<3;j++) origin0[j]=fVertex[j];
308   }
309   
310   ipa=0;
311   
312   // Generating fNpart HF-hadron pairs
313   fNprimaries = 0;
314   
315   while (ipa<2*fNpart) {
316     
317     GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
318     
319     GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
320     
321     if (fEnergy == 9 || fEnergy == -9) {      // boost particles from c.m.s. to ALICE lab frame
322       Double_t dyBoost = 0.47;
323       Double_t beta  = TMath::TanH(dyBoost);
324       Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
325       Double_t gb    = gamma * beta;
326       yq[0] += dyBoost;
327       yq[1] += dyBoost;
328       plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
329       plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
330       if (fEnergy == 9) {
331         yq[0] *= -1;
332         yq[1] *= -1;
333         plh[0] *= -1;
334         plh[1] *= -1;
335       }
336     }      
337
338     // Cuts from AliGenerator
339     
340     // Cut on theta
341     theta=TMath::ATan2(pth[0],plh[0]);
342     if (theta<fThetaMin || theta>fThetaMax) continue;
343     theta=TMath::ATan2(pth[1],plh[1]);
344     if (theta<fThetaMin || theta>fThetaMax) continue;
345     
346     // Cut on momentum
347     ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
348     if (ph[0]<fPMin || ph[0]>fPMax) continue;
349     ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
350     if (ph[1]<fPMin || ph[1]>fPMax) continue;
351     
352     // Add the quarks in the stack
353     
354     phiq[0] = Rndm()*k2PI;
355     if (Rndm() < 0.5) {
356       phiq[1] = phiq[0] + dphi*kDegrad; 
357     } else {
358       phiq[1] = phiq[0] - dphi*kDegrad; 
359     }    
360     if (phiq[1] > k2PI) phiq[1] -= k2PI;
361     if (phiq[1] < 0   ) phiq[1] += k2PI;
362     
363     // quarks pdg
364     iquark[0] = +fQuark;
365     iquark[1] = -fQuark;
366
367     // px and py
368     TVector2 qvect1 = TVector2();
369     TVector2 qvect2 = TVector2();
370     qvect1.SetMagPhi(ptq[0],phiq[0]);
371     qvect2.SetMagPhi(ptq[1],phiq[1]);
372     pq[0][0] = qvect1.Px();
373     pq[0][1] = qvect1.Py();
374     pq[1][0] = qvect2.Px();
375     pq[1][1] = qvect2.Py();
376
377     // pz
378     tanhy2 = TMath::TanH(yq[0]);
379     tanhy2 *= tanhy2;
380     pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
381     pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
382     tanhy2 = TMath::TanH(yq[1]);
383     tanhy2 *= tanhy2;
384     pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
385     pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
386
387     // Here we assume that  |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
388     // which is a good approximation for heavy flavors in Pythia
389     // ... moreover, same phi angles as for the quarks ...
390     
391     phih[0] = phiq[0];    
392     phih[1] = phiq[1];    
393
394     for (Int_t ihad = 0; ihad < 2; ihad++) {
395     while(1) {
396       //
397       // particle type 
398       Int_t iPart = ihadron[ihad];
399       fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
400       wgtp=fParentWeight;
401       wgtch=fChildWeight;
402       TParticlePDG *particle = pDataBase->GetParticle(iPart);
403       Float_t am = particle->Mass();
404       phi = phih[ihad];
405       pt  = pth[ihad];
406       pl  = plh[ihad];
407       ptot=TMath::Sqrt(pt*pt+pl*pl);
408
409       p[0]=pt*TMath::Cos(phi);
410       p[1]=pt*TMath::Sin(phi);
411       p[2]=pl;
412
413       if(fVertexSmear==kPerTrack) {
414         Rndm(random,6);
415         for (j=0;j<3;j++) {
416           origin0[j]=
417             fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
418             TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
419         }
420       }
421       
422       // Looking at fForceDecay : 
423       // if fForceDecay != none Primary particle decays using 
424       // AliPythia and children are tracked by GEANT
425       //
426       // if fForceDecay == none Primary particle is tracked by GEANT 
427       // (In the latest, make sure that GEANT actually does all the decays you want)      
428       //
429       
430       if (fForceDecay != kNoDecay) {
431         // Using lujet to decay particle
432         Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
433         TLorentzVector pmom(p[0], p[1], p[2], energy);
434         fDecayer->Decay(iPart,&pmom);
435         //
436         // select decay particles
437         Int_t np=fDecayer->ImportParticles(particles);
438         
439         //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
440         if (fGeometryAcceptance) 
441           if (!CheckAcceptanceGeometry(np,particles)) continue;
442         Int_t ncsel=0;
443         Int_t* pFlag      = new Int_t[np];
444         Int_t* pParent    = new Int_t[np];
445         Int_t* pSelected  = new Int_t[np];
446         Int_t* trackIt    = new Int_t[np];
447         
448         for (i=0; i<np; i++) {
449           pFlag[i]     =  0;
450           pSelected[i] =  0;
451           pParent[i]   = -1;
452         }
453         
454         if (np >1) {
455           TParticle* iparticle = 0;
456           Int_t ipF, ipL;
457           for (i = 1; i<np ; i++) {
458             trackIt[i] = 1;
459             iparticle = (TParticle *) particles->At(i);
460             Int_t kf = iparticle->GetPdgCode();
461             Int_t ks = iparticle->GetStatusCode();
462             // flagged particle
463             
464             if (pFlag[i] == 1) {
465               ipF = iparticle->GetFirstDaughter();
466               ipL = iparticle->GetLastDaughter();       
467               if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
468               continue;
469             }
470             
471             // flag decay products of particles with long life-time (c tau > .3 mum)                  
472             
473             if (ks != 1) { 
474               //TParticlePDG *particle = pDataBase->GetParticle(kf);
475               
476               Double_t lifeTime = fDecayer->GetLifetime(kf);
477               //Double_t mass     = particle->Mass();
478               //Double_t width    = particle->Width();
479               if (lifeTime > (Double_t) fMaxLifeTime) {
480                 ipF = iparticle->GetFirstDaughter();
481                 ipL = iparticle->GetLastDaughter();     
482                 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
483               } else {
484                 trackIt[i]     = 0;
485                 pSelected[i]   = 1;
486               }
487             } // ks==1 ?
488             //
489             // children
490             
491             if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[i])
492               {
493                 if (fCutOnChild) {
494                   pc[0]=iparticle->Px();
495                   pc[1]=iparticle->Py();
496                   pc[2]=iparticle->Pz();
497                   Bool_t  childok = KinematicSelection(iparticle, 1);
498                   if(childok) {
499                     pSelected[i]  = 1;
500                     ncsel++;
501                   } else {
502                     ncsel=-1;
503                     break;
504                   } // child kine cuts
505                 } else {
506                   pSelected[i]  = 1;
507                   ncsel++;
508                 } // if child selection
509               } // select muon
510           } // decay particle loop
511         } // if decay products
512         
513         Int_t iparent;
514         if ((fCutOnChild && ncsel >0) || !fCutOnChild){
515           ipa++;
516           //
517           // Parent
518           // quark
519           PushTrack(0, -1, iquark[ihad], pq[ihad], origin0, polar, 0, kPPrimary, nt, wgtp);
520           KeepTrack(nt);
521           ntq[ihad] = nt;
522           // hadron
523           PushTrack(0, ntq[ihad], iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
524           pParent[0] = nt;
525           KeepTrack(nt); 
526           fNprimaries++;
527           
528           //
529           // Decay Products
530           //              
531           for (i = 1; i < np; i++) {
532             if (pSelected[i]) {
533               TParticle* iparticle = (TParticle *) particles->At(i);
534               Int_t kf  = iparticle->GetPdgCode();
535               Int_t jpa = iparticle->GetFirstMother()-1;
536               
537               och[0] = origin0[0]+iparticle->Vx()/10;
538               och[1] = origin0[1]+iparticle->Vy()/10;
539               och[2] = origin0[2]+iparticle->Vz()/10;
540               pc[0]  = iparticle->Px();
541               pc[1]  = iparticle->Py();
542               pc[2]  = iparticle->Pz();
543               
544               if (jpa > -1) {
545                 iparent = pParent[jpa];
546               } else {
547                 iparent = -1;
548               }
549               
550               PushTrack(fTrackIt*trackIt[i], iparent, kf,
551                         pc, och, polar,
552                         0, kPDecay, nt, wgtch);
553               pParent[i] = nt;
554               KeepTrack(nt); 
555               fNprimaries++;
556             } // Selected
557           } // Particle loop 
558         }  // Decays by Lujet
559         particles->Clear();
560         if (pFlag)      delete[] pFlag;
561         if (pParent)    delete[] pParent;
562         if (pSelected)  delete[] pSelected;        
563         if (trackIt)    delete[] trackIt;
564       } // kinematic selection
565       else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
566         {
567           gAlice->GetMCApp()->
568             PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
569           ipa++; 
570           fNprimaries++;
571         }
572       break;
573     } // while(1) loop
574     } // hadron pair loop
575   }   // while (ipa<2*fNpart) loop
576
577   SetHighWaterMark(nt);
578   
579   AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
580   header->SetPrimaryVertex(fVertex);
581   header->SetNProduced(fNprimaries);
582   AddHeader(header);
583
584   delete particles;
585
586 }
587 //____________________________________________________________________________________
588 void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
589 {  
590 // Composition of a lower state charm hadron pair from a ccbar quark pair
591    Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};
592
593    Double_t id3, id4;
594    hProbHH->GetRandom2(id3, id4);
595    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
596    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
597
598    return;
599 }
600
601 void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
602 {  
603 // Composition of a lower state beauty hadron pair from a bbbar quark pair
604    // B-Bbar mixing will be done by Pythia at their decay point
605    Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};
606
607    Double_t id3, id4;
608    hProbHH->GetRandom2(id3, id4);
609    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
610    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
611
612    if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
613    if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;
614
615    return;
616 }
617
618 //____________________________________________________________________________________
619 Double_t AliGenCorrHF::ComputeIntegral(TFile* fG)       // needed by GetQuarkPair
620 {
621    // Read QQbar kinematical 5D grid's cell occupancy weights
622    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
623    TTree* tG = (TTree*) fG->Get("tGqq");
624    tG->GetBranch("cell")->SetAddress(&cell);
625    Int_t nbins = tG->GetEntries();
626
627    //   delete previously computed integral (if any)
628    if(fgIntegral) delete [] fgIntegral;
629
630    fgIntegral = new Double_t[nbins+1];
631    fgIntegral[0] = 0;
632    Int_t bin;
633    for(bin=0;bin<nbins;bin++) {
634      tG->GetEvent(bin);
635      fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
636    }
637    //   Normalize integral to 1
638    if (fgIntegral[nbins] == 0 ) {
639       return 0;
640    }
641    for (bin=1;bin<=nbins;bin++)  fgIntegral[bin] /= fgIntegral[nbins];
642
643    return fgIntegral[nbins];
644 }
645
646 //____________________________________________________________________________________
647 void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)              
648                                  // modification of ROOT's TH3::GetRandom3 for 5D
649 {
650    // Read QQbar kinematical 5D grid's cell coordinates
651    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
652    TTree* tG = (TTree*) fG->Get("tGqq");
653    tG->GetBranch("cell")->SetAddress(&cell);
654    Int_t nbins = tG->GetEntries();
655    Double_t rand[6];
656    gRandom->RndmArray(6,rand);
657    Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
658    tG->GetEvent(ibin);
659    y1   = fgy[cell[1]]  + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
660    y2   = fgy[cell[2]]  + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
661    pt1  = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
662    pt2  = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
663    dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
664 }
665
666 //____________________________________________________________________________________
667 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) 
668 {
669     // Generate a hadron pair
670     void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
671     fIpParaFunc = IpCharm;
672     Double_t mq = 1.2;              // c & b quark masses (used in AliPythia)
673     if (idq == 5) {
674       fIpParaFunc = IpBeauty;
675       mq = 4.75;
676     }
677     Double_t z11 = 0.;
678     Double_t z12 = 0.;
679     Double_t z21 = 0.;
680     Double_t z22 = 0.;
681     Double_t pz1, pz2, e1, e2, mh, ptemp, rand[2];
682     char tag[100]; 
683     TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
684     for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
685       snprintf(tag,100, "h2h_pt%d",ipt); 
686       h2h[ipt] = (TH2F*) fG->Get(tag); 
687       snprintf(tag,100, "h2s_pt%d",ipt); 
688       h2s[ipt] = (TH2F*) fG->Get(tag); 
689     }
690
691        if (y1*y2 < 0) {
692          for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
693            if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
694              h2h[ipt]->GetRandom2(z11, z21);
695            if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
696              h2h[ipt]->GetRandom2(z12, z22); 
697          }
698        }
699        else {
700          if (TMath::Abs(y1) > TMath::Abs(y2)) {
701            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
702              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
703                h2h[ipt]->GetRandom2(z11, z21);
704              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
705                h2s[ipt]->GetRandom2(z12, z22); 
706            }
707          }
708          else {
709            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
710              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
711                h2s[ipt]->GetRandom2(z11, z21);
712              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
713                h2h[ipt]->GetRandom2(z12, z22); 
714            }
715          }
716        }
717       gRandom->RndmArray(2,rand);
718       ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
719       pz1   = ptemp*TMath::SinH(y1); 
720       e1    = ptemp*TMath::CosH(y1); 
721       ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
722       pz2   = ptemp*TMath::SinH(y2); 
723       e2    = ptemp*TMath::CosH(y2); 
724
725       hProbHH = (TH2F*)fG->Get("hProbHH");
726       fIpParaFunc(hProbHH, id3, id4);
727       mh    = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
728       ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
729       if (idq==5) pt3   = pt1;                // an approximation at low pt, try better
730       else        pt3   = rand[0];            // pt3=pt1 gives less D-hadrons at low pt 
731       if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
732       if (pz1 > 0)   pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
733       else           pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
734       e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
735
736       mh    = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
737       ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
738       if (idq==5) pt4   = pt2;                // an approximation at low pt, try better
739       else        pt4   = rand[1];
740       if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
741       if (pz2 > 0)   pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
742       else           pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
743       e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
744
745       // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
746       Float_t ycorr = 0.2, y3, y4;
747       gRandom->RndmArray(2,rand);
748       y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
749       y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
750       if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
751         ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
752         y3  = 4*(1 - 2*rand[1]);
753         pz3 = ptemp*TMath::SinH(y3);
754         pz4 = pz3;
755       }
756 }
757