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