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