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