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