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