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