]>
Commit | Line | Data |
---|---|---|
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 | // | |
19 | // Generator using the TPythia interface (via AliPythia) | |
20 | // to generate pp collisions. | |
21 | // Using SetNuclei() also nuclear modifications to the structure functions | |
22 | // can be taken into account. This makes, of course, only sense for the | |
23 | // generation of the products of hard processes (heavy flavor, jets ...) | |
24 | // | |
25 | // andreas.morsch@cern.ch | |
26 | // | |
27 | ||
28 | #include <TDatabasePDG.h> | |
29 | #include <TParticle.h> | |
30 | #include <TPDGCode.h> | |
31 | #include <TSystem.h> | |
32 | #include <TTree.h> | |
33 | ||
34 | #include "AliConst.h" | |
35 | #include "AliDecayerPythia.h" | |
36 | #include "AliGenPythia.h" | |
37 | #include "AliHeader.h" | |
38 | #include "AliGenPythiaEventHeader.h" | |
39 | #include "AliPythia.h" | |
40 | #include "AliPythiaRndm.h" | |
41 | #include "AliRun.h" | |
42 | #include "AliStack.h" | |
43 | #include "AliRunLoader.h" | |
44 | #include "AliMC.h" | |
45 | ||
46 | ClassImp(AliGenPythia) | |
47 | ||
48 | AliGenPythia::AliGenPythia() | |
49 | :AliGenMC() | |
50 | { | |
51 | // Default Constructor | |
52 | fParticles = 0; | |
53 | fPythia = 0; | |
54 | fHeader = 0; | |
55 | fReadFromFile = 0; | |
56 | fDecayer = new AliDecayerPythia(); | |
57 | SetEventListRange(); | |
58 | SetJetPhiRange(); | |
59 | SetJetEtaRange(); | |
60 | SetJetEtRange(); | |
61 | SetGammaPhiRange(); | |
62 | SetGammaEtaRange(); | |
63 | SetPtKick(); | |
64 | SetQuench(); | |
65 | SetHadronisation(); | |
66 | fSetNuclei = kFALSE; | |
67 | if (!AliPythiaRndm::GetPythiaRandom()) | |
68 | AliPythiaRndm::SetPythiaRandom(GetRandom()); | |
69 | } | |
70 | ||
71 | AliGenPythia::AliGenPythia(Int_t npart) | |
72 | :AliGenMC(npart) | |
73 | { | |
74 | // default charm production at 5. 5 TeV | |
75 | // semimuonic decay | |
76 | // structure function GRVHO | |
77 | // | |
78 | fName = "Pythia"; | |
79 | fTitle= "Particle Generator using PYTHIA"; | |
80 | fXsection = 0.; | |
81 | fReadFromFile = 0; | |
82 | SetProcess(); | |
83 | SetStrucFunc(); | |
84 | SetForceDecay(); | |
85 | SetPtHard(); | |
86 | SetYHard(); | |
87 | SetEnergyCMS(); | |
88 | fDecayer = new AliDecayerPythia(); | |
89 | // Set random number generator | |
90 | if (!AliPythiaRndm::GetPythiaRandom()) | |
91 | AliPythiaRndm::SetPythiaRandom(GetRandom()); | |
92 | fFlavorSelect = 0; | |
93 | // Produced particles | |
94 | fParticles = new TClonesArray("TParticle",1000); | |
95 | fHeader = 0; | |
96 | SetEventListRange(); | |
97 | SetJetPhiRange(); | |
98 | SetJetEtaRange(); | |
99 | SetJetEtRange(); | |
100 | SetGammaPhiRange(); | |
101 | SetGammaEtaRange(); | |
102 | SetJetReconstructionMode(); | |
103 | SetQuench(); | |
104 | SetHadronisation(); | |
105 | SetPtKick(); | |
106 | // Options determining what to keep in the stack (Heavy flavour generation) | |
107 | fStackFillOpt = kFlavorSelection; // Keep particle with selected flavor | |
108 | fFeedDownOpt = kTRUE; // allow feed down from higher family | |
109 | // Fragmentation on/off | |
110 | fFragmentation = kTRUE; | |
111 | // Default counting mode | |
112 | fCountMode = kCountAll; | |
113 | // Pycel | |
114 | SetPycellParameters(); | |
115 | fSetNuclei = kFALSE; | |
116 | } | |
117 | ||
118 | AliGenPythia::AliGenPythia(const AliGenPythia & Pythia) | |
119 | :AliGenMC(Pythia) | |
120 | { | |
121 | // copy constructor | |
122 | Pythia.Copy(*this); | |
123 | } | |
124 | ||
125 | AliGenPythia::~AliGenPythia() | |
126 | { | |
127 | // Destructor | |
128 | } | |
129 | ||
130 | void AliGenPythia::SetPycellParameters(Float_t etamax, Int_t neta, Int_t nphi, | |
131 | Float_t thresh, Float_t etseed, Float_t minet, Float_t r) | |
132 | { | |
133 | // Set pycell parameters | |
134 | fPycellEtaMax = etamax; | |
135 | fPycellNEta = neta; | |
136 | fPycellNPhi = nphi; | |
137 | fPycellThreshold = thresh; | |
138 | fPycellEtSeed = etseed; | |
139 | fPycellMinEtJet = minet; | |
140 | fPycellMaxRadius = r; | |
141 | } | |
142 | ||
143 | ||
144 | ||
145 | void AliGenPythia::SetEventListRange(Int_t eventFirst, Int_t eventLast) | |
146 | { | |
147 | // Set a range of event numbers, for which a table | |
148 | // of generated particle will be printed | |
149 | fDebugEventFirst = eventFirst; | |
150 | fDebugEventLast = eventLast; | |
151 | if (fDebugEventLast==-1) fDebugEventLast=fDebugEventFirst; | |
152 | } | |
153 | ||
154 | void AliGenPythia::Init() | |
155 | { | |
156 | // Initialisation | |
157 | ||
158 | SetMC(AliPythia::Instance()); | |
159 | fPythia=(AliPythia*) fMCEvGen; | |
160 | ||
161 | // | |
162 | fParentWeight=1./Float_t(fNpart); | |
163 | // | |
164 | // Forward Paramters to the AliPythia object | |
165 | fDecayer->SetForceDecay(fForceDecay); | |
166 | fDecayer->Init(); | |
167 | ||
168 | ||
169 | fPythia->SetCKIN(3,fPtHardMin); | |
170 | fPythia->SetCKIN(4,fPtHardMax); | |
171 | fPythia->SetCKIN(7,fYHardMin); | |
172 | fPythia->SetCKIN(8,fYHardMax); | |
173 | ||
174 | if (fAProjectile > 0 && fATarget > 0) fPythia->SetNuclei(fAProjectile, fATarget); | |
175 | // Fragmentation? | |
176 | if (fFragmentation) { | |
177 | fPythia->SetMSTP(111,1); | |
178 | } else { | |
179 | fPythia->SetMSTP(111,0); | |
180 | } | |
181 | ||
182 | ||
183 | // initial state radiation | |
184 | fPythia->SetMSTP(61,fGinit); | |
185 | // final state radiation | |
186 | fPythia->SetMSTP(71,fGfinal); | |
187 | // pt - kick | |
188 | if (fPtKick > 0.) { | |
189 | fPythia->SetMSTP(91,1); | |
190 | fPythia->SetPARP(91,fPtKick); | |
191 | } else { | |
192 | fPythia->SetMSTP(91,0); | |
193 | } | |
194 | ||
195 | ||
196 | if (fReadFromFile) { | |
197 | fRL = AliRunLoader::Open(fFileName, "Partons"); | |
198 | fRL->LoadKinematics(); | |
199 | fRL->LoadHeader(); | |
200 | } else { | |
201 | fRL = 0x0; | |
202 | } | |
203 | ||
204 | ||
205 | // | |
206 | fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc); | |
207 | ||
208 | // Parent and Children Selection | |
209 | switch (fProcess) | |
210 | { | |
211 | case kPyCharm: | |
212 | case kPyCharmUnforced: | |
213 | case kPyCharmPbPbMNR: | |
214 | case kPyCharmppMNR: | |
215 | case kPyCharmpPbMNR: | |
216 | fParentSelect[0] = 411; | |
217 | fParentSelect[1] = 421; | |
218 | fParentSelect[2] = 431; | |
219 | fParentSelect[3] = 4122; | |
220 | fFlavorSelect = 4; | |
221 | break; | |
222 | case kPyD0PbPbMNR: | |
223 | case kPyD0pPbMNR: | |
224 | case kPyD0ppMNR: | |
225 | fParentSelect[0] = 421; | |
226 | fFlavorSelect = 4; | |
227 | break; | |
228 | case kPyDPlusPbPbMNR: | |
229 | case kPyDPluspPbMNR: | |
230 | case kPyDPlusppMNR: | |
231 | fParentSelect[0] = 411; | |
232 | fFlavorSelect = 4; | |
233 | break; | |
234 | case kPyBeauty: | |
235 | case kPyBeautyPbPbMNR: | |
236 | case kPyBeautypPbMNR: | |
237 | case kPyBeautyppMNR: | |
238 | fParentSelect[0]= 511; | |
239 | fParentSelect[1]= 521; | |
240 | fParentSelect[2]= 531; | |
241 | fParentSelect[3]= 5122; | |
242 | fParentSelect[4]= 5132; | |
243 | fParentSelect[5]= 5232; | |
244 | fParentSelect[6]= 5332; | |
245 | fFlavorSelect = 5; | |
246 | break; | |
247 | case kPyBeautyUnforced: | |
248 | fParentSelect[0] = 511; | |
249 | fParentSelect[1] = 521; | |
250 | fParentSelect[2] = 531; | |
251 | fParentSelect[3] = 5122; | |
252 | fParentSelect[4] = 5132; | |
253 | fParentSelect[5] = 5232; | |
254 | fParentSelect[6] = 5332; | |
255 | fFlavorSelect = 5; | |
256 | break; | |
257 | case kPyJpsiChi: | |
258 | case kPyJpsi: | |
259 | fParentSelect[0] = 443; | |
260 | break; | |
261 | case kPyMb: | |
262 | case kPyMbNonDiffr: | |
263 | case kPyJets: | |
264 | case kPyDirectGamma: | |
265 | break; | |
266 | case kPyW: | |
267 | break; | |
268 | } | |
269 | // | |
270 | // | |
271 | // JetFinder for Trigger | |
272 | // | |
273 | // Configure detector (EMCAL like) | |
274 | // | |
275 | fPythia->SetPARU(51, fPycellEtaMax); | |
276 | fPythia->SetMSTU(51, fPycellNEta); | |
277 | fPythia->SetMSTU(52, fPycellNPhi); | |
278 | // | |
279 | // Configure Jet Finder | |
280 | // | |
281 | fPythia->SetPARU(58, fPycellThreshold); | |
282 | fPythia->SetPARU(52, fPycellEtSeed); | |
283 | fPythia->SetPARU(53, fPycellMinEtJet); | |
284 | fPythia->SetPARU(54, fPycellMaxRadius); | |
285 | fPythia->SetMSTU(54, 2); | |
286 | // | |
287 | // This counts the total number of calls to Pyevnt() per run. | |
288 | fTrialsRun = 0; | |
289 | fQ = 0.; | |
290 | fX1 = 0.; | |
291 | fX2 = 0.; | |
292 | fNev = 0 ; | |
293 | // | |
294 | // | |
295 | // | |
296 | AliGenMC::Init(); | |
297 | // | |
298 | // | |
299 | // | |
300 | if (fSetNuclei) { | |
301 | fDyBoost = 0; | |
302 | Warning("Init","SetNuclei used. Use SetProjectile + SetTarget instead. fDyBoost has been reset to 0\n"); | |
303 | } | |
304 | ||
305 | if (fQuench) { | |
306 | fPythia->InitQuenching(0., 0.1, 0.6e6, 0); | |
307 | } | |
308 | ||
309 | } | |
310 | ||
311 | void AliGenPythia::Generate() | |
312 | { | |
313 | // Generate one event | |
314 | ||
315 | fDecayer->ForceDecay(); | |
316 | ||
317 | Float_t polar[3] = {0,0,0}; | |
318 | Float_t origin[3] = {0,0,0}; | |
319 | Float_t p[4]; | |
320 | // converts from mm/c to s | |
321 | const Float_t kconv=0.001/2.999792458e8; | |
322 | // | |
323 | Int_t nt=0; | |
324 | Int_t jev=0; | |
325 | Int_t j, kf; | |
326 | fTrials=0; | |
327 | ||
328 | ||
329 | // Set collision vertex position | |
330 | if (fVertexSmear == kPerEvent) Vertex(); | |
331 | ||
332 | // event loop | |
333 | while(1) | |
334 | { | |
335 | // | |
336 | // Produce event | |
337 | // | |
338 | // | |
339 | // Switch hadronisation off | |
340 | // | |
341 | fPythia->SetMSTJ(1, 0); | |
342 | // | |
343 | // Either produce new event or read partons from file | |
344 | // | |
345 | if (!fReadFromFile) { | |
346 | fPythia->Pyevnt(); | |
347 | fNpartons = fPythia->GetN(); | |
348 | } else { | |
349 | printf("Loading Event %d\n",AliRunLoader::GetRunLoader()->GetEventNumber()); | |
350 | fRL->GetEvent(AliRunLoader::GetRunLoader()->GetEventNumber()); | |
351 | fPythia->SetN(0); | |
352 | LoadEvent(fRL->Stack(), 0 , 1); | |
353 | fPythia->Pyedit(21); | |
354 | } | |
355 | ||
356 | // | |
357 | // Run quenching routine | |
358 | // | |
359 | if (fQuench == 1) { | |
360 | fPythia->Quench(); | |
361 | } else if (fQuench == 2){ | |
362 | fPythia->Pyquen(208., 0, 0.); | |
363 | } | |
364 | // | |
365 | // Switch hadronisation on | |
366 | // | |
367 | fPythia->SetMSTJ(1, 1); | |
368 | // | |
369 | // .. and perform hadronisation | |
370 | // printf("Calling hadronisation %d\n", fPythia->GetN()); | |
371 | fPythia->Pyexec(); | |
372 | fTrials++; | |
373 | fPythia->ImportParticles(fParticles,"All"); | |
374 | Boost(); | |
375 | // | |
376 | // | |
377 | // | |
378 | Int_t i; | |
379 | ||
380 | ||
381 | Int_t np = fParticles->GetEntriesFast(); | |
382 | ||
383 | if (np == 0 ) continue; | |
384 | // | |
385 | ||
386 | // | |
387 | Int_t* pParent = new Int_t[np]; | |
388 | Int_t* pSelected = new Int_t[np]; | |
389 | Int_t* trackIt = new Int_t[np]; | |
390 | for (i = 0; i < np; i++) { | |
391 | pParent[i] = -1; | |
392 | pSelected[i] = 0; | |
393 | trackIt[i] = 0; | |
394 | } | |
395 | ||
396 | Int_t nc = 0; // Total n. of selected particles | |
397 | Int_t nParents = 0; // Selected parents | |
398 | Int_t nTkbles = 0; // Trackable particles | |
399 | if (fProcess != kPyMb && fProcess != kPyJets && | |
400 | fProcess != kPyDirectGamma && | |
401 | fProcess != kPyMbNonDiffr && | |
402 | fProcess != kPyW) { | |
403 | ||
404 | for (i = 0; i < np; i++) { | |
405 | TParticle* iparticle = (TParticle *) fParticles->At(i); | |
406 | Int_t ks = iparticle->GetStatusCode(); | |
407 | kf = CheckPDGCode(iparticle->GetPdgCode()); | |
408 | // No initial state partons | |
409 | if (ks==21) continue; | |
410 | // | |
411 | // Heavy Flavor Selection | |
412 | // | |
413 | // quark ? | |
414 | kf = TMath::Abs(kf); | |
415 | Int_t kfl = kf; | |
416 | // Resonance | |
417 | if (kfl > 100000) kfl %= 100000; | |
418 | // meson ? | |
419 | if (kfl > 10) kfl/=100; | |
420 | // baryon | |
421 | if (kfl > 10) kfl/=10; | |
422 | Int_t ipa = iparticle->GetFirstMother()-1; | |
423 | Int_t kfMo = 0; | |
424 | ||
425 | if (ipa > -1) { | |
426 | TParticle * mother = (TParticle *) fParticles->At(ipa); | |
427 | kfMo = TMath::Abs(mother->GetPdgCode()); | |
428 | } | |
429 | // What to keep in Stack? | |
430 | Bool_t flavorOK = kFALSE; | |
431 | Bool_t selectOK = kFALSE; | |
432 | if (fFeedDownOpt) { | |
433 | if (kfl >= fFlavorSelect) flavorOK = kTRUE; | |
434 | } else { | |
435 | if (kfl > fFlavorSelect) { | |
436 | nc = -1; | |
437 | break; | |
438 | } | |
439 | if (kfl == fFlavorSelect) flavorOK = kTRUE; | |
440 | } | |
441 | switch (fStackFillOpt) { | |
442 | case kFlavorSelection: | |
443 | selectOK = kTRUE; | |
444 | break; | |
445 | case kParentSelection: | |
446 | if (ParentSelected(kf) || kf <= 10) selectOK = kTRUE; | |
447 | break; | |
448 | } | |
449 | if (flavorOK && selectOK) { | |
450 | // | |
451 | // Heavy flavor hadron or quark | |
452 | // | |
453 | // Kinematic seletion on final state heavy flavor mesons | |
454 | if (ParentSelected(kf) && !KinematicSelection(iparticle, 0)) | |
455 | { | |
456 | continue; | |
457 | } | |
458 | pSelected[i] = 1; | |
459 | if (ParentSelected(kf)) ++nParents; // Update parent count | |
460 | // printf("\n particle (HF) %d %d %d", i, pSelected[i], kf); | |
461 | } else { | |
462 | // Kinematic seletion on decay products | |
463 | if (fCutOnChild && ParentSelected(kfMo) && ChildSelected(kf) | |
464 | && !KinematicSelection(iparticle, 1)) | |
465 | { | |
466 | continue; | |
467 | } | |
468 | // | |
469 | // Decay products | |
470 | // Select if mother was selected and is not tracked | |
471 | ||
472 | if (pSelected[ipa] && | |
473 | !trackIt[ipa] && // mother will be tracked ? | |
474 | kfMo != 5 && // mother is b-quark, don't store fragments | |
475 | kfMo != 4 && // mother is c-quark, don't store fragments | |
476 | kf != 92) // don't store string | |
477 | { | |
478 | // | |
479 | // Semi-stable or de-selected: diselect decay products: | |
480 | // | |
481 | // | |
482 | if (pSelected[i] == -1 || fDecayer->GetLifetime(kf) > fMaxLifeTime) | |
483 | { | |
484 | Int_t ipF = iparticle->GetFirstDaughter(); | |
485 | Int_t ipL = iparticle->GetLastDaughter(); | |
486 | if (ipF > 0) for (j = ipF-1; j < ipL; j++) pSelected[j] = -1; | |
487 | } | |
488 | // printf("\n particle (decay) %d %d %d", i, pSelected[i], kf); | |
489 | pSelected[i] = (pSelected[i] == -1) ? 0 : 1; | |
490 | } | |
491 | } | |
492 | if (pSelected[i] == -1) pSelected[i] = 0; | |
493 | if (!pSelected[i]) continue; | |
494 | // Count quarks only if you did not include fragmentation | |
495 | if (fFragmentation && kf <= 10) continue; | |
496 | ||
497 | nc++; | |
498 | // Decision on tracking | |
499 | trackIt[i] = 0; | |
500 | // | |
501 | // Track final state particle | |
502 | if (ks == 1) trackIt[i] = 1; | |
503 | // Track semi-stable particles | |
504 | if ((ks == 1) || (fDecayer->GetLifetime(kf) > fMaxLifeTime)) trackIt[i] = 1; | |
505 | // Track particles selected by process if undecayed. | |
506 | if (fForceDecay == kNoDecay) { | |
507 | if (ParentSelected(kf)) trackIt[i] = 1; | |
508 | } else { | |
509 | if (ParentSelected(kf)) trackIt[i] = 0; | |
510 | } | |
511 | if (trackIt[i] == 1) ++nTkbles; // Update trackable counter | |
512 | // | |
513 | // | |
514 | ||
515 | } // particle selection loop | |
516 | if (nc > 0) { | |
517 | for (i = 0; i<np; i++) { | |
518 | if (!pSelected[i]) continue; | |
519 | TParticle * iparticle = (TParticle *) fParticles->At(i); | |
520 | kf = CheckPDGCode(iparticle->GetPdgCode()); | |
521 | Int_t ks = iparticle->GetStatusCode(); | |
522 | p[0] = iparticle->Px(); | |
523 | p[1] = iparticle->Py(); | |
524 | p[2] = iparticle->Pz(); | |
525 | p[3] = iparticle->Energy(); | |
526 | ||
527 | origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm] | |
528 | origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm] | |
529 | origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm] | |
530 | ||
531 | Float_t tof = kconv*iparticle->T(); | |
532 | Int_t ipa = iparticle->GetFirstMother()-1; | |
533 | Int_t iparent = (ipa > -1) ? pParent[ipa] : -1; | |
534 | ||
535 | PushTrack(fTrackIt*trackIt[i], iparent, kf, | |
536 | p[0], p[1], p[2], p[3], | |
537 | origin[0], origin[1], origin[2], tof, | |
538 | polar[0], polar[1], polar[2], | |
539 | kPPrimary, nt, 1., ks); | |
540 | pParent[i] = nt; | |
541 | KeepTrack(nt); | |
542 | } // PushTrack loop | |
543 | } | |
544 | } else { | |
545 | nc = GenerateMB(); | |
546 | } // mb ? | |
547 | ||
548 | if (pParent) delete[] pParent; | |
549 | if (pSelected) delete[] pSelected; | |
550 | if (trackIt) delete[] trackIt; | |
551 | ||
552 | if (nc > 0) { | |
553 | switch (fCountMode) { | |
554 | case kCountAll: | |
555 | // printf(" Count all \n"); | |
556 | jev += nc; | |
557 | break; | |
558 | case kCountParents: | |
559 | // printf(" Count parents \n"); | |
560 | jev += nParents; | |
561 | break; | |
562 | case kCountTrackables: | |
563 | // printf(" Count trackable \n"); | |
564 | jev += nTkbles; | |
565 | break; | |
566 | } | |
567 | if (jev >= fNpart || fNpart == -1) { | |
568 | fKineBias=Float_t(fNpart)/Float_t(fTrials); | |
569 | printf("\n Trials: %i %i %i\n",fTrials, fNpart, jev); | |
570 | ||
571 | fQ += fPythia->GetVINT(51); | |
572 | fX1 += fPythia->GetVINT(41); | |
573 | fX2 += fPythia->GetVINT(42); | |
574 | fTrialsRun += fTrials; | |
575 | fNev++; | |
576 | MakeHeader(); | |
577 | break; | |
578 | } | |
579 | } | |
580 | } // event loop | |
581 | SetHighWaterMark(nt); | |
582 | // adjust weight due to kinematic selection | |
583 | // AdjustWeights(); | |
584 | // get cross-section | |
585 | fXsection=fPythia->GetPARI(1); | |
586 | } | |
587 | ||
588 | Int_t AliGenPythia::GenerateMB() | |
589 | { | |
590 | // | |
591 | // Min Bias selection and other global selections | |
592 | // | |
593 | Int_t i, kf, nt, iparent; | |
594 | Int_t nc = 0; | |
595 | Float_t p[4]; | |
596 | Float_t polar[3] = {0,0,0}; | |
597 | Float_t origin[3] = {0,0,0}; | |
598 | // converts from mm/c to s | |
599 | const Float_t kconv=0.001/2.999792458e8; | |
600 | ||
601 | ||
602 | ||
603 | Int_t np = (fHadronisation) ? fParticles->GetEntriesFast() : fNpartons; | |
604 | ||
605 | ||
606 | Int_t* pParent = new Int_t[np]; | |
607 | for (i=0; i< np; i++) pParent[i] = -1; | |
608 | if (fProcess == kPyJets || fProcess == kPyDirectGamma) { | |
609 | TParticle* jet1 = (TParticle *) fParticles->At(6); | |
610 | TParticle* jet2 = (TParticle *) fParticles->At(7); | |
611 | if (!CheckTrigger(jet1, jet2)) return 0; | |
612 | } | |
613 | ||
614 | for (i = 0; i<np; i++) { | |
615 | Int_t trackIt = 0; | |
616 | TParticle * iparticle = (TParticle *) fParticles->At(i); | |
617 | kf = CheckPDGCode(iparticle->GetPdgCode()); | |
618 | Int_t ks = iparticle->GetStatusCode(); | |
619 | Int_t km = iparticle->GetFirstMother(); | |
620 | if ((ks == 1 && kf!=0 && KinematicSelection(iparticle, 0)) || | |
621 | (ks != 1) || | |
622 | (fProcess == kPyJets && ks == 21 && km == 0 && i>1)) { | |
623 | nc++; | |
624 | if (ks == 1) trackIt = 1; | |
625 | Int_t ipa = iparticle->GetFirstMother()-1; | |
626 | ||
627 | iparent = (ipa > -1) ? pParent[ipa] : -1; | |
628 | ||
629 | // | |
630 | // store track information | |
631 | p[0] = iparticle->Px(); | |
632 | p[1] = iparticle->Py(); | |
633 | p[2] = iparticle->Pz(); | |
634 | p[3] = iparticle->Energy(); | |
635 | ||
636 | ||
637 | origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm] | |
638 | origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm] | |
639 | origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm] | |
640 | ||
641 | Float_t tof=kconv*iparticle->T(); | |
642 | ||
643 | PushTrack(fTrackIt*trackIt, iparent, kf, | |
644 | p[0], p[1], p[2], p[3], | |
645 | origin[0], origin[1], origin[2], tof, | |
646 | polar[0], polar[1], polar[2], | |
647 | kPPrimary, nt, 1., ks); | |
648 | // | |
649 | // Special Treatment to store color-flow | |
650 | // | |
651 | if (ks == 3 || ks == 13 || ks == 14) { | |
652 | TParticle* particle = 0; | |
653 | if (fStack) { | |
654 | particle = fStack->Particle(nt); | |
655 | } else { | |
656 | particle = gAlice->Stack()->Particle(nt); | |
657 | } | |
658 | particle->SetFirstDaughter(fPythia->GetK(2, i)); | |
659 | particle->SetLastDaughter(fPythia->GetK(3, i)); | |
660 | } | |
661 | ||
662 | KeepTrack(nt); | |
663 | pParent[i] = nt; | |
664 | } // select particle | |
665 | } // particle loop | |
666 | ||
667 | if (pParent) delete[] pParent; | |
668 | ||
669 | printf("\n I've put %i particles on the stack \n",nc); | |
670 | return nc; | |
671 | } | |
672 | ||
673 | ||
674 | void AliGenPythia::FinishRun() | |
675 | { | |
676 | // Print x-section summary | |
677 | fPythia->Pystat(1); | |
678 | ||
679 | if (fNev > 0.) { | |
680 | fQ /= fNev; | |
681 | fX1 /= fNev; | |
682 | fX2 /= fNev; | |
683 | } | |
684 | ||
685 | printf("\nTotal number of Pyevnt() calls %d\n", fTrialsRun); | |
686 | printf("\nMean Q, x1, x2: %f %f %f\n", fQ, fX1, fX2); | |
687 | } | |
688 | ||
689 | void AliGenPythia::AdjustWeights() | |
690 | { | |
691 | // Adjust the weights after generation of all events | |
692 | // | |
693 | if (gAlice) { | |
694 | TParticle *part; | |
695 | Int_t ntrack=gAlice->GetMCApp()->GetNtrack(); | |
696 | for (Int_t i=0; i<ntrack; i++) { | |
697 | part= gAlice->GetMCApp()->Particle(i); | |
698 | part->SetWeight(part->GetWeight()*fKineBias); | |
699 | } | |
700 | } | |
701 | } | |
702 | ||
703 | void AliGenPythia::SetNuclei(Int_t a1, Int_t a2) | |
704 | { | |
705 | // Treat protons as inside nuclei with mass numbers a1 and a2 | |
706 | ||
707 | fAProjectile = a1; | |
708 | fATarget = a2; | |
709 | fSetNuclei = kTRUE; | |
710 | } | |
711 | ||
712 | ||
713 | void AliGenPythia::MakeHeader() | |
714 | { | |
715 | // Builds the event header, to be called after each event | |
716 | if (fHeader) delete fHeader; | |
717 | fHeader = new AliGenPythiaEventHeader("Pythia"); | |
718 | // | |
719 | // Event type | |
720 | ((AliGenPythiaEventHeader*) fHeader)->SetProcessType(fPythia->GetMSTI(1)); | |
721 | // | |
722 | // Number of trials | |
723 | ((AliGenPythiaEventHeader*) fHeader)->SetTrials(fTrials); | |
724 | // | |
725 | // Event Vertex | |
726 | fHeader->SetPrimaryVertex(fVertex); | |
727 | // | |
728 | // Jets that have triggered | |
729 | if (fProcess == kPyJets) | |
730 | { | |
731 | Int_t ntrig, njet; | |
732 | Float_t jets[4][10]; | |
733 | GetJets(njet, ntrig, jets); | |
734 | ||
735 | if (gAlice) { | |
736 | if (gAlice->GetEvNumber()>=fDebugEventFirst && | |
737 | gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(1); | |
738 | } | |
739 | ||
740 | for (Int_t i = 0; i < ntrig; i++) { | |
741 | ((AliGenPythiaEventHeader*) fHeader)->AddJet(jets[0][i], jets[1][i], jets[2][i], | |
742 | jets[3][i]); | |
743 | } | |
744 | } | |
745 | // | |
746 | // Copy relevant information from external header, if present. | |
747 | // | |
748 | Float_t uqJet[4]; | |
749 | ||
750 | if (fRL) { | |
751 | AliGenPythiaEventHeader* exHeader = (AliGenPythiaEventHeader*) (fRL->GetHeader()->GenEventHeader()); | |
752 | for (Int_t i = 0; i < exHeader->NTriggerJets(); i++) | |
753 | { | |
754 | printf("Adding Jet %d %d \n", i, exHeader->NTriggerJets()); | |
755 | ||
756 | ||
757 | exHeader->TriggerJet(i, uqJet); | |
758 | ((AliGenPythiaEventHeader*) fHeader)->AddUQJet(uqJet[0], uqJet[1], uqJet[2], uqJet[3]); | |
759 | } | |
760 | } | |
761 | // | |
762 | // Store quenching parameters | |
763 | // | |
764 | if (fQuench){ | |
765 | Double_t z[4]; | |
766 | Double_t xp, yp; | |
767 | ||
768 | fPythia->GetQuenchingParameters(xp, yp, z); | |
769 | ||
770 | ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp); | |
771 | ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z); | |
772 | } | |
773 | ||
774 | // | |
775 | // Pass header to RunLoader | |
776 | // | |
777 | AliRunLoader::GetRunLoader()->GetHeader()->SetGenEventHeader(fHeader); | |
778 | } | |
779 | ||
780 | ||
781 | Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2) | |
782 | { | |
783 | // Check the kinematic trigger condition | |
784 | // | |
785 | Double_t eta[2]; | |
786 | eta[0] = jet1->Eta(); | |
787 | eta[1] = jet2->Eta(); | |
788 | Double_t phi[2]; | |
789 | phi[0] = jet1->Phi(); | |
790 | phi[1] = jet2->Phi(); | |
791 | Int_t pdg[2]; | |
792 | pdg[0] = jet1->GetPdgCode(); | |
793 | pdg[1] = jet2->GetPdgCode(); | |
794 | Bool_t triggered = kFALSE; | |
795 | ||
796 | if (fProcess == kPyJets) { | |
797 | Int_t njets = 0; | |
798 | Int_t ntrig = 0; | |
799 | Float_t jets[4][10]; | |
800 | // | |
801 | // Use Pythia clustering on parton level to determine jet axis | |
802 | // | |
803 | GetJets(njets, ntrig, jets); | |
804 | ||
805 | if (ntrig) triggered = kTRUE; | |
806 | // | |
807 | } else { | |
808 | Int_t ij = 0; | |
809 | Int_t ig = 1; | |
810 | if (pdg[0] == kGamma) { | |
811 | ij = 1; | |
812 | ig = 0; | |
813 | } | |
814 | //Check eta range first... | |
815 | if ((eta[ij] < fEtaMaxJet && eta[ij] > fEtaMinJet) && | |
816 | (eta[ig] < fEtaMaxGamma && eta[ig] > fEtaMinGamma)) | |
817 | { | |
818 | //Eta is okay, now check phi range | |
819 | if ((phi[ij] < fPhiMaxJet && phi[ij] > fPhiMinJet) && | |
820 | (phi[ig] < fPhiMaxGamma && phi[ig] > fPhiMinGamma)) | |
821 | { | |
822 | triggered = kTRUE; | |
823 | } | |
824 | } | |
825 | } | |
826 | return triggered; | |
827 | } | |
828 | ||
829 | AliGenPythia& AliGenPythia::operator=(const AliGenPythia& rhs) | |
830 | { | |
831 | // Assignment operator | |
832 | rhs.Copy(*this); | |
833 | return *this; | |
834 | } | |
835 | ||
836 | void AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr) | |
837 | { | |
838 | // | |
839 | // Load event into Pythia Common Block | |
840 | // | |
841 | ||
842 | Int_t npart = stack -> GetNprimary(); | |
843 | Int_t n0 = 0; | |
844 | ||
845 | if (!flag) { | |
846 | (fPythia->GetPyjets())->N = npart; | |
847 | } else { | |
848 | n0 = (fPythia->GetPyjets())->N; | |
849 | (fPythia->GetPyjets())->N = n0 + npart; | |
850 | } | |
851 | ||
852 | ||
853 | for (Int_t part = 0; part < npart; part++) { | |
854 | TParticle *MPart = stack->Particle(part); | |
855 | ||
856 | Int_t kf = MPart->GetPdgCode(); | |
857 | Int_t ks = MPart->GetStatusCode(); | |
858 | Int_t idf = MPart->GetFirstDaughter(); | |
859 | Int_t idl = MPart->GetLastDaughter(); | |
860 | ||
861 | if (reHadr) { | |
862 | if (ks == 11 || ks == 12) { | |
863 | ks -= 10; | |
864 | idf = -1; | |
865 | idl = -1; | |
866 | } | |
867 | } | |
868 | ||
869 | Float_t px = MPart->Px(); | |
870 | Float_t py = MPart->Py(); | |
871 | Float_t pz = MPart->Pz(); | |
872 | Float_t e = MPart->Energy(); | |
873 | Float_t m = MPart->GetCalcMass(); | |
874 | ||
875 | ||
876 | (fPythia->GetPyjets())->P[0][part+n0] = px; | |
877 | (fPythia->GetPyjets())->P[1][part+n0] = py; | |
878 | (fPythia->GetPyjets())->P[2][part+n0] = pz; | |
879 | (fPythia->GetPyjets())->P[3][part+n0] = e; | |
880 | (fPythia->GetPyjets())->P[4][part+n0] = m; | |
881 | ||
882 | (fPythia->GetPyjets())->K[1][part+n0] = kf; | |
883 | (fPythia->GetPyjets())->K[0][part+n0] = ks; | |
884 | (fPythia->GetPyjets())->K[3][part+n0] = idf + 1; | |
885 | (fPythia->GetPyjets())->K[4][part+n0] = idl + 1; | |
886 | (fPythia->GetPyjets())->K[2][part+n0] = MPart->GetFirstMother() + 1; | |
887 | } | |
888 | } | |
889 | ||
890 | ||
891 | void AliGenPythia::RecJetsUA1(Int_t& njets, Float_t jets [4][50]) | |
892 | { | |
893 | // | |
894 | // Calls the Pythia jet finding algorithm to find jets in the current event | |
895 | // | |
896 | // | |
897 | // | |
898 | // Save jets | |
899 | Int_t n = fPythia->GetN(); | |
900 | ||
901 | // | |
902 | // Run Jet Finder | |
903 | fPythia->Pycell(njets); | |
904 | Int_t i; | |
905 | for (i = 0; i < njets; i++) { | |
906 | Float_t px = (fPythia->GetPyjets())->P[0][n+i]; | |
907 | Float_t py = (fPythia->GetPyjets())->P[1][n+i]; | |
908 | Float_t pz = (fPythia->GetPyjets())->P[2][n+i]; | |
909 | Float_t e = (fPythia->GetPyjets())->P[3][n+i]; | |
910 | ||
911 | jets[0][i] = px; | |
912 | jets[1][i] = py; | |
913 | jets[2][i] = pz; | |
914 | jets[3][i] = e; | |
915 | } | |
916 | } | |
917 | ||
918 | ||
919 | ||
920 | void AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10]) | |
921 | { | |
922 | // | |
923 | // Calls the Pythia clustering algorithm to find jets in the current event | |
924 | // | |
925 | Int_t n = fPythia->GetN(); | |
926 | nJets = 0; | |
927 | nJetsTrig = 0; | |
928 | if (fJetReconstruction == kCluster) { | |
929 | // | |
930 | // Configure cluster algorithm | |
931 | // | |
932 | fPythia->SetPARU(43, 2.); | |
933 | fPythia->SetMSTU(41, 1); | |
934 | // | |
935 | // Call cluster algorithm | |
936 | // | |
937 | fPythia->Pyclus(nJets); | |
938 | // | |
939 | // Loading jets from common block | |
940 | // | |
941 | } else { | |
942 | ||
943 | // | |
944 | // Run Jet Finder | |
945 | fPythia->Pycell(nJets); | |
946 | } | |
947 | ||
948 | Int_t i; | |
949 | for (i = 0; i < nJets; i++) { | |
950 | Float_t px = (fPythia->GetPyjets())->P[0][n+i]; | |
951 | Float_t py = (fPythia->GetPyjets())->P[1][n+i]; | |
952 | Float_t pz = (fPythia->GetPyjets())->P[2][n+i]; | |
953 | Float_t e = (fPythia->GetPyjets())->P[3][n+i]; | |
954 | Float_t pt = TMath::Sqrt(px * px + py * py); | |
955 | Float_t phi = TMath::Pi() + TMath::ATan2(-py, -px); | |
956 | Float_t theta = TMath::ATan2(pt,pz); | |
957 | Float_t et = e * TMath::Sin(theta); | |
958 | Float_t eta = -TMath::Log(TMath::Tan(theta / 2.)); | |
959 | if ( | |
960 | eta > fEtaMinJet && eta < fEtaMaxJet && | |
961 | phi > fPhiMinJet && phi < fPhiMaxJet && | |
962 | et > fEtMinJet && et < fEtMaxJet | |
963 | ) | |
964 | { | |
965 | jets[0][nJetsTrig] = px; | |
966 | jets[1][nJetsTrig] = py; | |
967 | jets[2][nJetsTrig] = pz; | |
968 | jets[3][nJetsTrig] = e; | |
969 | nJetsTrig++; | |
970 | // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg); | |
971 | } else { | |
972 | // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg); | |
973 | } | |
974 | } | |
975 | } | |
976 | ||
977 | ||
978 | #ifdef never | |
979 | void AliGenPythia::Streamer(TBuffer &R__b) | |
980 | { | |
981 | // Stream an object of class AliGenPythia. | |
982 | ||
983 | if (R__b.IsReading()) { | |
984 | Version_t R__v = R__b.ReadVersion(); if (R__v) { } | |
985 | AliGenerator::Streamer(R__b); | |
986 | R__b >> (Int_t&)fProcess; | |
987 | R__b >> (Int_t&)fStrucFunc; | |
988 | R__b >> (Int_t&)fForceDecay; | |
989 | R__b >> fEnergyCMS; | |
990 | R__b >> fKineBias; | |
991 | R__b >> fTrials; | |
992 | fParentSelect.Streamer(R__b); | |
993 | fChildSelect.Streamer(R__b); | |
994 | R__b >> fXsection; | |
995 | // (AliPythia::Instance())->Streamer(R__b); | |
996 | R__b >> fPtHardMin; | |
997 | R__b >> fPtHardMax; | |
998 | // if (fDecayer) fDecayer->Streamer(R__b); | |
999 | } else { | |
1000 | R__b.WriteVersion(AliGenPythia::IsA()); | |
1001 | AliGenerator::Streamer(R__b); | |
1002 | R__b << (Int_t)fProcess; | |
1003 | R__b << (Int_t)fStrucFunc; | |
1004 | R__b << (Int_t)fForceDecay; | |
1005 | R__b << fEnergyCMS; | |
1006 | R__b << fKineBias; | |
1007 | R__b << fTrials; | |
1008 | fParentSelect.Streamer(R__b); | |
1009 | fChildSelect.Streamer(R__b); | |
1010 | R__b << fXsection; | |
1011 | // R__b << fPythia; | |
1012 | R__b << fPtHardMin; | |
1013 | R__b << fPtHardMax; | |
1014 | // fDecayer->Streamer(R__b); | |
1015 | } | |
1016 | } | |
1017 | #endif | |
1018 | ||
1019 | ||
1020 |