1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 // ~2.10^6 Pythia6.214 events generated with kCharmppMNRwmi & kBeautyppMNRwmi,
24 // CTEQ5L PDF and Pt_hard = 2.76 GeV/c for p-p collisions at 7, 10 and 14 TeV,
25 // and with kCharmppMNR (Pt_hard = 2.10 GeV/c) & kBeautyppMNR (Pt_hard = 2.75 GeV/c),
26 // CTEQ4L PDF for Pb-Pb at 3.94 TeV, for p-Pb & Pb-p at 8.8 TeV.
27 // Decays are performed by Pythia.
28 // Author: S. Grigoryan, LPC Clermont-Fd & YerPhI, Smbat.Grigoryan@cern.ch
29 // July 07: added quarks in the stack (B. Vulpescu)
30 // April 09: added energy choice between 10 and 14 TeV (S. Grigoryan)
31 // Sept 09: added hadron pair composition probabilities via 2D histo (X.M. Zhang)
32 // Oct 09: added energy choice between 7, 10, 14 TeV (for p-p), 4 TeV (for Pb-Pb),
33 // 9 TeV (for p-Pb) and -9 TeV (for Pb-p) (S. Grigoryan)
34 // April 10: removed "static" from definition of some variables (B. Vulpescu)
35 // May 11: added Flag for transportation of background particles while using
36 // SetForceDecay() function (L. Manceau)
37 // June 11: added modifications allowing the setting of cuts on HF-hadron children.
38 // Quarks, hadrons and decay particles are loaded in the stack outside the loop
39 // of HF-hadrons, when the cuts on their children are satisfied (L. Manceau)
40 // Oct 11: added Pb-Pb at 2.76 TeV (S. Grigoryan)
41 // June 12: added p-Pb & Pb-p at 5 TeV (S. Grigoryan)
43 //-------------------------------------------------------------------------
44 // How it works (for the given flavor and p-p energy):
46 // 1) Reads QQbar kinematical grid (TTree) from the Input file and generates
47 // quark pairs according to the weights of the cells.
48 // It is a 5D grid in y1,y2,pt1,pt2 and deltaphi, with occupancy weights
49 // of the cells obtained from Pythia (see details in GetQuarkPair).
50 // 2) Reads "soft" and "hard" fragmentation functions (12 2D-histograms each,
51 // for 12 pt bins) from the Input file, applies to quarks and produces hadrons
52 // (only lower states, with proportions of species obtained from Pythia).
53 // Fragmentation functions are the same for all hadron species and depend
54 // on 2 variables - light cone energy-momentum fractions:
55 // z1=(E_H + Pz_H)/(E_Q + Pz_Q), z2=(E_H - Pz_H)/(E_Q - Pz_Q).
56 // "soft" & "hard" FFs correspond to "slower" & "faster" quark of a pair
57 // (see details in GetHadronPair). Fragmentation does not depend on p-p energy.
58 // 3) Decays the hadrons and saves all the particles in the event stack in the
59 // following order: HF hadron from Q, then its decay products, then HF hadron
60 // from Qbar, then its decay productes, then next HF hadon pair (if any)
61 // in the same way, etc...
62 // 4) It is fast, e.g., generates the same number of events with a beauty pair
63 // ~15 times faster than AliGenPythia with kBeautyppMNRwmi (w/o tracking)
65 // An Input file for each quark flavor and p-p energy is in EVGEN/dataCorrHF/
66 // One can use also user-defined Input files.
68 // More details could be found in my presentation at DiMuonNet Workshop, Dec 2006:
69 // http://www-dapnia.cea.fr/Sphn/Alice/DiMuonNet.
71 //-------------------------------------------------------------------------
74 // add the following typical lines in Config.C
76 if (!strcmp(option,"corr")) {
77 // An example for correlated charm or beauty hadron pair production at 14 TeV
79 // AliGenCorrHF *gener = new AliGenCorrHF(1, 4, 14); // for charm, 1 pair per event
80 AliGenCorrHF *gener = new AliGenCorrHF(1, 5, 14); // for beauty, 1 pair per event
82 gener->SetMomentumRange(0,9999);
83 gener->SetCutOnChild(0); // 1/0 means cuts on children enable/disable
84 gener->SetChildThetaRange(171.0,178.0);
85 gener->SetOrigin(0,0,0); //vertex position
86 gener->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
87 gener->SetForceDecay(kSemiMuonic);
88 gener->SetSelectAll(kTRUE); //Force the transport of all particles.
89 //Necessary while using a different
90 //option than kAll for SetForceDecay
91 gener->SetTrackingFlag(1); //1: Decay during transport,
92 //0: No Decay during transport
96 // and in aliroot do e.g. gAlice->Run(10,"Config.C") to produce 10 events.
97 // One can include AliGenCorrHF in an AliGenCocktail generator.
98 //--------------------------------------------------------------------------
100 #include <Riostream.h>
102 #include <TClonesArray.h>
103 #include <TDatabasePDG.h>
106 #include <TLorentzVector.h>
108 #include <TParticle.h>
109 #include <TParticlePDG.h>
113 #include <TVirtualMC.h>
114 #include <TVector3.h>
116 #include "AliGenCorrHF.h"
118 #include "AliConst.h"
119 #include "AliDecayer.h"
122 #include "AliGenEventHeader.h"
124 ClassImp(AliGenCorrHF)
128 <img src="picts/AliGenCorrHF.gif">
132 Double_t AliGenCorrHF::fgdph[19] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100, 110, 120, 130, 140, 150, 160, 170, 180};
133 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};
134 Double_t AliGenCorrHF::fgpt[51] = {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, 14, 15, 16, 17, 18, 19, 20.1, 21.5, 23, 24.5, 26, 27.5, 29.1, 31, 33, 35, 37, 39.2, 42, 45, 48, 51, 55.2, 60, 65, 71, 81, 100};
135 Int_t AliGenCorrHF::fgnptbins = 12;
136 Double_t AliGenCorrHF::fgptbmin[12] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9};
137 Double_t AliGenCorrHF::fgptbmax[12] = {0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9, 100};
139 //____________________________________________________________
140 AliGenCorrHF::AliGenCorrHF():
151 // Default constructor
154 //____________________________________________________________
155 AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t idquark, Int_t energy):
167 // Constructor using particle number, quark type, energy & default InputFile
171 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP7PythiaMNRwmi.root";
172 else if (fEnergy == 10)
173 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP10PythiaMNRwmi.root";
174 else if (fEnergy == 14)
175 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP14PythiaMNRwmi.root";
176 else if (fEnergy == 3)
177 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb276PythiaMNR.root";
178 else if (fEnergy == 4)
179 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
180 else if (fEnergy == 5 || fEnergy == -5)
181 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb5PythiaMNR.root";
182 else if (fEnergy == 9 || fEnergy == -9)
183 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb88PythiaMNR.root";
184 else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
189 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP7PythiaMNRwmi.root";
190 else if (fEnergy == 10)
191 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP10PythiaMNRwmi.root";
192 else if (fEnergy == 14)
193 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP14PythiaMNRwmi.root";
194 else if (fEnergy == 3)
195 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb276PythiaMNR.root";
196 else if (fEnergy == 4)
197 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
198 else if (fEnergy == 5 || fEnergy == -5)
199 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb5PythiaMNR.root";
200 else if (fEnergy == 9 || fEnergy == -9)
201 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb88PythiaMNR.root";
202 else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
205 fTitle= "Generator for correlated pairs of HF hadrons";
208 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
211 SetChildMomentumRange();
214 SetChildThetaRange();
217 //___________________________________________________________________
218 AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t idquark, Int_t energy):
230 // Constructor using particle number, quark type, energy & user-defined InputFile
232 if (fQuark != 5) fQuark = 4;
233 fName = "UserDefined";
234 fTitle= "Generator for correlated pairs of HF hadrons";
237 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
240 SetChildMomentumRange();
243 SetChildThetaRange();
246 //____________________________________________________________
247 AliGenCorrHF::~AliGenCorrHF()
253 //____________________________________________________________
254 void AliGenCorrHF::Init()
257 AliInfo(Form("Number of HF-hadron pairs = %d",fNpart));
258 AliInfo(Form(" QQbar kinematics and fragm. functions from: %s",fFileName.Data() ));
259 fFile = TFile::Open(fFileName.Data());
260 if(!fFile->IsOpen()){
261 AliError(Form("Could not open file %s",fFileName.Data() ));
264 ComputeIntegral(fFile);
266 fParentWeight = 1./fNpart; // fNpart is number of HF-hadron pairs
268 // particle decay related initialization
270 if (gMC) fDecayer = gMC->GetDecayer();
271 fDecayer->SetForceDecay(fForceDecay);
277 //____________________________________________________________
278 void AliGenCorrHF::Generate()
281 // Generate fNpart of correlated HF hadron pairs per event
282 // in the the desired theta and momentum windows (phi = 0 - 2pi).
285 // Reinitialize decayer
287 fDecayer->SetForceDecay(fForceDecay);
290 Float_t polar[2][3]; // Polarisation of the parent particle (for GEANT tracking)
291 Float_t origin0[2][3]; // Origin of the generated parent particle (for GEANT tracking)
292 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
293 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
294 Float_t p[2][3]; // Momenta
295 Int_t nt, i, j, ihad, ipa, ipa0, ipa1, ihadron[2], iquark[2];
296 Float_t wgtp[2], wgtch[2], random[6];
297 Float_t pq[2][3], pc[3]; // Momenta of the two quarks
298 Double_t tanhy2, qm = 0;
300 Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2], phiq[2];
302 Int_t** pSelected = new Int_t* [2];
303 Int_t** trackIt = new Int_t* [2];
305 for (i=0; i<2; i++) {
314 for (j=0; j<3; j++) polar[i][j]=0;
317 // same quarks mass as in the fragmentation functions
318 if (fQuark == 4) qm = 1.20;
321 TClonesArray *particleshad1 = new TClonesArray("TParticle",1000);
322 TClonesArray *particleshad2 = new TClonesArray("TParticle",1000);
324 TList *particleslist = new TList();
325 particleslist->Add(particleshad1);
326 particleslist->Add(particleshad2);
328 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
330 // Calculating vertex position per event
331 if (fVertexSmear==kPerEvent) {
334 for (j=0;j<3;j++) origin0[i][j]=fVertex[j];
339 for (j=0;j<3;j++) origin0[i][j]=fOrigin[j];
347 // Generating fNpart HF-hadron pairs
350 while (ipa<2*fNpart) {
352 GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
354 GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
356 // Boost particles from c.m.s. to ALICE lab frame for p-Pb & Pb-p collisions
357 if (fEnergy == 5 || fEnergy == -5 || fEnergy == 9 || fEnergy == -9) {
358 Double_t dyBoost = 0.47;
359 Double_t beta = TMath::TanH(dyBoost);
360 Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
361 Double_t gb = gamma * beta;
364 plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
365 plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
366 if (fEnergy == 5 || fEnergy == 9) {
374 // Cuts from AliGenerator
377 theta=TMath::ATan2(pth[0],plh[0]);
378 if (theta<fThetaMin || theta>fThetaMax) continue;
379 theta=TMath::ATan2(pth[1],plh[1]);
380 if (theta<fThetaMin || theta>fThetaMax) continue;
383 ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
384 if (ph[0]<fPMin || ph[0]>fPMax) continue;
385 ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
386 if (ph[1]<fPMin || ph[1]>fPMax) continue;
388 // Add the quarks in the stack
390 phiq[0] = Rndm()*k2PI;
392 phiq[1] = phiq[0] + dphi*kDegrad;
394 phiq[1] = phiq[0] - dphi*kDegrad;
396 if (phiq[1] > k2PI) phiq[1] -= k2PI;
397 if (phiq[1] < 0 ) phiq[1] += k2PI;
404 TVector2 qvect1 = TVector2();
405 TVector2 qvect2 = TVector2();
406 qvect1.SetMagPhi(ptq[0],phiq[0]);
407 qvect2.SetMagPhi(ptq[1],phiq[1]);
408 pq[0][0] = qvect1.Px();
409 pq[0][1] = qvect1.Py();
410 pq[1][0] = qvect2.Px();
411 pq[1][1] = qvect2.Py();
414 tanhy2 = TMath::TanH(yq[0]);
416 pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
417 pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
418 tanhy2 = TMath::TanH(yq[1]);
420 pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
421 pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
423 // Here we assume that |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
424 // which is a good approximation for heavy flavors in Pythia
425 // ... moreover, same phi angles as for the quarks ...
432 for (ihad = 0; ihad < 2; ihad++) {
438 fChildWeight=(fDecayer->GetPartialBranchingRatio(ihadron[ihad]))*fParentWeight;
439 wgtp[ihad]=fParentWeight;
440 wgtch[ihad]=fChildWeight;
441 TParticlePDG *particle = pDataBase->GetParticle(ihadron[ihad]);
442 Float_t am = particle->Mass();
446 ptot=TMath::Sqrt(pt*pt+pl*pl);
448 p[ihad][0]=pt*TMath::Cos(phi);
449 p[ihad][1]=pt*TMath::Sin(phi);
452 if(fVertexSmear==kPerTrack) {
456 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
457 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
461 // Looking at fForceDecay :
462 // if fForceDecay != none Primary particle decays using
463 // AliPythia and children are tracked by GEANT
465 // if fForceDecay == none Primary particle is tracked by GEANT
466 // (In the latest, make sure that GEANT actually does all the decays you want)
468 if (fForceDecay != kNoDecay) {
469 // Using lujet to decay particle
470 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
471 TLorentzVector pmom(p[ihad][0], p[ihad][1], p[ihad][2], energy);
472 fDecayer->Decay(ihadron[ihad],&pmom);
474 // select decay particles
476 np[ihad]=fDecayer->ImportParticles((TClonesArray *)particleslist->At(ihad));
478 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
480 if (fGeometryAcceptance)
481 if (!CheckAcceptanceGeometry(np[ihad],(TClonesArray*)particleslist->At(ihad))) continue;
483 trackIt[ihad] = new Int_t [np[ihad]];
484 pSelected[ihad] = new Int_t [np[ihad]];
485 Int_t* pFlag = new Int_t [np[ihad]];
487 for (i=0; i<np[ihad]; i++) {
489 pSelected[ihad][i] = 0;
493 TParticle* iparticle = 0;
496 for (i = 1; i<np[ihad] ; i++) {
497 trackIt[ihad][i] = 1;
499 (TParticle *) ((TClonesArray *) particleslist->At(ihad))->At(i);
500 Int_t kf = iparticle->GetPdgCode();
501 Int_t ks = iparticle->GetStatusCode();
504 ipF = iparticle->GetFirstDaughter();
505 ipL = iparticle->GetLastDaughter();
506 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
510 // flag decay products of particles with long life-time (c tau > .3 mum)
512 Double_t lifeTime = fDecayer->GetLifetime(kf);
513 if (lifeTime > (Double_t) fMaxLifeTime) {
514 ipF = iparticle->GetFirstDaughter();
515 ipL = iparticle->GetLastDaughter();
516 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
518 trackIt[ihad][i] = 0;
519 pSelected[ihad][i] = 1;
524 if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[ihad][i])
527 pc[0]=iparticle->Px();
528 pc[1]=iparticle->Py();
529 pc[2]=iparticle->Pz();
530 //printf("px %f py %f pz %f\n",pc[0],pc[1],pc[2]);
531 Bool_t childok = KinematicSelection(iparticle, 1);
533 pSelected[ihad][i] = 1;
540 pSelected[ihad][i] = 1;
542 } // if child selection
544 } // decay particle loop
545 } // if decay products
547 if ((fCutOnChild && ncsel[ihad] >0) || !fCutOnChild) ipa1++;
549 if (pFlag) delete[] pFlag;
551 } // kinematic selection
552 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
555 PushTrack(fTrackIt,-1,ihadron[ihad],p[ihad],origin0[ihad],polar[ihad],0,kPPrimary,nt,wgtp[ihad]);
564 if (pSelected[ihad]) delete pSelected[ihad];
565 if (trackIt[ihad]) delete trackIt[ihad];
566 particleshad1->Clear();
567 particleshad2->Clear();
569 }//go out of loop and generate new pair if at least one hadron is rejected
570 } // hadron pair loop
575 if(fForceDecay != kNoDecay){
576 for(ihad=0;ihad<2;ihad++){
578 //load tracks in the stack if both hadrons of the pair accepted
579 LoadTracks(iquark[ihad],pq[ihad],ihadron[ihad],p[ihad],np[ihad],
580 (TClonesArray *)particleslist->At(ihad),origin0[ihad],
581 polar[ihad],wgtp[ihad],wgtch[ihad],nt,ncsel[ihad],
582 pSelected[ihad],trackIt[ihad]);
584 if (pSelected[ihad]) delete pSelected[ihad];
585 if (trackIt[ihad]) delete trackIt[ihad];
588 particleshad1->Clear();
589 particleshad2->Clear();
592 } // while (ipa<2*fNpart) loop
594 SetHighWaterMark(nt);
596 AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
597 header->SetPrimaryVertex(fVertex);
598 header->SetNProduced(fNprimaries);
602 delete particleshad1;
603 delete particleshad2;
604 delete particleslist;
609 //____________________________________________________________________________________
610 void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
612 // Composition of a lower state charm hadron pair from a ccbar quark pair
613 Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};
616 hProbHH->GetRandom2(id3, id4);
617 pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
618 pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
623 void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
625 // Composition of a lower state beauty hadron pair from a bbbar quark pair
626 // B-Bbar mixing will be done by Pythia at their decay point
627 Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};
630 hProbHH->GetRandom2(id3, id4);
631 pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
632 pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
634 if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
635 if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;
640 //____________________________________________________________________________________
641 Double_t AliGenCorrHF::ComputeIntegral(TFile* fG) // needed by GetQuarkPair
643 // Read QQbar kinematical 5D grid's cell occupancy weights
644 Int_t cell[6]; // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
645 TTree* tG = (TTree*) fG->Get("tGqq");
646 tG->GetBranch("cell")->SetAddress(&cell);
647 Int_t nbins = tG->GetEntries();
649 // delete previously computed integral (if any)
650 if(fgIntegral) delete [] fgIntegral;
652 fgIntegral = new Double_t[nbins+1];
655 for(bin=0;bin<nbins;bin++) {
657 fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
659 // Normalize integral to 1
660 if (fgIntegral[nbins] == 0 ) {
663 for (bin=1;bin<=nbins;bin++) fgIntegral[bin] /= fgIntegral[nbins];
665 return fgIntegral[nbins];
669 //____________________________________________________________________________________
670 void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)
671 // modification of ROOT's TH3::GetRandom3 for 5D
673 // Read QQbar kinematical 5D grid's cell coordinates
674 Int_t cell[6]; // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
675 TTree* tG = (TTree*) fG->Get("tGqq");
676 tG->GetBranch("cell")->SetAddress(&cell);
677 Int_t nbins = tG->GetEntries();
679 gRandom->RndmArray(6,rand);
680 Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
682 y1 = fgy[cell[1]] + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
683 y2 = fgy[cell[2]] + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
684 pt1 = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
685 pt2 = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
686 dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
689 //____________________________________________________________________________________
690 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)
692 // Generate a hadron pair
693 void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
694 fIpParaFunc = IpCharm;
695 Double_t mq = 1.2; // c & b quark masses (used in AliPythia)
697 fIpParaFunc = IpBeauty;
704 Double_t pz1, pz2, e1, e2, mh, ptemp, rand[2];
706 TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
707 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
708 snprintf(tag,100, "h2h_pt%d",ipt);
709 h2h[ipt] = (TH2F*) fG->Get(tag);
710 snprintf(tag,100, "h2s_pt%d",ipt);
711 h2s[ipt] = (TH2F*) fG->Get(tag);
715 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
716 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
717 h2h[ipt]->GetRandom2(z11, z21);
718 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
719 h2h[ipt]->GetRandom2(z12, z22);
723 if (TMath::Abs(y1) > TMath::Abs(y2)) {
724 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
725 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
726 h2h[ipt]->GetRandom2(z11, z21);
727 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
728 h2s[ipt]->GetRandom2(z12, z22);
732 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
733 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
734 h2s[ipt]->GetRandom2(z11, z21);
735 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
736 h2h[ipt]->GetRandom2(z12, z22);
740 gRandom->RndmArray(2,rand);
741 ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
742 pz1 = ptemp*TMath::SinH(y1);
743 e1 = ptemp*TMath::CosH(y1);
744 ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
745 pz2 = ptemp*TMath::SinH(y2);
746 e2 = ptemp*TMath::CosH(y2);
748 hProbHH = (TH2F*)fG->Get("hProbHH");
749 fIpParaFunc(hProbHH, id3, id4);
750 mh = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
751 ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
752 if (idq==5) pt3 = pt1; // an approximation at low pt, try better
753 else pt3 = rand[0]; // pt3=pt1 gives less D-hadrons at low pt
754 if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
755 if (pz1 > 0) pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
756 else pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
757 e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
759 mh = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
760 ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
761 if (idq==5) pt4 = pt2; // an approximation at low pt, try better
763 if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
764 if (pz2 > 0) pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
765 else pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
766 e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
768 // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
769 Float_t ycorr = 0.2, y3, y4;
770 gRandom->RndmArray(2,rand);
771 y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
772 y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
773 if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
774 ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
775 y3 = 4*(1 - 2*rand[1]);
776 pz3 = ptemp*TMath::SinH(y3);
782 //____________________________________________________________________________________
783 void AliGenCorrHF::LoadTracks(Int_t iquark, Float_t *pq,
784 Int_t iPart, Float_t *p,
785 Int_t np, TClonesArray *particles,
786 Float_t *origin0, Float_t *polar,
787 Float_t wgtp, Float_t wgtch,
788 Int_t &nt, Int_t ncsel, Int_t *pSelected,
792 Int_t* pParent = new Int_t[np];
793 Float_t pc[3], och[3];
796 for(i=0;i<np;i++) pParent[i]=-1;
798 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
801 PushTrack(0, -1, iquark, pq, origin0, polar, 0, kPPrimary, nt, wgtp);
805 PushTrack(0, ntq, iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
811 for (i = 1; i < np; i++) {
814 TParticle* iparticle = (TParticle *) particles->At(i);
815 Int_t kf = iparticle->GetPdgCode();
816 Int_t jpa = iparticle->GetFirstMother()-1;
817 // RS: note, the conversion mm->cm is done now in the decayer. The time is ignored here!
818 och[0] = origin0[0]+iparticle->Vx();
819 och[1] = origin0[1]+iparticle->Vy();
820 och[2] = origin0[2]+iparticle->Vz();
821 pc[0] = iparticle->Px();
822 pc[1] = iparticle->Py();
823 pc[2] = iparticle->Pz();
826 iparent = pParent[jpa];
831 PushTrack(fTrackIt*trackIt[i], iparent, kf,
833 0, kPDecay, nt, wgtch);
841 if (pParent) delete[] pParent;