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