]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenCorrHF.cxx
Updates by Constantin Loizides.
[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         fDecayer(0),
136         fgIntegral(0)
137 {
138 // Default constructor
139 }
140
141 //____________________________________________________________
142 AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t idquark, Int_t energy):
143     AliGenMC(npart),
144     fFileName(0),
145     fFile(0),
146     fQuark(idquark),
147     fEnergy(energy),
148     fBias(0.),
149     fTrials(0),
150     //    fDecayer(new AliDecayerPythia())
151     fDecayer(0),
152     fgIntegral(0)
153 {
154 // Constructor using particle number, quark type, energy & default InputFile
155 //
156     if (fQuark == 5) {
157       if (fEnergy == 7)
158            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP7PythiaMNRwmi.root";
159       else if (fEnergy == 10)
160            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP10PythiaMNRwmi.root";
161       else if (fEnergy == 14)
162            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP14PythiaMNRwmi.root";
163       else if (fEnergy == 4)
164            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
165       else if (fEnergy == 9 || fEnergy == -9)
166            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb88PythiaMNR.root";
167       else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
168     }
169     else {
170       fQuark = 4;
171       if (fEnergy == 7)
172            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP7PythiaMNRwmi.root";
173       else if (fEnergy == 10)
174            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP10PythiaMNRwmi.root";
175       else if (fEnergy == 14)
176            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP14PythiaMNRwmi.root";
177       else if (fEnergy == 4)
178            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
179       else if (fEnergy == 9 || fEnergy == -9)
180            fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb88PythiaMNR.root";
181       else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
182     }
183     fName = "Default";
184     fTitle= "Generator for correlated pairs of HF hadrons";
185       
186     fChildSelect.Set(5);
187     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
188     SetForceDecay();
189     SetCutOnChild();
190     SetChildMomentumRange();
191     SetChildPtRange();
192     SetChildPhiRange();
193     SetChildThetaRange(); 
194 }
195
196 //___________________________________________________________________
197 AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t idquark, Int_t energy):
198     AliGenMC(npart),
199     fFileName(tname),
200     fFile(0),
201     fQuark(idquark),
202     fEnergy(energy),
203     fBias(0.),
204     fTrials(0),
205     //    fDecayer(new AliDecayerPythia())
206     fDecayer(0),
207     fgIntegral(0)
208 {
209 // Constructor using particle number, quark type, energy & user-defined InputFile
210 //
211     if (fQuark != 5) fQuark = 4;
212     fName = "UserDefined";
213     fTitle= "Generator for correlated pairs of HF hadrons";
214       
215     fChildSelect.Set(5);
216     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
217     SetForceDecay();
218     SetCutOnChild();
219     SetChildMomentumRange();
220     SetChildPtRange();
221     SetChildPhiRange();
222     SetChildThetaRange(); 
223 }
224
225 //____________________________________________________________
226 AliGenCorrHF::~AliGenCorrHF()
227 {
228 // Destructor
229   delete fFile;
230 }
231
232 //____________________________________________________________
233 void AliGenCorrHF::Init()
234 {
235 // Initialisation
236
237   AliInfo(Form(" QQbar kinematics and fragm. functions from:  %s",fFileName.Data() )); 
238     fFile = TFile::Open(fFileName.Data());
239     if(!fFile->IsOpen()){
240       AliError(Form("Could not open file %s",fFileName.Data() ));
241     }
242
243     ComputeIntegral(fFile);
244     
245     fParentWeight = 1./fNpart;   // fNpart is number of HF-hadron pairs
246
247 // particle decay related initialization
248
249     if (gMC) fDecayer = gMC->GetDecayer();
250     fDecayer->SetForceDecay(fForceDecay);
251     fDecayer->Init();
252
253 //
254     AliGenMC::Init();
255 }
256 //____________________________________________________________
257 void AliGenCorrHF::Generate()
258 {
259 //
260 // Generate fNpart of correlated HF hadron pairs per event
261 // in the the desired theta and momentum windows (phi = 0 - 2pi).
262 //
263
264 //  Reinitialize decayer
265
266   fDecayer->SetForceDecay(fForceDecay);
267   fDecayer->Init();
268   
269   //
270   Float_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
271   Float_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
272   Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
273   Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
274   Float_t p[3], pc[3], och[3];// Momentum, polarisation and origin of the children particles from lujet
275   Int_t nt, i, j, ipa, ihadron[2], iquark[2];
276   Float_t  wgtp, wgtch, random[6];
277   Float_t pq[2][3];           // Momenta of the two quarks
278   Int_t ntq[2] = {-1, -1};
279   Double_t tanhy2, qm = 0;
280
281   Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2], phiq[2];
282   for (i=0; i<2; i++) { 
283     ptq[i]     =0; 
284     yq[i]      =0; 
285     pth[i]     =0; 
286     plh[i]     =0;
287     phih[i]    =0; 
288     phiq[i]    =0;
289     ihadron[i] =0; 
290     iquark[i]  =0;
291   }
292
293   // same quarks mass as in the fragmentation functions
294   if (fQuark == 4) qm = 1.20;
295   else             qm = 4.75;
296
297   TClonesArray *particles = new TClonesArray("TParticle",1000);
298   
299   TDatabasePDG *pDataBase = TDatabasePDG::Instance();
300   //
301   
302   // Calculating vertex position per event
303   for (j=0;j<3;j++) origin0[j]=fOrigin[j];
304   if (fVertexSmear==kPerEvent) {
305     Vertex();
306     for (j=0;j<3;j++) origin0[j]=fVertex[j];
307   }
308   
309   ipa=0;
310   
311   // Generating fNpart HF-hadron pairs
312   fNprimaries = 0;
313   
314   while (ipa<2*fNpart) {
315     
316     GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
317     
318     GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
319     
320     if (fEnergy == 9 || fEnergy == -9) {      // boost particles from c.m.s. to ALICE lab frame
321       Double_t dyBoost = 0.47;
322       Double_t beta  = TMath::TanH(dyBoost);
323       Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
324       Double_t gb    = gamma * beta;
325       yq[0] += dyBoost;
326       yq[1] += dyBoost;
327       plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
328       plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
329       if (fEnergy == 9) {
330         yq[0] *= -1;
331         yq[1] *= -1;
332         plh[0] *= -1;
333         plh[1] *= -1;
334       }
335     }      
336
337     // Cuts from AliGenerator
338     
339     // Cut on theta
340     theta=TMath::ATan2(pth[0],plh[0]);
341     if (theta<fThetaMin || theta>fThetaMax) continue;
342     theta=TMath::ATan2(pth[1],plh[1]);
343     if (theta<fThetaMin || theta>fThetaMax) continue;
344     
345     // Cut on momentum
346     ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
347     if (ph[0]<fPMin || ph[0]>fPMax) continue;
348     ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
349     if (ph[1]<fPMin || ph[1]>fPMax) continue;
350     
351     // Add the quarks in the stack
352     
353     phiq[0] = Rndm()*k2PI;
354     if (Rndm() < 0.5) {
355       phiq[1] = phiq[0] + dphi*kDegrad; 
356     } else {
357       phiq[1] = phiq[0] - dphi*kDegrad; 
358     }    
359     if (phiq[1] > k2PI) phiq[1] -= k2PI;
360     if (phiq[1] < 0   ) phiq[1] += k2PI;
361     
362     // quarks pdg
363     iquark[0] = +fQuark;
364     iquark[1] = -fQuark;
365
366     // px and py
367     TVector2 qvect1 = TVector2();
368     TVector2 qvect2 = TVector2();
369     qvect1.SetMagPhi(ptq[0],phiq[0]);
370     qvect2.SetMagPhi(ptq[1],phiq[1]);
371     pq[0][0] = qvect1.Px();
372     pq[0][1] = qvect1.Py();
373     pq[1][0] = qvect2.Px();
374     pq[1][1] = qvect2.Py();
375
376     // pz
377     tanhy2 = TMath::TanH(yq[0]);
378     tanhy2 *= tanhy2;
379     pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
380     pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
381     tanhy2 = TMath::TanH(yq[1]);
382     tanhy2 *= tanhy2;
383     pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
384     pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
385
386     // Here we assume that  |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
387     // which is a good approximation for heavy flavors in Pythia
388     // ... moreover, same phi angles as for the quarks ...
389     
390     phih[0] = phiq[0];    
391     phih[1] = phiq[1];    
392
393     for (Int_t ihad = 0; ihad < 2; ihad++) {
394     while(1) {
395       //
396       // particle type 
397       Int_t iPart = ihadron[ihad];
398       fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
399       wgtp=fParentWeight;
400       wgtch=fChildWeight;
401       TParticlePDG *particle = pDataBase->GetParticle(iPart);
402       Float_t am = particle->Mass();
403       phi = phih[ihad];
404       pt  = pth[ihad];
405       pl  = plh[ihad];
406       ptot=TMath::Sqrt(pt*pt+pl*pl);
407
408       p[0]=pt*TMath::Cos(phi);
409       p[1]=pt*TMath::Sin(phi);
410       p[2]=pl;
411
412       if(fVertexSmear==kPerTrack) {
413         Rndm(random,6);
414         for (j=0;j<3;j++) {
415           origin0[j]=
416             fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
417             TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
418         }
419       }
420       
421       // Looking at fForceDecay : 
422       // if fForceDecay != none Primary particle decays using 
423       // AliPythia and children are tracked by GEANT
424       //
425       // if fForceDecay == none Primary particle is tracked by GEANT 
426       // (In the latest, make sure that GEANT actually does all the decays you want)      
427       //
428       
429       if (fForceDecay != kNoDecay) {
430         // Using lujet to decay particle
431         Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
432         TLorentzVector pmom(p[0], p[1], p[2], energy);
433         fDecayer->Decay(iPart,&pmom);
434         //
435         // select decay particles
436         Int_t np=fDecayer->ImportParticles(particles);
437         
438         //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
439         if (fGeometryAcceptance) 
440           if (!CheckAcceptanceGeometry(np,particles)) continue;
441         Int_t ncsel=0;
442         Int_t* pFlag      = new Int_t[np];
443         Int_t* pParent    = new Int_t[np];
444         Int_t* pSelected  = new Int_t[np];
445         Int_t* trackIt    = new Int_t[np];
446         
447         for (i=0; i<np; i++) {
448           pFlag[i]     =  0;
449           pSelected[i] =  0;
450           pParent[i]   = -1;
451         }
452         
453         if (np >1) {
454           TParticle* iparticle =  (TParticle *) particles->At(0);
455           Int_t ipF, ipL;
456           for (i = 1; i<np ; i++) {
457             trackIt[i] = 1;
458             iparticle = (TParticle *) particles->At(i);
459             Int_t kf = iparticle->GetPdgCode();
460             Int_t ks = iparticle->GetStatusCode();
461             // flagged particle
462             
463             if (pFlag[i] == 1) {
464               ipF = iparticle->GetFirstDaughter();
465               ipL = iparticle->GetLastDaughter();       
466               if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
467               continue;
468             }
469             
470             // flag decay products of particles with long life-time (c tau > .3 mum)                  
471             
472             if (ks != 1) { 
473               //TParticlePDG *particle = pDataBase->GetParticle(kf);
474               
475               Double_t lifeTime = fDecayer->GetLifetime(kf);
476               //Double_t mass     = particle->Mass();
477               //Double_t width    = particle->Width();
478               if (lifeTime > (Double_t) fMaxLifeTime) {
479                 ipF = iparticle->GetFirstDaughter();
480                 ipL = iparticle->GetLastDaughter();     
481                 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
482               } else {
483                 trackIt[i]     = 0;
484                 pSelected[i]   = 1;
485               }
486             } // ks==1 ?
487             //
488             // children
489             
490             if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll) && trackIt[i])
491               {
492                 if (fCutOnChild) {
493                   pc[0]=iparticle->Px();
494                   pc[1]=iparticle->Py();
495                   pc[2]=iparticle->Pz();
496                   Bool_t  childok = KinematicSelection(iparticle, 1);
497                   if(childok) {
498                     pSelected[i]  = 1;
499                     ncsel++;
500                   } else {
501                     ncsel=-1;
502                     break;
503                   } // child kine cuts
504                 } else {
505                   pSelected[i]  = 1;
506                   ncsel++;
507                 } // if child selection
508               } // select muon
509           } // decay particle loop
510         } // if decay products
511         
512         Int_t iparent;
513         if ((fCutOnChild && ncsel >0) || !fCutOnChild){
514           ipa++;
515           //
516           // Parent
517           // quark
518           PushTrack(0, -1, iquark[ihad], pq[ihad], origin0, polar, 0, kPPrimary, nt, wgtp);
519           KeepTrack(nt);
520           ntq[ihad] = nt;
521           // hadron
522           PushTrack(0, ntq[ihad], iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
523           pParent[0] = nt;
524           KeepTrack(nt); 
525           fNprimaries++;
526           
527           //
528           // Decay Products
529           //              
530           for (i = 1; i < np; i++) {
531             if (pSelected[i]) {
532               TParticle* iparticle = (TParticle *) particles->At(i);
533               Int_t kf  = iparticle->GetPdgCode();
534               Int_t jpa = iparticle->GetFirstMother()-1;
535               
536               och[0] = origin0[0]+iparticle->Vx()/10;
537               och[1] = origin0[1]+iparticle->Vy()/10;
538               och[2] = origin0[2]+iparticle->Vz()/10;
539               pc[0]  = iparticle->Px();
540               pc[1]  = iparticle->Py();
541               pc[2]  = iparticle->Pz();
542               
543               if (jpa > -1) {
544                 iparent = pParent[jpa];
545               } else {
546                 iparent = -1;
547               }
548               
549               PushTrack(fTrackIt*trackIt[i], iparent, kf,
550                         pc, och, polar,
551                         0, kPDecay, nt, wgtch);
552               pParent[i] = nt;
553               KeepTrack(nt); 
554               fNprimaries++;
555             } // Selected
556           } // Particle loop 
557         }  // Decays by Lujet
558         particles->Clear();
559         if (pFlag)      delete[] pFlag;
560         if (pParent)    delete[] pParent;
561         if (pSelected)  delete[] pSelected;        
562         if (trackIt)    delete[] trackIt;
563       } // kinematic selection
564       else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
565         {
566           gAlice->GetMCApp()->
567             PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
568           ipa++; 
569           fNprimaries++;
570         }
571       break;
572     } // while(1) loop
573     } // hadron pair loop
574   }   // while (ipa<2*fNpart) loop
575
576   SetHighWaterMark(nt);
577   
578   AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
579   header->SetPrimaryVertex(fVertex);
580   header->SetNProduced(fNprimaries);
581   AddHeader(header);
582
583   delete particles;
584
585 }
586 //____________________________________________________________________________________
587 void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
588 {  
589 // Composition of a lower state charm hadron pair from a ccbar quark pair
590    Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};
591
592    Double_t id3, id4;
593    hProbHH->GetRandom2(id3, id4);
594    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
595    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
596
597    return;
598 }
599
600 void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
601 {  
602 // Composition of a lower state beauty hadron pair from a bbbar quark pair
603    // B-Bbar mixing will be done by Pythia at their decay point
604    Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};
605
606    Double_t id3, id4;
607    hProbHH->GetRandom2(id3, id4);
608    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
609    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
610
611    if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
612    if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;
613
614    return;
615 }
616
617 //____________________________________________________________________________________
618 Double_t AliGenCorrHF::ComputeIntegral(TFile* fG)       // needed by GetQuarkPair
619 {
620    // Read QQbar kinematical 5D grid's cell occupancy weights
621    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
622    TTree* tG = (TTree*) fG->Get("tGqq");
623    tG->GetBranch("cell")->SetAddress(&cell);
624    Int_t nbins = tG->GetEntries();
625
626    //   delete previously computed integral (if any)
627    if(fgIntegral) delete [] fgIntegral;
628
629    fgIntegral = new Double_t[nbins+1];
630    fgIntegral[0] = 0;
631    Int_t bin;
632    for(bin=0;bin<nbins;bin++) {
633      tG->GetEvent(bin);
634      fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
635    }
636    //   Normalize integral to 1
637    if (fgIntegral[nbins] == 0 ) {
638       return 0;
639    }
640    for (bin=1;bin<=nbins;bin++)  fgIntegral[bin] /= fgIntegral[nbins];
641
642    return fgIntegral[nbins];
643 }
644
645 //____________________________________________________________________________________
646 void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)              
647                                  // modification of ROOT's TH3::GetRandom3 for 5D
648 {
649    // Read QQbar kinematical 5D grid's cell coordinates
650    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
651    TTree* tG = (TTree*) fG->Get("tGqq");
652    tG->GetBranch("cell")->SetAddress(&cell);
653    Int_t nbins = tG->GetEntries();
654    Double_t rand[6];
655    gRandom->RndmArray(6,rand);
656    Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
657    tG->GetEvent(ibin);
658    y1   = fgy[cell[1]]  + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
659    y2   = fgy[cell[2]]  + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
660    pt1  = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
661    pt2  = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
662    dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
663 }
664
665 //____________________________________________________________________________________
666 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) 
667 {
668     // Generate a hadron pair
669     void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
670     fIpParaFunc = IpCharm;
671     Double_t mq = 1.2;              // c & b quark masses (used in AliPythia)
672     if (idq == 5) {
673       fIpParaFunc = IpBeauty;
674       mq = 4.75;
675     }
676     Double_t z11, z12, z21, z22, pz1, pz2, e1, e2, mh, ptemp, rand[2];
677     char tag[100]; 
678     TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
679     for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
680       sprintf(tag,"h2h_pt%d",ipt); 
681       h2h[ipt] = (TH2F*) fG->Get(tag); 
682       sprintf(tag,"h2s_pt%d",ipt); 
683       h2s[ipt] = (TH2F*) fG->Get(tag); 
684     }
685
686        if (y1*y2 < 0) {
687          for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
688            if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
689              h2h[ipt]->GetRandom2(z11, z21);
690            if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
691              h2h[ipt]->GetRandom2(z12, z22); 
692          }
693        }
694        else {
695          if (TMath::Abs(y1) > TMath::Abs(y2)) {
696            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
697              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
698                h2h[ipt]->GetRandom2(z11, z21);
699              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
700                h2s[ipt]->GetRandom2(z12, z22); 
701            }
702          }
703          else {
704            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
705              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
706                h2s[ipt]->GetRandom2(z11, z21);
707              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
708                h2h[ipt]->GetRandom2(z12, z22); 
709            }
710          }
711        }
712       gRandom->RndmArray(2,rand);
713       ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
714       pz1   = ptemp*TMath::SinH(y1); 
715       e1    = ptemp*TMath::CosH(y1); 
716       ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
717       pz2   = ptemp*TMath::SinH(y2); 
718       e2    = ptemp*TMath::CosH(y2); 
719
720       hProbHH = (TH2F*)fG->Get("hProbHH");
721       fIpParaFunc(hProbHH, id3, id4);
722       mh    = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
723       ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
724       if (idq==5) pt3   = pt1;                // an approximation at low pt, try better
725       else        pt3   = rand[0];            // pt3=pt1 gives less D-hadrons at low pt 
726       if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
727       if (pz1 > 0)   pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
728       else           pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
729       e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
730
731       mh    = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
732       ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
733       if (idq==5) pt4   = pt2;                // an approximation at low pt, try better
734       else        pt4   = rand[1];
735       if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
736       if (pz2 > 0)   pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
737       else           pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
738       e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
739
740       // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
741       Float_t ycorr = 0.2, y3, y4;
742       gRandom->RndmArray(2,rand);
743       y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
744       y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
745       if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
746         ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
747         y3  = 4*(1 - 2*rand[1]);
748         pz3 = ptemp*TMath::SinH(y3);
749         pz4 = pz3;
750       }
751 }
752