]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenCorrHF.cxx
Make revisions needed for compilation from the main AliRoot tree
[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 //
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
97 #include "AliGenCorrHF.h"
98 #include "AliLog.h"
99 #include "AliConst.h"
100 #include "AliDecayer.h"
101 #include "AliMC.h"
102 #include "AliRun.h"
103
104 ClassImp(AliGenCorrHF)
105
106   //Begin_Html
107   /*
108     <img src="picts/AliGenCorrHF.gif">
109   */
110   //End_Html
111
112 Double_t AliGenCorrHF::fgdph[19] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100, 110, 120, 130, 140, 150, 160, 170, 180};
113 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};
114 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};
115 Int_t AliGenCorrHF::fgnptbins = 12;
116 Double_t AliGenCorrHF::fgptbmin[12] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9};
117 Double_t AliGenCorrHF::fgptbmax[12] = {0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9, 100};
118
119 Double_t* AliGenCorrHF::fgIntegral = 0;
120
121 //____________________________________________________________
122     AliGenCorrHF::AliGenCorrHF():
123         fFileName(0),
124         fFile(0),
125         fQuark(0),
126         fBias(0.),
127         fTrials(0),
128         fDecayer(0)
129 {
130 // Default constructor
131 }
132
133 //____________________________________________________________
134 AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t param):
135     AliGenMC(npart),
136     fFileName(0),
137     fFile(0),
138     fQuark(param),
139     fBias(0.),
140     fTrials(0),
141     //    fDecayer(new AliDecayerPythia())
142     fDecayer(0)
143 {
144 // Constructor using number of particles, quark type & default InputFile
145 //
146     if (fQuark != 5) fQuark = 4;
147     fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmppMNRwmiCorr100K.root";
148     if (fQuark == 5) fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyppMNRwmiCorr100K.root";
149
150     fName = "Default";
151     fTitle= "Generator for correlated pairs of HF hadrons";
152       
153     fChildSelect.Set(5);
154     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
155     SetForceDecay();
156     SetCutOnChild();
157     SetChildMomentumRange();
158     SetChildPtRange();
159     SetChildPhiRange();
160     SetChildThetaRange(); 
161 }
162
163 //___________________________________________________________________
164 AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t param):
165     AliGenMC(npart),
166     fFileName(tname),
167     fFile(0),
168     fQuark(param),
169     fBias(0.),
170     fTrials(0),
171     //    fDecayer(new AliDecayerPythia())
172     fDecayer(0)
173 {
174 // Constructor using number of particles, quark type & user-defined InputFile
175 //
176     if (fQuark != 5) fQuark = 4;
177     fName = "UserDefined";
178     fTitle= "Generator for correlated pairs of HF hadrons";
179       
180     fChildSelect.Set(5);
181     for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
182     SetForceDecay();
183     SetCutOnChild();
184     SetChildMomentumRange();
185     SetChildPtRange();
186     SetChildPhiRange();
187     SetChildThetaRange(); 
188 }
189
190 //____________________________________________________________
191 AliGenCorrHF::~AliGenCorrHF()
192 {
193 // Destructor
194   delete fFile;
195 }
196
197 //____________________________________________________________
198 void AliGenCorrHF::Init()
199 {
200 // Initialisation
201
202   AliInfo(Form(" QQbar kinematics and fragm. functions from:  %s",fFileName.Data() )); 
203     fFile = TFile::Open(fFileName.Data());
204     if(!fFile->IsOpen()){
205       AliError(Form("Could not open file %s",fFileName.Data() ));
206     }
207
208     ComputeIntegral(fFile);
209
210     fParentWeight = 1./fNpart;   // fNpart is number of HF-hadron pairs
211
212 // particle decay related initialization
213
214     if (gMC) fDecayer = gMC->GetDecayer();
215     fDecayer->SetForceDecay(fForceDecay);
216     fDecayer->Init();
217
218 //
219     AliGenMC::Init();
220 }
221
222 //____________________________________________________________
223 void AliGenCorrHF::Generate()
224 {
225 //
226 // Generate fNpart of correlated HF hadron pairs per event
227 // in the the desired theta and momentum windows (phi = 0 - 2pi). 
228 // Gaussian smearing on the vertex is done if selected. 
229 // The decay of heavy hadrons is done using lujet, 
230 //    and the childern particle are tracked by GEANT
231 // However, light mesons are directly tracked by GEANT 
232 // setting fForceDecay = nodecay (SetForceDecay(nodecay)) 
233 //
234
235
236   Float_t polar[3]= {0,0,0};  // Polarisation of the parent particle (for GEANT tracking)
237   Float_t origin0[3];         // Origin of the generated parent particle (for GEANT tracking)
238   Float_t pt, pl, ptot;       // Transverse, logitudinal and total momenta of the parent particle
239   Float_t phi, theta;         // Phi and theta spherical angles of the parent particle momentum
240   Float_t p[3], pc[3], och[3];// Momentum, polarisation and origin of the children particles from lujet
241
242
243   Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2];
244   Int_t i, j, ipair, ihadron[2];
245   for (i=0; i<2; i++) { 
246     ptq[i] =0; 
247     yq[i]  =0; 
248     pth[i] =0; 
249     plh[i] =0; 
250     ihadron[i] =0; 
251   }
252
253   static TClonesArray *particles;
254   //
255   if(!particles) particles = new TClonesArray("TParticle",1000);
256   
257   TDatabasePDG* pDataBase = TDatabasePDG::Instance();
258  
259 // Calculating vertex position per event
260   for (j=0;j<3;j++) origin0[j]=fOrigin[j];
261   if(fVertexSmear==kPerEvent) {
262       Vertex();
263       for (j=0;j<3;j++) origin0[j]=fVertex[j];
264   }
265   
266   Float_t wgtp, wgtch, random[6];
267   Int_t ipap = 0;
268   Int_t nt   = 0;
269
270 // Generating fNpart HF-hadron pairs per event
271   while(ipap<fNpart) {
272
273     while(1) {
274
275       GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
276
277       GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
278
279 // Here we assume that  |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
280 // which is a good approximation for heavy flavors in Pythia
281
282       /* // doesn't work if PhiMax < k2PI or PhiMin > 0, since dphi = 0 - 180 
283       phih[0] = fPhiMin + Rndm()*(fPhiMax-fPhiMin);
284       phih[1] = phih[0] + dphi*kDegrad;
285       if (phih[0] > fPhiMax/2.) phih[1] = phih[0] - dphi*kDegrad;
286       */
287       phih[0] = Rndm()*k2PI;
288       phih[1] = phih[0] + dphi*kDegrad; 
289       if (phih[0] > TMath::Pi()) phih[1] = phih[0] - dphi*kDegrad;
290
291 // Cut on theta
292       theta=TMath::ATan2(pth[0],plh[0]);
293       if(theta<fThetaMin || theta>fThetaMax) continue;
294       theta=TMath::ATan2(pth[1],plh[1]);
295       if(theta<fThetaMin || theta>fThetaMax) continue;
296
297 // Cut on momentum
298       ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
299       if (ph[0]<fPMin || ph[0]>fPMax) continue;
300       ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
301       if (ph[1]<fPMin || ph[1]>fPMax) continue;
302
303 // Common origin for particles of the HF-hadron pair
304       if(fVertexSmear==kPerTrack) {
305           Rndm(random,6);
306           for (j=0;j<3;j++) {
307               origin0[j]=
308                 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
309                 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
310           }
311       }
312
313       Int_t np1=0, kf1[100], select1[100], iparent1[100], trackIt1[100];
314       Float_t  wgtch1=0, p1[3], pc1[100][3], och1[100][3];
315
316       for (j=0; j<3; j++)   p1[j] = 0;
317       for (i=0; i<100; i++) {
318         kf1[i]      = 0;  
319         select1[i]  = 0;  
320         iparent1[i] = 0;  
321         trackIt1[i] = 0;
322         for (j=0; j<3; j++) {
323           pc1[i][j]  = 0;
324           och1[i][j] = 0;
325         }
326       }
327
328 //
329 // Loop over particles of the HF-hadron pair
330       Int_t nhadron = 0;
331       for (ipair=0;ipair<2;ipair++) {
332         phi  = phih[ipair];
333         pl   = plh[ipair];
334         pt   = pth[ipair];
335         ptot = ph[ipair];
336 //
337 // particle type 
338           Int_t iPart = ihadron[ipair];
339           Float_t am  = pDataBase->GetParticle(iPart)->Mass();
340           fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
341
342           wgtp  = fParentWeight;
343           wgtch = fChildWeight;
344
345 //
346           p[0]=pt*TMath::Cos(phi);
347           p[1]=pt*TMath::Sin(phi);
348           p[2]=pl;
349           
350 // Looking at fForceDecay : 
351 // if fForceDecay != none Primary particle decays using 
352 // AliPythia and children are tracked by GEANT
353 //
354 // if fForceDecay == none Primary particle is tracked by GEANT 
355 // (In the latest, make sure that GEANT actually does all the decays you want)    
356 //
357
358           if (fForceDecay != kNoDecay) {
359 // Using lujet to decay particle
360               Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
361               TLorentzVector pmom(p[0], p[1], p[2], energy);
362               fDecayer->Decay(iPart,&pmom);
363 //
364 // select decay particles
365               Int_t np=fDecayer->ImportParticles(particles);
366
367 //  Selecting  GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
368               if (fGeometryAcceptance)
369                 if (!CheckAcceptanceGeometry(np,particles)) break;
370               Int_t ncsel=0;
371               Int_t* pParent    = new Int_t[np];
372               Int_t* pSelected  = new Int_t[np];
373               Int_t* trackIt    = new Int_t[np];
374
375               for (i=0; i<np; i++) {
376                   pSelected[i] =  0;
377                   pParent[i]   = -1;
378               }
379               
380               if (np >1) {
381                   TParticle* iparticle =  (TParticle *) particles->At(0);
382                   for (i=1; i<np; i++) {
383                       trackIt[i] = 1;
384                       iparticle = (TParticle *) particles->At(i);
385                       Int_t kf = iparticle->GetPdgCode();
386                       Int_t ks = iparticle->GetStatusCode();
387
388 // particles with long life-time (c tau > .3 mum)
389                       if (ks != 1) { 
390                           Double_t lifeTime = fDecayer->GetLifetime(kf);
391                           if (lifeTime <= (Double_t) fMaxLifeTime) {
392                               trackIt[i]     = 0;
393                               pSelected[i]   = 1;
394                           }
395                       } // ks==1 ?
396 //
397 // children, discard neutrinos
398                       if (TMath::Abs(kf) == 12 || TMath::Abs(kf) == 14) continue;
399                       if (trackIt[i])
400                       {
401                           if (fCutOnChild) {
402                               pc[0]=iparticle->Px();
403                               pc[1]=iparticle->Py();
404                               pc[2]=iparticle->Pz();
405                               Bool_t  childok = KinematicSelection(iparticle, 1);
406                               if(childok) {
407                                   pSelected[i]  = 1;
408                                   ncsel++;
409                               } else {
410                                   ncsel=-1;
411                                   break;
412                               } // child kine cuts
413                           } else {
414                               pSelected[i]  = 1;
415                               ncsel++;
416                           } // if child selection
417                       } // select muon
418                   } // decay particle loop
419               } // if decay products
420               
421               Int_t iparent;
422               if ((fCutOnChild && ncsel >0) || !fCutOnChild){
423
424                   nhadron++;
425 //
426 // Parents and Decay Products
427                   if (ipair == 0) {
428                       np1    = np;
429                       wgtch1 = wgtch;
430                       p1[0] = p[0]; p1[1] = p[1]; p1[2] = p[2];
431                   } else {
432                       ipap++;
433                       PushTrack(0, -1, ihadron[0], p1, origin0, polar, 0, 
434                                 kPPrimary, nt, wgtp);
435                       KeepTrack(nt);
436                       for (i = 1; i < np1; i++) {
437                         if (select1[i]) {
438                           for (j=0; j<3; j++) {
439                               och[j] = och1[i][j];
440                               pc[j]  = pc1[i][j];
441                           }
442                           PushTrack(fTrackIt*trackIt1[i], iparent1[i], kf1[i], pc, och,
443                                     polar, 0, kPDecay, nt, wgtch1);
444                           KeepTrack(nt);
445                           }
446                       }
447                       PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
448                       KeepTrack(nt);
449                   } 
450                   pParent[0] = nt;
451 //
452 // Decay Products
453                   Int_t ntcount = 0;
454                   for (i = 1; i < np; i++) {
455                     if (pSelected[i]) {
456                           TParticle* iparticle = (TParticle *) particles->At(i);
457                           Int_t kf  = iparticle->GetPdgCode();
458                           Int_t ipa = iparticle->GetFirstMother()-1;
459
460                           och[0] = origin0[0]+iparticle->Vx()/10;
461                           och[1] = origin0[1]+iparticle->Vy()/10;
462                           och[2] = origin0[2]+iparticle->Vz()/10;
463                           pc[0]  = iparticle->Px();
464                           pc[1]  = iparticle->Py();
465                           pc[2]  = iparticle->Pz();
466                           
467                           if (ipa > -1) {
468                               iparent = pParent[ipa];
469                           } else {
470                               iparent = -1;
471                           }
472
473                           if (ipair == 0) {
474                               kf1[i]      = kf;  
475                               select1[i]  = pSelected[i];
476                               iparent1[i] = iparent;
477                               trackIt1[i] = trackIt[i];
478                               for (j=0; j<3; j++) {
479                                   och1[i][j] = och[j];
480                                   pc1[i][j]  = pc[j];
481                               }
482                               ntcount++;
483                           } else {
484                             PushTrack(fTrackIt*trackIt[i], iparent, kf, pc, och,
485                                       polar, 0, kPDecay, nt, wgtch);
486                             KeepTrack(nt); 
487                           } 
488                           pParent[i] = nt + ntcount;
489                     } // Selected
490                   } // Particle loop 
491               }  // Decays by Lujet
492               particles->Clear();
493               if (pParent)    delete[] pParent;
494               if (pSelected)  delete[] pSelected;          
495               if (trackIt)    delete[] trackIt;
496           } // kinematic selection
497           else  // nodecay option, so parent will be tracked by GEANT
498           {
499             nhadron++;
500             if (ipair == 0) {
501                 p1[0] = p[0]; p1[1] = p[1]; p1[2] = p[2];
502             } else {
503                 ipap++;
504                 gAlice->GetMCApp()->
505                 PushTrack(fTrackIt,-1,ihadron[0],p1,origin0,polar,0,kPPrimary,nt,wgtp);
506                 gAlice->GetMCApp()->
507                 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
508             }
509           }
510           if (nhadron == 0) break;
511       } // ipair loop
512       if (nhadron != 2) continue;
513       break;
514     } // while(1)
515     nt++;
516   } // while(ipa<fNpart) --> event loop
517
518   SetHighWaterMark(nt);
519 }
520
521 //____________________________________________________________________________________
522 Int_t AliGenCorrHF::IpCharm(TRandom* ran)
523 {  
524 // Composition of lower state charm hadrons, containing a c-quark
525     Float_t random;
526     Int_t ip;            // +- 411,421,431,4122,4132,4232,4332
527     random = ran->Rndm();
528 //  Rates from Pythia6.214 using 100Kevents with kPyCharmppMNRwmi at 14 TeV.   
529   
530     if (random < 0.6027) {                       
531         ip=421;
532     } else if (random < 0.7962) {
533         ip=411;
534     } else if (random < 0.9127) {
535         ip=431;        
536     } else if (random < 0.9899) {
537         ip=4122;        
538     } else if (random < 0.9948) {
539         ip=4132;        
540     } else if (random < 0.9999) {
541         ip=4232;        
542     } else {
543         ip=4332;
544     }
545     
546     return ip;
547 }
548
549 Int_t AliGenCorrHF::IpBeauty(TRandom* ran)
550 {  
551 // Composition of lower state beauty hadrons, containing a b-quark
552     Float_t random;
553     Int_t ip;            // +- 511,521,531,5122,5132,5232,5332
554     random = ran->Rndm(); 
555 //  Rates from Pythia6.214 using 100Kevents with kPyBeautyppMNRwmi at 14 TeV.   
556                         // B-Bbar mixing will be done by Pythia at the decay point
557  if (random < 0.3965) {                       
558         ip=-511;
559     } else if (random < 0.7930) {
560         ip=-521;
561     } else if (random < 0.9112) {
562         ip=-531;        
563     } else if (random < 0.9887) {
564         ip=5122;        
565     } else if (random < 0.9943) {
566         ip=5132;        
567     } else if (random < 0.9999) {
568         ip=5232;        
569     } else {
570         ip=5332;
571     }
572     
573   return ip;
574 }
575
576 //____________________________________________________________________________________
577 Double_t AliGenCorrHF::ComputeIntegral(TFile* fG)       // needed by GetQuarkPair
578 {
579    // Read QQbar kinematical 5D grid's cell occupancy weights
580    Int_t* cell = new Int_t[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
581    TTree* tG = (TTree*) fG->Get("tGqq");
582    tG->GetBranch("cell")->SetAddress(cell);
583    Int_t nbins = tG->GetEntries();
584
585    //   delete previously computed integral (if any)
586    if(fgIntegral) delete [] fgIntegral;
587
588    fgIntegral = new Double_t[nbins+1];
589    fgIntegral[0] = 0;
590    Int_t bin;
591    for(bin=0;bin<nbins;bin++) {
592      tG->GetEvent(bin);
593      fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
594    }
595    //   Normalize integral to 1
596    if (fgIntegral[nbins] == 0 ) {
597       return 0;
598    }
599    for (bin=1;bin<=nbins;bin++)  fgIntegral[bin] /= fgIntegral[nbins];
600
601    return fgIntegral[nbins];
602 }
603
604 //____________________________________________________________________________________
605 void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)              
606                                  // modification of ROOT's TH3::GetRandom3 for 5D
607 {
608    // Read QQbar kinematical 5D grid's cell coordinates
609    Int_t* cell = new Int_t[6];           // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
610    TTree* tG = (TTree*) fG->Get("tGqq");
611    tG->GetBranch("cell")->SetAddress(cell);
612    Int_t nbins = tG->GetEntries();
613    Double_t rand[6];
614    gRandom->RndmArray(6,rand);
615    Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
616    tG->GetEvent(ibin);
617    y1   = fgy[cell[1]]  + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
618    y2   = fgy[cell[2]]  + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
619    pt1  = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
620    pt2  = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
621    dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
622 }
623
624 //____________________________________________________________________________________
625 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) 
626 {
627     // Generate a hadron pair
628     Int_t (*fIpParaFunc )(TRandom*);//Pointer to particle type parametrisation function
629     fIpParaFunc = IpCharm;
630     Double_t mq = 1.2;              // c & b quark masses (used in AliPythia)
631     if (idq == 5) {
632       fIpParaFunc = IpBeauty;
633       mq = 4.75;
634     }
635     Double_t z11, z12, z21, z22, pz1, pz2, e1, e2, mh, ptemp, rand[2];
636     char tag[100]; 
637     TH2F *h2h[12], *h2s[12];      // hard & soft Fragmentation Functions
638     for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
639       sprintf(tag,"h2h_pt%d",ipt); 
640       h2h[ipt] = (TH2F*) fG->Get(tag); 
641       sprintf(tag,"h2s_pt%d",ipt); 
642       h2s[ipt] = (TH2F*) fG->Get(tag); 
643     }
644
645        if (y1*y2 < 0) {
646          for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
647            if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
648              h2h[ipt]->GetRandom2(z11, z21);
649            if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
650              h2h[ipt]->GetRandom2(z12, z22); 
651          }
652        }
653        else {
654          if (TMath::Abs(y1) > TMath::Abs(y2)) {
655            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
656              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
657                h2h[ipt]->GetRandom2(z11, z21);
658              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
659                h2s[ipt]->GetRandom2(z12, z22); 
660            }
661          }
662          else {
663            for (Int_t ipt = 0; ipt<fgnptbins; ipt++) { 
664              if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt]) 
665                h2s[ipt]->GetRandom2(z11, z21);
666              if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt]) 
667                h2h[ipt]->GetRandom2(z12, z22); 
668            }
669          }
670        }
671       gRandom->RndmArray(2,rand);
672       ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
673       pz1   = ptemp*TMath::SinH(y1); 
674       e1    = ptemp*TMath::CosH(y1); 
675       ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
676       pz2   = ptemp*TMath::SinH(y2); 
677       e2    = ptemp*TMath::CosH(y2); 
678
679       id3   = fIpParaFunc(gRandom);
680       mh    = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
681       ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
682       pt3   = (idq-3)*rand[0];                // some smearing at low pt, try better
683       if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
684       if (pz1 > 0)   pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
685       else           pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
686       e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
687
688       id4   = - fIpParaFunc(gRandom);
689       mh    = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
690       ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
691       pt4   = (idq-3)*rand[1];                // some smearing at low pt, try better
692       if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
693       if (pz2 > 0)   pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
694       else           pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
695       e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
696
697       // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
698       Float_t ycorr = 0.2, y3, y4;
699       gRandom->RndmArray(2,rand);
700       y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
701       y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
702       if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
703         ptemp = TMath::Sqrt(e1*e1 - pz3*pz3);
704         y3  = 4*(1 - 2*rand[1]);
705         pz3 = ptemp*TMath::SinH(y3);
706         pz4 = pz3;
707       }
708 }