]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenCorrHF.cxx
Corrections.
[u/mrichter/AliRoot.git] / EVGEN / AliGenCorrHF.cxx
CommitLineData
2c890605 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
b33adf51 23// ~2.10^6 Pythia6.214 events generated with kCharmppMNRwmi & kBeautyppMNRwmi,
9fd56238 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.
b33adf51 27// Decays are performed by Pythia.
2c890605 28// Author: S. Grigoryan, LPC Clermont-Fd & YerPhI, Smbat.Grigoryan@cern.ch
b44c3901 29// July 07: added quarks in the stack (B. Vulpescu)
b33adf51 30// April 09: added energy choice between 10 and 14 TeV (S. Grigoryan)
9fd56238 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)
5cc3bcd2 34// April 10: removed "static" from definition of some variables (B. Vulpescu)
826fd892 35// May 11: Added Flag for transportation Background While using SetForceDecay() function (L. Manceau)
36// June 11 modifications allowing the setting of cuts on Children added
37//(L. Manceau)
38// Quarks, hadrons and decayed particles are loaded in the stack outside the
39// main loop on heavy hadrons.
40// The particles are loaded only when a pair containing
41// two heavy hadrons given children wich statify cut conditions
42// are tagged in the main loop
43//
2c890605 44//-------------------------------------------------------------------------
b33adf51 45// How it works (for the given flavor and p-p energy):
2c890605 46//
5cc3bcd2 47// 1) Reads QQbar kinematical grid (TTree) from the Input file and generates
2c890605 48// quark pairs according to the weights of the cells.
49// It is a 5D grid in y1,y2,pt1,pt2 and deltaphi, with occupancy weights
50// of the cells obtained from Pythia (see details in GetQuarkPair).
51// 2) Reads "soft" and "hard" fragmentation functions (12 2D-histograms each,
52// for 12 pt bins) from the Input file, applies to quarks and produces hadrons
53// (only lower states, with proportions of species obtained from Pythia).
54// Fragmentation functions are the same for all hadron species and depend
55// on 2 variables - light cone energy-momentum fractions:
56// z1=(E_H + Pz_H)/(E_Q + Pz_Q), z2=(E_H - Pz_H)/(E_Q - Pz_Q).
57// "soft" & "hard" FFs correspond to "slower" & "faster" quark of a pair
b33adf51 58// (see details in GetHadronPair). Fragmentation does not depend on p-p energy.
2c890605 59// 3) Decays the hadrons and saves all the particles in the event stack in the
60// following order: HF hadron from Q, then its decay products, then HF hadron
61// from Qbar, then its decay productes, then next HF hadon pair (if any)
62// in the same way, etc...
63// 4) It is fast, e.g., generates the same number of events with a beauty pair
64// ~15 times faster than AliGenPythia with kBeautyppMNRwmi (w/o tracking)
65//
b33adf51 66// An Input file for each quark flavor and p-p energy is in EVGEN/dataCorrHF/
2c890605 67// One can use also user-defined Input files.
68//
69// More details could be found in my presentation at DiMuonNet Workshop, Dec 2006:
b33adf51 70// http://www-dapnia.cea.fr/Sphn/Alice/DiMuonNet.
2c890605 71//
72//-------------------------------------------------------------------------
73// How to use it:
74//
75// add the following typical lines in Config.C
76/*
77 if (!strcmp(option,"corr")) {
9fd56238 78 // An example for correlated charm or beauty hadron pair production at 14 TeV
2c890605 79
b33adf51 80 // AliGenCorrHF *gener = new AliGenCorrHF(1, 4, 14); // for charm, 1 pair per event
81 AliGenCorrHF *gener = new AliGenCorrHF(1, 5, 14); // for beauty, 1 pair per event
2c890605 82
83 gener->SetMomentumRange(0,9999);
9fd56238 84 gener->SetCutOnChild(0); // 1/0 means cuts on children enable/disable
2c890605 85 gener->SetChildThetaRange(171.0,178.0);
86 gener->SetOrigin(0,0,0); //vertex position
87 gener->SetSigma(0,0,0); //Sigma in (X,Y,Z) (cm) on IP position
88 gener->SetForceDecay(kSemiMuonic);
89 gener->SetTrackingFlag(0);
90 gener->Init();
91}
92*/
93// and in aliroot do e.g. gAlice->Run(10,"Config.C") to produce 10 events.
94// One can include AliGenCorrHF in an AliGenCocktail generator.
95//--------------------------------------------------------------------------
96
7ca4655f 97#include <Riostream.h>
98#include <TCanvas.h>
99#include <TClonesArray.h>
100#include <TDatabasePDG.h>
2c890605 101#include <TFile.h>
2c890605 102#include <TH2F.h>
2c890605 103#include <TLorentzVector.h>
7ca4655f 104#include <TMath.h>
2c890605 105#include <TParticle.h>
106#include <TParticlePDG.h>
7ca4655f 107#include <TROOT.h>
108#include <TRandom.h>
109#include <TTree.h>
2c890605 110#include <TVirtualMC.h>
b44c3901 111#include <TVector3.h>
2c890605 112
113#include "AliGenCorrHF.h"
114#include "AliLog.h"
115#include "AliConst.h"
116#include "AliDecayer.h"
117#include "AliMC.h"
118#include "AliRun.h"
84555c93 119#include "AliGenEventHeader.h"
2c890605 120
121ClassImp(AliGenCorrHF)
122
123 //Begin_Html
124 /*
125 <img src="picts/AliGenCorrHF.gif">
126 */
127 //End_Html
128
129Double_t AliGenCorrHF::fgdph[19] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100, 110, 120, 130, 140, 150, 160, 170, 180};
130Double_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};
b33adf51 131Double_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};
2c890605 132Int_t AliGenCorrHF::fgnptbins = 12;
133Double_t AliGenCorrHF::fgptbmin[12] = {0, 0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9};
134Double_t AliGenCorrHF::fgptbmax[12] = {0.5, 1, 1.5, 2, 2.5, 3, 4, 5, 6, 7, 9, 100};
135
2c890605 136//____________________________________________________________
137 AliGenCorrHF::AliGenCorrHF():
138 fFileName(0),
139 fFile(0),
140 fQuark(0),
b33adf51 141 fEnergy(0),
2c890605 142 fBias(0.),
143 fTrials(0),
18e09c20 144 fSelectAll(kFALSE),
5cc3bcd2 145 fDecayer(0),
146 fgIntegral(0)
2c890605 147{
148// Default constructor
149}
150
151//____________________________________________________________
b33adf51 152AliGenCorrHF::AliGenCorrHF(Int_t npart, Int_t idquark, Int_t energy):
2c890605 153 AliGenMC(npart),
154 fFileName(0),
155 fFile(0),
b33adf51 156 fQuark(idquark),
157 fEnergy(energy),
2c890605 158 fBias(0.),
159 fTrials(0),
18e09c20 160 fSelectAll(kFALSE),
5cc3bcd2 161 fDecayer(0),
162 fgIntegral(0)
2c890605 163{
b33adf51 164// Constructor using particle number, quark type, energy & default InputFile
2c890605 165//
b33adf51 166 if (fQuark == 5) {
9fd56238 167 if (fEnergy == 7)
168 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP7PythiaMNRwmi.root";
169 else if (fEnergy == 10)
170 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP10PythiaMNRwmi.root";
171 else if (fEnergy == 14)
172 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPP14PythiaMNRwmi.root";
173 else if (fEnergy == 4)
174 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
175 else if (fEnergy == 9 || fEnergy == -9)
176 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPPb88PythiaMNR.root";
177 else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/BeautyPbPb394PythiaMNR.root";
b33adf51 178 }
179 else {
180 fQuark = 4;
9fd56238 181 if (fEnergy == 7)
182 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP7PythiaMNRwmi.root";
183 else if (fEnergy == 10)
184 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP10PythiaMNRwmi.root";
185 else if (fEnergy == 14)
186 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPP14PythiaMNRwmi.root";
187 else if (fEnergy == 4)
188 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
189 else if (fEnergy == 9 || fEnergy == -9)
190 fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPPb88PythiaMNR.root";
191 else fFileName = "$ALICE_ROOT/EVGEN/dataCorrHF/CharmPbPb394PythiaMNR.root";
b33adf51 192 }
2c890605 193 fName = "Default";
194 fTitle= "Generator for correlated pairs of HF hadrons";
195
196 fChildSelect.Set(5);
197 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
198 SetForceDecay();
199 SetCutOnChild();
200 SetChildMomentumRange();
201 SetChildPtRange();
202 SetChildPhiRange();
203 SetChildThetaRange();
204}
205
206//___________________________________________________________________
b33adf51 207AliGenCorrHF::AliGenCorrHF(char* tname, Int_t npart, Int_t idquark, Int_t energy):
2c890605 208 AliGenMC(npart),
209 fFileName(tname),
210 fFile(0),
b33adf51 211 fQuark(idquark),
212 fEnergy(energy),
2c890605 213 fBias(0.),
214 fTrials(0),
18e09c20 215 fSelectAll(kFALSE),
5cc3bcd2 216 fDecayer(0),
217 fgIntegral(0)
2c890605 218{
b33adf51 219// Constructor using particle number, quark type, energy & user-defined InputFile
2c890605 220//
221 if (fQuark != 5) fQuark = 4;
222 fName = "UserDefined";
223 fTitle= "Generator for correlated pairs of HF hadrons";
224
225 fChildSelect.Set(5);
226 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
227 SetForceDecay();
228 SetCutOnChild();
229 SetChildMomentumRange();
230 SetChildPtRange();
231 SetChildPhiRange();
232 SetChildThetaRange();
233}
234
2c890605 235//____________________________________________________________
236AliGenCorrHF::~AliGenCorrHF()
237{
238// Destructor
239 delete fFile;
240}
241
242//____________________________________________________________
243void AliGenCorrHF::Init()
244{
245// Initialisation
826fd892 246 AliInfo(Form("Number of HF-hadron pairs = %d",fNpart));
2c890605 247 AliInfo(Form(" QQbar kinematics and fragm. functions from: %s",fFileName.Data() ));
248 fFile = TFile::Open(fFileName.Data());
249 if(!fFile->IsOpen()){
250 AliError(Form("Could not open file %s",fFileName.Data() ));
251 }
252
253 ComputeIntegral(fFile);
b44c3901 254
2c890605 255 fParentWeight = 1./fNpart; // fNpart is number of HF-hadron pairs
256
257// particle decay related initialization
258
259 if (gMC) fDecayer = gMC->GetDecayer();
260 fDecayer->SetForceDecay(fForceDecay);
261 fDecayer->Init();
262
263//
264 AliGenMC::Init();
265}
2c890605 266//____________________________________________________________
267void AliGenCorrHF::Generate()
268{
269//
cb11fd27 270// Generate fNpart of correlated HF hadron pairs per event
271// in the the desired theta and momentum windows (phi = 0 - 2pi).
84555c93 272//
cb11fd27 273
84555c93 274// Reinitialize decayer
2c890605 275
84555c93 276 fDecayer->SetForceDecay(fForceDecay);
277 fDecayer->Init();
826fd892 278
279 Float_t polar[2][3]; // Polarisation of the parent particle (for GEANT tracking)
280 Float_t origin0[2][3]; // Origin of the generated parent particle (for GEANT tracking)
2c890605 281 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
282 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
826fd892 283 Float_t p[2][3]; // Momenta
284 Int_t nt, i, j, ihad, ipa, ipa0, ipa1, ihadron[2], iquark[2];
285 Float_t wgtp[2], wgtch[2], random[6];
286 Float_t pq[2][3], pc[3]; // Momenta of the two quarks
cb11fd27 287 Double_t tanhy2, qm = 0;
826fd892 288 Int_t np[2];
b44c3901 289 Double_t dphi=0, ptq[2], yq[2], pth[2], plh[2], ph[2], phih[2], phiq[2];
826fd892 290 Int_t ncsel[2];
291 Int_t** pSelected = new Int_t* [2];
292 Int_t** trackIt = new Int_t* [2];
293
2c890605 294 for (i=0; i<2; i++) {
b44c3901 295 ptq[i] =0;
296 yq[i] =0;
297 pth[i] =0;
298 plh[i] =0;
299 phih[i] =0;
300 phiq[i] =0;
2c890605 301 ihadron[i] =0;
84555c93 302 iquark[i] =0;
826fd892 303 for (j=0; j<3; j++) polar[i][j]=0;
2c890605 304 }
305
cb11fd27 306 // same quarks mass as in the fragmentation functions
307 if (fQuark == 4) qm = 1.20;
308 else qm = 4.75;
2c890605 309
826fd892 310 TClonesArray *particleshad1 = new TClonesArray("TParticle",1000);
311 TClonesArray *particleshad2 = new TClonesArray("TParticle",1000);
84555c93 312
826fd892 313 TList *particleslist = new TList();
314 particleslist->Add(particleshad1);
315 particleslist->Add(particleshad2);
316
317 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
318
84555c93 319 // Calculating vertex position per event
826fd892 320 for (i=0;i<2;i++){
321 for (j=0;j<3;j++) origin0[i][j]=fOrigin[j];
322 if (fVertexSmear==kPerEvent) {
323 Vertex();
324 for (j=0;j<3;j++) origin0[i][j]=fVertex[j];
325 }
2c890605 326 }
327
826fd892 328 ipa = 0;
329 ipa1 = 0;
330 ipa0 = 0;
84555c93 331
9fd56238 332 // Generating fNpart HF-hadron pairs
84555c93 333 fNprimaries = 0;
826fd892 334
84555c93 335 while (ipa<2*fNpart) {
826fd892 336
84555c93 337 GetQuarkPair(fFile, fgIntegral, yq[0], yq[1], ptq[0], ptq[1], dphi);
338
339 GetHadronPair(fFile, fQuark, yq[0], yq[1], ptq[0], ptq[1], ihadron[0], ihadron[1], plh[0], plh[1], pth[0], pth[1]);
340
9fd56238 341 if (fEnergy == 9 || fEnergy == -9) { // boost particles from c.m.s. to ALICE lab frame
342 Double_t dyBoost = 0.47;
343 Double_t beta = TMath::TanH(dyBoost);
344 Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
345 Double_t gb = gamma * beta;
346 yq[0] += dyBoost;
347 yq[1] += dyBoost;
348 plh[0] = gb * TMath::Sqrt(plh[0]*plh[0] + pth[0]*pth[0]) + gamma * plh[0];
349 plh[1] = gb * TMath::Sqrt(plh[1]*plh[1] + pth[1]*pth[1]) + gamma * plh[1];
350 if (fEnergy == 9) {
351 yq[0] *= -1;
352 yq[1] *= -1;
353 plh[0] *= -1;
354 plh[1] *= -1;
355 }
356 }
826fd892 357
84555c93 358 // Cuts from AliGenerator
359
360 // Cut on theta
361 theta=TMath::ATan2(pth[0],plh[0]);
b33adf51 362 if (theta<fThetaMin || theta>fThetaMax) continue;
84555c93 363 theta=TMath::ATan2(pth[1],plh[1]);
b33adf51 364 if (theta<fThetaMin || theta>fThetaMax) continue;
84555c93 365
366 // Cut on momentum
367 ph[0]=TMath::Sqrt(pth[0]*pth[0]+plh[0]*plh[0]);
368 if (ph[0]<fPMin || ph[0]>fPMax) continue;
369 ph[1]=TMath::Sqrt(pth[1]*pth[1]+plh[1]*plh[1]);
370 if (ph[1]<fPMin || ph[1]>fPMax) continue;
371
372 // Add the quarks in the stack
373
374 phiq[0] = Rndm()*k2PI;
375 if (Rndm() < 0.5) {
b44c3901 376 phiq[1] = phiq[0] + dphi*kDegrad;
84555c93 377 } else {
378 phiq[1] = phiq[0] - dphi*kDegrad;
379 }
380 if (phiq[1] > k2PI) phiq[1] -= k2PI;
381 if (phiq[1] < 0 ) phiq[1] += k2PI;
382
383 // quarks pdg
384 iquark[0] = +fQuark;
385 iquark[1] = -fQuark;
826fd892 386
84555c93 387 // px and py
388 TVector2 qvect1 = TVector2();
389 TVector2 qvect2 = TVector2();
390 qvect1.SetMagPhi(ptq[0],phiq[0]);
391 qvect2.SetMagPhi(ptq[1],phiq[1]);
392 pq[0][0] = qvect1.Px();
393 pq[0][1] = qvect1.Py();
394 pq[1][0] = qvect2.Px();
395 pq[1][1] = qvect2.Py();
826fd892 396
84555c93 397 // pz
cb11fd27 398 tanhy2 = TMath::TanH(yq[0]);
399 tanhy2 *= tanhy2;
400 pq[0][2] = TMath::Sqrt((ptq[0]*ptq[0]+qm*qm)*tanhy2/(1-tanhy2));
401 pq[0][2] = TMath::Sign((Double_t)pq[0][2],yq[0]);
402 tanhy2 = TMath::TanH(yq[1]);
403 tanhy2 *= tanhy2;
404 pq[1][2] = TMath::Sqrt((ptq[1]*ptq[1]+qm*qm)*tanhy2/(1-tanhy2));
405 pq[1][2] = TMath::Sign((Double_t)pq[1][2],yq[1]);
826fd892 406
84555c93 407 // Here we assume that |phi_H1 - phi_H2| = |phi_Q1 - phi_Q2| = dphi
408 // which is a good approximation for heavy flavors in Pythia
409 // ... moreover, same phi angles as for the quarks ...
410
411 phih[0] = phiq[0];
412 phih[1] = phiq[1];
826fd892 413
414 ipa1 = 0;
b44c3901 415
826fd892 416 for (ihad = 0; ihad < 2; ihad++) {
417 while(1) {
418
419 ipa0=ipa1;
84555c93 420
826fd892 421 // particle type
422 fChildWeight=(fDecayer->GetPartialBranchingRatio(ihadron[ihad]))*fParentWeight;
423 wgtp[ihad]=fParentWeight;
424 wgtch[ihad]=fChildWeight;
425 TParticlePDG *particle = pDataBase->GetParticle(ihadron[ihad]);
426 Float_t am = particle->Mass();
427 phi = phih[ihad];
428 pt = pth[ihad];
429 pl = plh[ihad];
430 ptot=TMath::Sqrt(pt*pt+pl*pl);
84555c93 431
826fd892 432 p[ihad][0]=pt*TMath::Cos(phi);
433 p[ihad][1]=pt*TMath::Sin(phi);
434 p[ihad][2]=pl;
435
436 if(fVertexSmear==kPerTrack) {
437 Rndm(random,6);
438 for (j=0;j<3;j++) {
439 origin0[ihad][j]=
440 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
441 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
442 }
84555c93 443 }
444
826fd892 445 // Looking at fForceDecay :
446 // if fForceDecay != none Primary particle decays using
447 // AliPythia and children are tracked by GEANT
448 //
449 // if fForceDecay == none Primary particle is tracked by GEANT
450 // (In the latest, make sure that GEANT actually does all the decays you want)
451
452 if (fForceDecay != kNoDecay) {
453 // Using lujet to decay particle
454 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
455 TLorentzVector pmom(p[ihad][0], p[ihad][1], p[ihad][2], energy);
456 fDecayer->Decay(ihadron[ihad],&pmom);
457
458 // select decay particles
459
460 np[ihad]=fDecayer->ImportParticles((TClonesArray *)particleslist->At(ihad));
461
462 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
463
464 if (fGeometryAcceptance)
465 if (!CheckAcceptanceGeometry(np[ihad],(TClonesArray*)particleslist->At(ihad))) continue;
466
467 trackIt[ihad] = new Int_t [np[ihad]];
468 pSelected[ihad] = new Int_t [np[ihad]];
469 Int_t* pFlag = new Int_t [np[ihad]];
470
471 for (i=0; i<np[ihad]; i++) {
472 pFlag[i] = 0;
473 pSelected[ihad][i] = 0;
474 }
475
476 if (np[ihad] >1) {
477 TParticle* iparticle = 0;
478 Int_t ipF, ipL;
84555c93 479
826fd892 480 for (i = 1; i<np[ihad] ; i++) {
481 trackIt[ihad][i] = 1;
482 iparticle =
483 (TParticle *) ((TClonesArray *) particleslist->At(ihad))->At(i);
484 Int_t kf = iparticle->GetPdgCode();
485 Int_t ks = iparticle->GetStatusCode();
486 // flagged particle
487 if (pFlag[i] == 1) {
84555c93 488 ipF = iparticle->GetFirstDaughter();
489 ipL = iparticle->GetLastDaughter();
490 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
826fd892 491 continue;
2c890605 492 }
493
826fd892 494 // flag decay products of particles with long life-time (c tau > .3 mum)
495 if (ks != 1) {
496 Double_t lifeTime = fDecayer->GetLifetime(kf);
497 if (lifeTime > (Double_t) fMaxLifeTime) {
498 ipF = iparticle->GetFirstDaughter();
499 ipL = iparticle->GetLastDaughter();
500 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
501 } else {
502 trackIt[ihad][i] = 0;
503 pSelected[ihad][i] = 1;
504 }
505 } // ks==1 ?
506 //
507 // children
508 if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[ihad][i])
509 {
510 if (fCutOnChild) {
511 pc[0]=iparticle->Px();
512 pc[1]=iparticle->Py();
513 pc[2]=iparticle->Pz();
514 //printf("px %f py %f pz %f\n",pc[0],pc[1],pc[2]);
515 Bool_t childok = KinematicSelection(iparticle, 1);
516 if(childok) {
517 pSelected[ihad][i] = 1;
518 ncsel[ihad]++;
519 } else {
520 ncsel[ihad]=-1;
521 break;
522 } // child kine cuts
523 } else {
524 pSelected[ihad][i] = 1;
525 ncsel[ihad]++;
526 } // if child selection
527 } // select muon
528 } // decay particle loop
529 } // if decay products
530
531 if ((fCutOnChild && ncsel[ihad] >0) || !fCutOnChild) ipa1++;
532
533 if (pFlag) delete[] pFlag;
534
535 } // kinematic selection
536 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
537 {
538 gAlice->GetMCApp()->
539 PushTrack(fTrackIt,-1,ihadron[ihad],p[ihad],origin0[ihad],polar[ihad],0,kPPrimary,nt,wgtp[ihad]);
540 ipa1++;
541 fNprimaries++;
542
543 }
544 break;
545 } // while(1) loop
546 if (ipa1<ipa0+1){
547 ipa1=0;
548 if (pSelected[ihad]) delete pSelected[ihad];
549 if (trackIt[ihad]) delete trackIt[ihad];
550 particleshad1->Clear();
551 particleshad2->Clear();
552 break;
553 }//go out of loop and generate new pair if at least one hadron is rejected
84555c93 554 } // hadron pair loop
826fd892 555 if(ipa1==2){
556
557 ipa=ipa+ipa1;
558
559 if(fForceDecay != kNoDecay){
560 for(ihad=0;ihad<2;ihad++){
561
562 //load tracks in the stack if both hadrons of the pair accepted
563 LoadTracks(iquark[ihad],pq[ihad],ihadron[ihad],p[ihad],np[ihad],
564 (TClonesArray *)particleslist->At(ihad),origin0[ihad],
565 polar[ihad],wgtp[ihad],wgtch[ihad],nt,ncsel[ihad],
566 pSelected[ihad],trackIt[ihad]);
567
568 if (pSelected[ihad]) delete pSelected[ihad];
569 if (trackIt[ihad]) delete trackIt[ihad];
9fd56238 570
826fd892 571 }
572 particleshad1->Clear();
573 particleshad2->Clear();
574 }
575 }
576 } // while (ipa<2*fNpart) loop
577
2c890605 578 SetHighWaterMark(nt);
84555c93 579
cb11fd27 580 AliGenEventHeader* header = new AliGenEventHeader("CorrHF");
84555c93 581 header->SetPrimaryVertex(fVertex);
582 header->SetNProduced(fNprimaries);
583 AddHeader(header);
826fd892 584
585
586 delete particleshad1;
587 delete particleshad2;
588 delete particleslist;
589
00fce64a 590 delete[] pSelected;
591 delete[] trackIt;
84555c93 592}
2c890605 593//____________________________________________________________________________________
9fd56238 594void AliGenCorrHF::IpCharm(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
2c890605 595{
9fd56238 596// Composition of a lower state charm hadron pair from a ccbar quark pair
597 Int_t pdgH[] = {411, 421, 431, 4122, 4132, 4232, 4332};
598
599 Double_t id3, id4;
600 hProbHH->GetRandom2(id3, id4);
601 pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
602 pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
603
604 return;
2c890605 605}
606
9fd56238 607void AliGenCorrHF::IpBeauty(TH2F *hProbHH, Int_t &pdg3, Int_t &pdg4)
2c890605 608{
9fd56238 609// Composition of a lower state beauty hadron pair from a bbbar quark pair
610 // B-Bbar mixing will be done by Pythia at their decay point
611 Int_t pdgH[] = {511, 521, 531, 5122, 5132, 5232, 5332};
612
613 Double_t id3, id4;
614 hProbHH->GetRandom2(id3, id4);
615 pdg3 = pdgH[(Int_t)TMath::Floor(id3)];
616 pdg4 = -1*pdgH[(Int_t)TMath::Floor(id4)];
617
618 if ( (pdg3== 511) || (pdg3== 521) || (pdg3== 531) ) pdg3 *= -1;
619 if ( (pdg4==-511) || (pdg4==-521) || (pdg4==-531) ) pdg4 *= -1;
620
621 return;
2c890605 622}
623
624//____________________________________________________________________________________
625Double_t AliGenCorrHF::ComputeIntegral(TFile* fG) // needed by GetQuarkPair
626{
627 // Read QQbar kinematical 5D grid's cell occupancy weights
4c111067 628 Int_t cell[6]; // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
2c890605 629 TTree* tG = (TTree*) fG->Get("tGqq");
4c111067 630 tG->GetBranch("cell")->SetAddress(&cell);
2c890605 631 Int_t nbins = tG->GetEntries();
632
633 // delete previously computed integral (if any)
634 if(fgIntegral) delete [] fgIntegral;
635
636 fgIntegral = new Double_t[nbins+1];
637 fgIntegral[0] = 0;
638 Int_t bin;
639 for(bin=0;bin<nbins;bin++) {
640 tG->GetEvent(bin);
641 fgIntegral[bin+1] = fgIntegral[bin] + cell[0];
642 }
643 // Normalize integral to 1
644 if (fgIntegral[nbins] == 0 ) {
645 return 0;
646 }
647 for (bin=1;bin<=nbins;bin++) fgIntegral[bin] /= fgIntegral[nbins];
648
649 return fgIntegral[nbins];
650}
651
826fd892 652
2c890605 653//____________________________________________________________________________________
654void AliGenCorrHF::GetQuarkPair(TFile* fG, Double_t* fInt, Double_t &y1, Double_t &y2, Double_t &pt1, Double_t &pt2, Double_t &dphi)
655 // modification of ROOT's TH3::GetRandom3 for 5D
656{
657 // Read QQbar kinematical 5D grid's cell coordinates
4c111067 658 Int_t cell[6]; // cell[6]={wght,iy1,iy2,ipt1,ipt2,idph}
2c890605 659 TTree* tG = (TTree*) fG->Get("tGqq");
4c111067 660 tG->GetBranch("cell")->SetAddress(&cell);
2c890605 661 Int_t nbins = tG->GetEntries();
662 Double_t rand[6];
663 gRandom->RndmArray(6,rand);
664 Int_t ibin = TMath::BinarySearch(nbins,fInt,rand[0]);
665 tG->GetEvent(ibin);
666 y1 = fgy[cell[1]] + (fgy[cell[1]+1]-fgy[cell[1]])*rand[1];
667 y2 = fgy[cell[2]] + (fgy[cell[2]+1]-fgy[cell[2]])*rand[2];
668 pt1 = fgpt[cell[3]] + (fgpt[cell[3]+1]-fgpt[cell[3]])*rand[3];
669 pt2 = fgpt[cell[4]] + (fgpt[cell[4]+1]-fgpt[cell[4]])*rand[4];
670 dphi = fgdph[cell[5]]+ (fgdph[cell[5]+1]-fgdph[cell[5]])*rand[5];
671}
672
673//____________________________________________________________________________________
674void 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)
675{
676 // Generate a hadron pair
9fd56238 677 void (*fIpParaFunc)(TH2F *, Int_t &, Int_t &);//Pointer to hadron pair composition function
2c890605 678 fIpParaFunc = IpCharm;
679 Double_t mq = 1.2; // c & b quark masses (used in AliPythia)
680 if (idq == 5) {
681 fIpParaFunc = IpBeauty;
682 mq = 4.75;
683 }
08bffa4d 684 Double_t z11 = 0.;
685 Double_t z12 = 0.;
686 Double_t z21 = 0.;
687 Double_t z22 = 0.;
688 Double_t pz1, pz2, e1, e2, mh, ptemp, rand[2];
2c890605 689 char tag[100];
9fd56238 690 TH2F *h2h[12], *h2s[12], *hProbHH; // hard & soft fragmentation and HH-probability functions
2c890605 691 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
09298f0c 692 snprintf(tag,100, "h2h_pt%d",ipt);
2c890605 693 h2h[ipt] = (TH2F*) fG->Get(tag);
09298f0c 694 snprintf(tag,100, "h2s_pt%d",ipt);
2c890605 695 h2s[ipt] = (TH2F*) fG->Get(tag);
696 }
697
698 if (y1*y2 < 0) {
699 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
700 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
701 h2h[ipt]->GetRandom2(z11, z21);
702 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
703 h2h[ipt]->GetRandom2(z12, z22);
704 }
705 }
706 else {
707 if (TMath::Abs(y1) > TMath::Abs(y2)) {
708 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
709 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
710 h2h[ipt]->GetRandom2(z11, z21);
711 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
712 h2s[ipt]->GetRandom2(z12, z22);
713 }
714 }
715 else {
716 for (Int_t ipt = 0; ipt<fgnptbins; ipt++) {
717 if(pt1 >= fgptbmin[ipt] && pt1 < fgptbmax[ipt])
718 h2s[ipt]->GetRandom2(z11, z21);
719 if(pt2 >= fgptbmin[ipt] && pt2 < fgptbmax[ipt])
720 h2h[ipt]->GetRandom2(z12, z22);
721 }
722 }
723 }
724 gRandom->RndmArray(2,rand);
725 ptemp = TMath::Sqrt(pt1*pt1 + mq*mq);
726 pz1 = ptemp*TMath::SinH(y1);
727 e1 = ptemp*TMath::CosH(y1);
728 ptemp = TMath::Sqrt(pt2*pt2 + mq*mq);
729 pz2 = ptemp*TMath::SinH(y2);
730 e2 = ptemp*TMath::CosH(y2);
731
9fd56238 732 hProbHH = (TH2F*)fG->Get("hProbHH");
733 fIpParaFunc(hProbHH, id3, id4);
2c890605 734 mh = TDatabasePDG::Instance()->GetParticle(id3)->Mass();
735 ptemp = z11*z21*(e1*e1-pz1*pz1) - mh*mh;
b33adf51 736 if (idq==5) pt3 = pt1; // an approximation at low pt, try better
9fd56238 737 else pt3 = rand[0]; // pt3=pt1 gives less D-hadrons at low pt
2c890605 738 if (ptemp > 0) pt3 = TMath::Sqrt(ptemp);
739 if (pz1 > 0) pz3 = (z11*(e1 + pz1) - z21*(e1 - pz1)) / 2;
740 else pz3 = (z21*(e1 + pz1) - z11*(e1 - pz1)) / 2;
741 e1 = TMath::Sqrt(pz3*pz3 + pt3*pt3 + mh*mh);
742
2c890605 743 mh = TDatabasePDG::Instance()->GetParticle(id4)->Mass();
744 ptemp = z12*z22*(e2*e2-pz2*pz2) - mh*mh;
b33adf51 745 if (idq==5) pt4 = pt2; // an approximation at low pt, try better
746 else pt4 = rand[1];
2c890605 747 if (ptemp > 0) pt4 = TMath::Sqrt(ptemp);
748 if (pz2 > 0) pz4 = (z12*(e2 + pz2) - z22*(e2 - pz2)) / 2;
749 else pz4 = (z22*(e2 + pz2) - z12*(e2 - pz2)) / 2;
750 e2 = TMath::Sqrt(pz4*pz4 + pt4*pt4 + mh*mh);
751
752 // small corr. instead of using Frag. Func. depending on yQ (in addition to ptQ)
753 Float_t ycorr = 0.2, y3, y4;
754 gRandom->RndmArray(2,rand);
755 y3 = 0.5 * TMath::Log((e1 + pz3 + 1.e-13)/(e1 - pz3 + 1.e-13));
756 y4 = 0.5 * TMath::Log((e2 + pz4 + 1.e-13)/(e2 - pz4 + 1.e-13));
757 if(TMath::Abs(y3)<ycorr && TMath::Abs(y4)<ycorr && rand[0]>0.5) {
bfd20868 758 ptemp = TMath::Sqrt((e1-pz3)*(e1+pz3));
2c890605 759 y3 = 4*(1 - 2*rand[1]);
760 pz3 = ptemp*TMath::SinH(y3);
761 pz4 = pz3;
762 }
763}
5cc3bcd2 764
826fd892 765
766//____________________________________________________________________________________
767void AliGenCorrHF::LoadTracks(Int_t iquark, Float_t *pq,
768 Int_t iPart, Float_t *p,
769 Int_t np, TClonesArray *particles,
770 Float_t *origin0, Float_t *polar,
771 Float_t wgtp, Float_t wgtch,
772 Int_t &nt, Int_t ncsel, Int_t *pSelected,
773 Int_t *trackIt){
774 Int_t i;
775 Int_t ntq=-1;
776 Int_t* pParent = new Int_t[np];
777 Float_t pc[3], och[3];
778 Int_t iparent;
779
780 for(i=0;i<np;i++) pParent[i]=-1;
781
782 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
783 // Parents
784 // quark
785 PushTrack(0, -1, iquark, pq, origin0, polar, 0, kPPrimary, nt, wgtp);
786 KeepTrack(nt);
787 ntq = nt;
788 // hadron
789 PushTrack(0, ntq, iPart, p, origin0, polar, 0, kPDecay, nt, wgtp);
790 pParent[0] = nt;
791 KeepTrack(nt);
792 fNprimaries++;
793
794 // Decay Products
795 for (i = 1; i < np; i++) {
796 if (pSelected[i]) {
797
798 TParticle* iparticle = (TParticle *) particles->At(i);
799 Int_t kf = iparticle->GetPdgCode();
800 Int_t jpa = iparticle->GetFirstMother()-1;
801
802 och[0] = origin0[0]+iparticle->Vx()/10;
803 och[1] = origin0[1]+iparticle->Vy()/10;
804 och[2] = origin0[2]+iparticle->Vz()/10;
805 pc[0] = iparticle->Px();
806 pc[1] = iparticle->Py();
807 pc[2] = iparticle->Pz();
808
809 if (jpa > -1) {
810 iparent = pParent[jpa];
811 } else {
812 iparent = -1;
813 }
814
815 PushTrack(fTrackIt*trackIt[i], iparent, kf,
816 pc, och, polar,
817 0, kPDecay, nt, wgtch);
818 pParent[i] = nt;
819 KeepTrack(nt);
820 fNprimaries++;
821
822 } // Selected
823 } // Particle loop
824 }
825 if (pParent) delete[] pParent;
826
827 return;
828}
829