]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenCorrHF.cxx
Flattened rapidity distributions.
[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   if (fVertexSmear==kPerEvent) {
323     Vertex();
324     for (i=0;i<2;i++){
325       for (j=0;j<3;j++) origin0[i][j]=fVertex[j];
326     }
327   }
328   else {
329     for (i=0;i<2;i++){
330       for (j=0;j<3;j++) origin0[i][j]=fOrigin[j];
331     }
332   }
333   
334   ipa  = 0;
335   ipa1 = 0;
336   ipa0 = 0;
337   
338   // Generating fNpart HF-hadron pairs
339   fNprimaries = 0;
340  
341   while (ipa<2*fNpart) {
342
343     GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
344     
345     GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
346     
347     if (fEnergy == 9 || fEnergy == -9) {      // boost particles from c.m.s. to ALICE lab frame
348       Double_t dyBoost = 0.47;
349       Double_t beta  = TMath::TanH(dyBoost);
350       Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
351       Double_t gb    = gamma * beta;
352       yq[0] += dyBoost;
353       yq[1] += dyBoost;
354       plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
355       plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
356       if (fEnergy == 9) {
357         yq[0] *= -1;
358         yq[1] *= -1;
359         plh[0] *= -1;
360         plh[1] *= -1;
361       }
362     }      
363     
364     // Cuts from AliGenerator
365     
366     // Cut on theta
367     theta=TMath::ATan2(pth[0],plh[0]);
368     if (theta<fThetaMin || theta>fThetaMax) continue;
369     theta=TMath::ATan2(pth[1],plh[1]);
370     if (theta<fThetaMin || theta>fThetaMax) continue;
371     
372     // Cut on momentum
373     ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
374     if (ph[0]<fPMin || ph[0]>fPMax) continue;
375     ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
376     if (ph[1]<fPMin || ph[1]>fPMax) continue;
377     
378     // Add the quarks in the stack
379     
380     phiq[0] = Rndm()*k2PI;
381     if (Rndm() < 0.5) {
382       phiq[1] = phiq[0] + dphi*kDegrad; 
383     } else {
384       phiq[1] = phiq[0] - dphi*kDegrad; 
385     }    
386     if (phiq[1] > k2PI) phiq[1] -= k2PI;
387     if (phiq[1] < 0   ) phiq[1] += k2PI;
388     
389     // quarks pdg
390     iquark[0] = +fQuark;
391     iquark[1] = -fQuark;
392     
393     // px and py
394     TVector2 qvect1 = TVector2();
395     TVector2 qvect2 = TVector2();
396     qvect1.SetMagPhi(ptq[0],phiq[0]);
397     qvect2.SetMagPhi(ptq[1],phiq[1]);
398     pq[0][0] = qvect1.Px();
399     pq[0][1] = qvect1.Py();
400     pq[1][0] = qvect2.Px();
401     pq[1][1] = qvect2.Py();
402     
403     // pz
404     tanhy2 = TMath::TanH(yq[0]);
405     tanhy2 *= tanhy2;
406     pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
407     pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
408     tanhy2 = TMath::TanH(yq[1]);
409     tanhy2 *= tanhy2;
410     pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
411     pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
412     
413     // Here we assume that  |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
414     // which is a good approximation for heavy flavors in Pythia
415     // ... moreover, same phi angles as for the quarks ...
416     
417     phih[0] = phiq[0];    
418     phih[1] = phiq[1];    
419     
420     ipa1 = 0;
421
422     for (ihad = 0; ihad < 2; ihad++) {
423       while(1) {
424         
425         ipa0=ipa1;
426         
427         // particle type 
428         fChildWeight=(fDecayer->GetPartialBranchingRatio(ihadron[ihad]))*fParentWeight;
429         wgtp[ihad]=fParentWeight;
430         wgtch[ihad]=fChildWeight;
431         TParticlePDG *particle = pDataBase->GetParticle(ihadron[ihad]);
432         Float_t am = particle->Mass();
433         phi = phih[ihad];
434         pt  = pth[ihad];
435         pl  = plh[ihad];
436         ptot=TMath::Sqrt(pt*pt+pl*pl);
437         
438         p[ihad][0]=pt*TMath::Cos(phi);
439         p[ihad][1]=pt*TMath::Sin(phi);
440         p[ihad][2]=pl;
441         
442         if(fVertexSmear==kPerTrack) {
443           Rndm(random,6);
444           for (j=0;j<3;j++) {
445             origin0[ihad][j]=
446               fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
447               TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
448           }
449         }
450         
451         // Looking at fForceDecay : 
452         // if fForceDecay != none Primary particle decays using 
453         // AliPythia and children are tracked by GEANT
454         //
455         // if fForceDecay == none Primary particle is tracked by GEANT 
456         // (In the latest, make sure that GEANT actually does all the decays you want)    
457
458         if (fForceDecay != kNoDecay) {
459           // Using lujet to decay particle
460           Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
461           TLorentzVector pmom(p[ihad][0], p[ihad][1], p[ihad][2], energy);
462           fDecayer->Decay(ihadron[ihad],&pmom);
463           
464           // select decay particles
465          
466           np[ihad]=fDecayer->ImportParticles((TClonesArray *)particleslist->At(ihad));
467
468           //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
469           
470           if (fGeometryAcceptance) 
471             if (!CheckAcceptanceGeometry(np[ihad],(TClonesArray*)particleslist->At(ihad))) continue;
472           
473           trackIt[ihad]     = new Int_t [np[ihad]];
474           pSelected[ihad]   = new Int_t [np[ihad]];
475           Int_t* pFlag      = new Int_t [np[ihad]];
476           
477           for (i=0; i<np[ihad]; i++) {
478             pFlag[i]     =  0;
479             pSelected[ihad][i] =  0;
480           }
481
482           if (np[ihad] >1) {
483             TParticle* iparticle = 0;
484             Int_t ipF, ipL;
485             
486             for (i = 1; i<np[ihad] ; i++) {
487               trackIt[ihad][i] = 1;
488               iparticle = 
489                 (TParticle *) ((TClonesArray *) particleslist->At(ihad))->At(i);
490               Int_t kf = iparticle->GetPdgCode();
491               Int_t ks = iparticle->GetStatusCode();
492               // flagged particle
493               if (pFlag[i] == 1) {
494                 ipF = iparticle->GetFirstDaughter();
495                 ipL = iparticle->GetLastDaughter();     
496                 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
497                 continue;
498               }
499               
500               // flag decay products of particles with long life-time (c tau > .3 mum)
501               if (ks != 1) { 
502                 Double_t lifeTime = fDecayer->GetLifetime(kf);
503                 if (lifeTime > (Double_t) fMaxLifeTime) {
504                   ipF = iparticle->GetFirstDaughter();
505                   ipL = iparticle->GetLastDaughter();   
506                   if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
507                 } else {
508                   trackIt[ihad][i]     = 0;
509                   pSelected[ihad][i]   = 1;
510                 }
511               } // ks==1 ?
512               //
513               // children
514               if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[ihad][i])
515                 {      
516                   if (fCutOnChild) {
517                     pc[0]=iparticle->Px();
518                     pc[1]=iparticle->Py();
519                     pc[2]=iparticle->Pz();
520                     //printf("px %f py %f pz %f\n",pc[0],pc[1],pc[2]);
521                     Bool_t  childok = KinematicSelection(iparticle, 1);
522                     if(childok) {
523                       pSelected[ihad][i]  = 1;
524                       ncsel[ihad]++;
525                     } else {
526                       ncsel[ihad]=-1;
527                       break;
528                     } // child kine cuts
529                   } else {
530                     pSelected[ihad][i]  = 1;
531                     ncsel[ihad]++;
532                   } // if child selection
533                 } // select muon
534             } // decay particle loop
535           } // if decay products
536
537           if ((fCutOnChild && ncsel[ihad] >0) || !fCutOnChild) ipa1++;
538
539           if (pFlag) delete[] pFlag;
540
541         } // kinematic selection
542         else  // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
543           {
544             gAlice->GetMCApp()->
545               PushTrack(fTrackIt,-1,ihadron[ihad],p[ihad],origin0[ihad],polar[ihad],0,kPPrimary,nt,wgtp[ihad]);
546             ipa1++; 
547             fNprimaries++;
548             
549           }
550         break;
551       } // while(1) loop
552       if (ipa1<ipa0+1){
553         ipa1=0; 
554         if (pSelected[ihad]) delete pSelected[ihad];
555         if (trackIt[ihad])   delete trackIt[ihad];
556         particleshad1->Clear();
557         particleshad2->Clear();
558         break;
559       }//go out of loop and generate new pair if at least one hadron is rejected
560     } // hadron pair loop
561     if(ipa1==2){ 
562  
563       ipa=ipa+ipa1;
564  
565       if(fForceDecay != kNoDecay){
566         for(ihad=0;ihad<2;ihad++){
567
568           //load tracks in the stack if both hadrons of the pair accepted
569           LoadTracks(iquark[ihad],pq[ihad],ihadron[ihad],p[ihad],np[ihad],
570                      (TClonesArray *)particleslist->At(ihad),origin0[ihad],
571                      polar[ihad],wgtp[ihad],wgtch[ihad],nt,ncsel[ihad],
572                      pSelected[ihad],trackIt[ihad]);
573
574           if (pSelected[ihad]) delete pSelected[ihad];
575           if (trackIt[ihad])   delete trackIt[ihad];
576
577         }
578           particleshad1->Clear();
579           particleshad2->Clear();
580       }
581     }
582   }   // while (ipa<2*fNpart) loop
583   
584   SetHighWaterMark(nt);
585   
586   AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
587   header->SetPrimaryVertex(fVertex);
588   header->SetNProduced(fNprimaries);
589   AddHeader(header);
590   
591   
592   delete particleshad1;
593   delete particleshad2;
594   delete particleslist;
595  
596   delete[] pSelected;
597   delete[] trackIt;
598 }
599 //____________________________________________________________________________________
600 void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
601 {  
602 // Composition of a lower state charm hadron pair from a ccbar quark pair
603    Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};
604
605    Double_t id3, id4;
606    hProbHH->GetRandom2(id3, id4);
607    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
608    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
609
610    return;
611 }
612
613 void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
614 {  
615 // Composition of a lower state beauty hadron pair from a bbbar quark pair
616    // B-Bbar mixing will be done by Pythia at their decay point
617    Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};
618
619    Double_t id3, id4;
620    hProbHH->GetRandom2(id3, id4);
621    pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
622    pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
623
624    if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
625    if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;
626
627    return;
628 }
629
630 //____________________________________________________________________________________
631 Double_t AliGenCorrHF::ComputeIntegral(TFile* fG)       // needed by GetQuarkPair
632 {
633    // Read QQbar kinematical 5D grid's cell occupancy weights
634    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
635    TTree* tG = (TTree*) fG->Get("tGqq");
636    tG->GetBranch("cell")->SetAddress(&cell);
637    Int_t nbins = tG->GetEntries();
638
639    //   delete previously computed integral (if any)
640    if(fgIntegral) delete [] fgIntegral;
641
642    fgIntegral = new Double_t[nbins+1];
643    fgIntegral[0] = 0;
644    Int_t bin;
645    for(bin=0;bin<nbins;bin++) {
646      tG->GetEvent(bin);
647      fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
648    }
649    //   Normalize integral to 1
650    if (fgIntegral[nbins] == 0 ) {
651       return 0;
652    }
653    for (bin=1;bin<=nbins;bin++)  fgIntegral[bin] /= fgIntegral[nbins];
654
655    return fgIntegral[nbins];
656 }
657
658
659 //____________________________________________________________________________________
660 void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)              
661                                  // modification of ROOT's TH3::GetRandom3 for 5D
662 {
663    // Read QQbar kinematical 5D grid's cell coordinates
664    Int_t cell[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
665    TTree* tG = (TTree*) fG->Get("tGqq");
666    tG->GetBranch("cell")->SetAddress(&cell);
667    Int_t nbins = tG->GetEntries();
668    Double_t rand[6];
669    gRandom->RndmArray(6,rand);
670    Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
671    tG->GetEvent(ibin);
672    y1   = fgy[cell[1]]  + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
673    y2   = fgy[cell[2]]  + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
674    pt1  = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
675    pt2  = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
676    dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
677 }
678
679 //____________________________________________________________________________________
680 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) 
681 {
682     // Generate a hadron pair
683     void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
684     fIpParaFunc = IpCharm;
685     Double_t mq = 1.2;              // c & b quark masses (used in AliPythia)
686     if (idq == 5) {
687       fIpParaFunc = IpBeauty;
688       mq = 4.75;
689     }
690     Double_t z11 = 0.;
691     Double_t z12 = 0.;
692     Double_t z21 = 0.;
693     Double_t z22 = 0.;
694     Double_t pz1, pz2, e1, e2, mh, ptemp, rand[2];
695     char tag[100]; 
696     TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
697     for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
698       snprintf(tag,100, "h2h_pt%d",ipt); 
699       h2h[ipt] = (TH2F*) fG->Get(tag); 
700       snprintf(tag,100, "h2s_pt%d",ipt); 
701       h2s[ipt] = (TH2F*) fG->Get(tag); 
702     }
703
704        if (y1*y2 < 0) {
705          for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
706            if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
707              h2h[ipt]->GetRandom2(z11, z21);
708            if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
709              h2h[ipt]->GetRandom2(z12, z22); 
710          }
711        }
712        else {
713          if (TMath::Abs(y1) > TMath::Abs(y2)) {
714            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
715              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
716                h2h[ipt]->GetRandom2(z11, z21);
717              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
718                h2s[ipt]->GetRandom2(z12, z22); 
719            }
720          }
721          else {
722            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
723              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
724                h2s[ipt]->GetRandom2(z11, z21);
725              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
726                h2h[ipt]->GetRandom2(z12, z22); 
727            }
728          }
729        }
730       gRandom->RndmArray(2,rand);
731       ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
732       pz1   = ptemp*TMath::SinH(y1); 
733       e1    = ptemp*TMath::CosH(y1); 
734       ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
735       pz2   = ptemp*TMath::SinH(y2); 
736       e2    = ptemp*TMath::CosH(y2); 
737
738       hProbHH = (TH2F*)fG->Get("hProbHH");
739       fIpParaFunc(hProbHH, id3, id4);
740       mh    = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
741       ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
742       if (idq==5) pt3   = pt1;                // an approximation at low pt, try better
743       else        pt3   = rand[0];            // pt3=pt1 gives less D-hadrons at low pt 
744       if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
745       if (pz1 > 0)   pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
746       else           pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
747       e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
748
749       mh    = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
750       ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
751       if (idq==5) pt4   = pt2;                // an approximation at low pt, try better
752       else        pt4   = rand[1];
753       if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
754       if (pz2 > 0)   pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
755       else           pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
756       e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
757
758       // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
759       Float_t ycorr = 0.2, y3, y4;
760       gRandom->RndmArray(2,rand);
761       y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
762       y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
763       if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
764         ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
765         y3  = 4*(1 - 2*rand[1]);
766         pz3 = ptemp*TMath::SinH(y3);
767         pz4 = pz3;
768       }
769 }
770
771                 
772 //____________________________________________________________________________________
773 void AliGenCorrHF::LoadTracks(Int_t iquark, Float_t *pq, 
774                               Int_t iPart, Float_t *p, 
775                               Int_t np, TClonesArray *particles,
776                               Float_t *origin0, Float_t *polar, 
777                               Float_t wgtp, Float_t wgtch,
778                               Int_t &nt, Int_t ncsel, Int_t *pSelected, 
779                               Int_t *trackIt){
780   Int_t i; 
781   Int_t ntq=-1;
782   Int_t* pParent = new Int_t[np];
783   Float_t pc[3], och[3];
784   Int_t iparent;
785
786   for(i=0;i<np;i++) pParent[i]=-1;
787
788   if ((fCutOnChild && ncsel >0) || !fCutOnChild){  
789     // Parents
790     // quark
791     PushTrack(0, -1, iquark, pq, origin0, polar, 0, kPPrimary, nt, wgtp);
792     KeepTrack(nt);
793     ntq = nt;
794     // hadron
795     PushTrack(0, ntq, iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
796     pParent[0] = nt;
797     KeepTrack(nt); 
798     fNprimaries++;
799
800     // Decay Products  
801     for (i = 1; i < np; i++) {
802       if (pSelected[i]) {
803
804         TParticle* iparticle = (TParticle *) particles->At(i);
805         Int_t kf  = iparticle->GetPdgCode();
806         Int_t jpa = iparticle->GetFirstMother()-1;
807
808         och[0] = origin0[0]+iparticle->Vx()/10;
809         och[1] = origin0[1]+iparticle->Vy()/10;
810         och[2] = origin0[2]+iparticle->Vz()/10;
811         pc[0]  = iparticle->Px();
812         pc[1]  = iparticle->Py();
813         pc[2]  = iparticle->Pz();
814         
815         if (jpa > -1) {
816           iparent = pParent[jpa];
817         } else {
818           iparent = -1;
819         }
820         
821         PushTrack(fTrackIt*trackIt[i], iparent, kf,
822                   pc, och, polar,
823                   0, kPDecay, nt, wgtch);
824         pParent[i] = nt;
825         KeepTrack(nt); 
826         fNprimaries++;
827
828       } // Selected
829     } // Particle loop
830   }
831   if (pParent) delete[] pParent;
832  
833   return;
834 }
835