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