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 // Implementation of AliDecayer using Pythia8
19 // Author: andreas.morsch@cern.ch
22 #include <TLorentzVector.h>
23 #include "AliTPythia8.h"
24 #include "AliDecayerPythia8.h"
25 #include "ParticleData.h"
27 ClassImp(AliDecayerPythia8)
29 Bool_t AliDecayerPythia8::fgInit = kFALSE;
31 AliDecayerPythia8::AliDecayerPythia8():
33 fPythia8(AliTPythia8::Instance()),
39 fPythia8->Pythia8()->readString("SoftQCD:elastic = on");
40 fPythia8->Pythia8()->init();
43 //___________________________________________________________________________
44 void AliDecayerPythia8::Decay(Int_t pdg, TLorentzVector* p)
46 // Decay a single particle
48 AppendParticle(pdg, p);
49 Int_t idPart = fPythia8->Pythia8()->event[0].id();
50 fPythia8->Pythia8()->particleData.mayDecay(idPart,kTRUE);
51 fPythia8->Pythia8()->moreDecays();
52 if (fDebug > 0) fPythia8->EventListing();
55 //___________________________________________________________________________
56 Int_t AliDecayerPythia8::ImportParticles(TClonesArray *particles)
58 //import the decay products into particles array
59 return (fPythia8->ImportParticles(particles, "All"));
63 void AliDecayerPythia8::Init()
69 // fPythia->SetDecayTable();
72 // Switch on heavy flavor decays
75 Int_t heavy[14] = {411, 421, 431, 4122, 4132, 4232, 4332, 511, 521, 531, 5122, 5132, 5232, 5332};
76 // fPythia->ResetDecayTable();
77 for (j=0; j < 14; j++) {
78 if (fDecay == kNoDecayHeavy) {
79 fPythia8->ReadString(Form("%d:onMode = off", heavy[j]));
81 fPythia8->ReadString(Form("%d:onMode = on", heavy[j]));
86 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
88 if (fDecay != kNeutralPion) {
89 fPythia8->ReadString("111:onMode = off");
91 fPythia8->ReadString("111:onMode = on");
94 fPythia8->ReadString("310:onMode = off");
95 fPythia8->ReadString("3112:onMode = off");
96 fPythia8->ReadString("3212:onMode = off");
97 fPythia8->ReadString("3222:onMode = off");
98 fPythia8->ReadString("3312:onMode = off");
99 fPythia8->ReadString("3322:onMode = off");
100 fPythia8->ReadString("3334:onMode = off");
101 // .. Force decay channels
105 void AliDecayerPythia8::ForceDecay()
108 // Force a particle decay mode
109 // Switch heavy flavour production off if requested
110 if (!fHeavyFlavour) SwitchOffHeavyFlavour();
112 Decay_t decay = fDecay;
113 fPythia8->ReadString("HadronLevel:Decay = on");
115 if (decay == kNoDecayHeavy) return;
123 fPythia8->ReadString("511:onMode = off");
124 fPythia8->ReadString("511:onIfAny = 13 443 100443");
126 fPythia8->ReadString("521:onMode = off");
127 fPythia8->ReadString("521:onIfAny = 13 443 100443");
129 fPythia8->ReadString("531:onMode = off");
130 fPythia8->ReadString("531:onIfAny = 13 443 100443");
132 fPythia8->ReadString("5122:onMode = off");
133 fPythia8->ReadString("5122:onIfAny = 13 443 100443");
135 fPythia8->ReadString("5132:onMode = off");
136 fPythia8->ReadString("5132:onIfAny = 13 443 100443");
138 fPythia8->ReadString("5232:onMode = off");
139 fPythia8->ReadString("5232:onIfAny = 13 443 100443");
141 fPythia8->ReadString("5332:onMode = off");
142 fPythia8->ReadString("5332:onIfAny = 13 443 100443");
144 fPythia8->ReadString("100443:onMode = off");
145 fPythia8->ReadString("100443:onIfAny = 443");
147 fPythia8->ReadString("443:onMode = off");
148 fPythia8->ReadString("443:onIfAll = 13 13");
150 fPythia8->ReadString("411:onMode = off");
151 fPythia8->ReadString("411:onIfAll = 13");
153 fPythia8->ReadString("421:onMode = off");
154 fPythia8->ReadString("421:onIfAll = 13");
156 fPythia8->ReadString("431:onMode = off");
157 fPythia8->ReadString("431:onIfAll = 13");
159 fPythia8->ReadString("4122:onMode = off");
160 fPythia8->ReadString("4122:onIfAll = 13");
162 fPythia8->ReadString("4132:onMode = off");
163 fPythia8->ReadString("4132:onIfAll = 13");
165 fPythia8->ReadString("4232:onMode = off");
166 fPythia8->ReadString("4232:onIfAll = 13");
168 fPythia8->ReadString("4332:onMode = off");
169 fPythia8->ReadString("4332:onIfAll = 13");
172 case kChiToJpsiGammaToMuonMuon:
173 // Chi_1c -> J/Psi Gamma
174 fPythia8->ReadString("20443:onMode = off");
175 fPythia8->ReadString("20443:onIfAll = 443 22");
176 // Chi_2c -> J/Psi Gamma
177 fPythia8->ReadString("445:onMode = off");
178 fPythia8->ReadString("445:onIfAll = 443 22");
180 fPythia8->ReadString("443:onMode = off");
181 fPythia8->ReadString("443:onIfAll = 13 13");
183 case kChiToJpsiGammaToElectronElectron:
184 // Chi_1c -> J/Psi Gamma
185 fPythia8->ReadString("20443:onMode = off");
186 fPythia8->ReadString("20443:onIfAll = 443 22");
187 // Chi_2c -> J/Psi Gamma
188 fPythia8->ReadString("445:onMode = off");
189 fPythia8->ReadString("445:onIfAll = 443 22");
191 fPythia8->ReadString("443:onMode = off");
192 fPythia8->ReadString("443:onIfAll = 11 11");
197 fPythia8->ReadString("511:onMode = off");
198 fPythia8->ReadString("511:onIfAny = 13");
200 fPythia8->ReadString("521:onMode = off");
201 fPythia8->ReadString("521:onIfAny = 13");
203 fPythia8->ReadString("531:onMode = off");
204 fPythia8->ReadString("531:onIfAny = 13");
206 fPythia8->ReadString("5122:onMode = off");
207 fPythia8->ReadString("5122:onIfAny = 13");
209 fPythia8->ReadString("5132:onMode = off");
210 fPythia8->ReadString("5132:onIfAny = 13");
212 fPythia8->ReadString("5232:onMode = off");
213 fPythia8->ReadString("5232:onIfAny = 13");
215 fPythia8->ReadString("5332:onMode = off");
216 fPythia8->ReadString("5332:onIfAny = 13");
220 fPythia8->ReadString("411:onMode = off");
221 fPythia8->ReadString("411:onIfAll = 13");
223 fPythia8->ReadString("421:onMode = off");
224 fPythia8->ReadString("421:onIfAll = 13");
226 fPythia8->ReadString("431:onMode = off");
227 fPythia8->ReadString("431:onIfAll = 13");
229 fPythia8->ReadString("4122:onMode = off");
230 fPythia8->ReadString("4122:onIfAll = 13");
232 fPythia8->ReadString("4132:onMode = off");
233 fPythia8->ReadString("4132:onIfAll = 13");
235 fPythia8->ReadString("4232:onMode = off");
236 fPythia8->ReadString("4232:onIfAll = 13");
238 fPythia8->ReadString("4332:onMode = off");
239 fPythia8->ReadString("4332:onIfAll = 13");
241 fPythia8->ReadString("511:onMode = off");
242 fPythia8->ReadString("511:onIfAny = 13");
244 fPythia8->ReadString("521:onMode = off");
245 fPythia8->ReadString("521:onIfAny = 13");
247 fPythia8->ReadString("531:onMode = off");
248 fPythia8->ReadString("531:onIfAny = 13");
250 fPythia8->ReadString("5122:onMode = off");
251 fPythia8->ReadString("5122:onIfAny = 13");
253 fPythia8->ReadString("5132:onMode = off");
254 fPythia8->ReadString("5132:onIfAny = 13");
256 fPythia8->ReadString("5232:onMode = off");
257 fPythia8->ReadString("5232:onIfAny = 13");
259 fPythia8->ReadString("5332:onMode = off");
260 fPythia8->ReadString("5332:onIfAny = 13");
265 fPythia8->ReadString("443:onMode = off");
266 fPythia8->ReadString("443:onIfAll = 13 13");
270 fPythia8->ReadString("113:onMode = off");
271 fPythia8->ReadString("113:onIfAll = 13 13");
273 fPythia8->ReadString("221:onMode = off");
274 fPythia8->ReadString("221:onIfAll = 13 13");
276 fPythia8->ReadString("223:onMode = off");
277 fPythia8->ReadString("223:onIfAll = 13 13");
279 fPythia8->ReadString("333:onMode = off");
280 fPythia8->ReadString("333:onIfAll = 13 13");
282 fPythia8->ReadString("443:onMode = off");
283 fPythia8->ReadString("443:onIfAll = 13 13");
285 fPythia8->ReadString("100443:onMode = off");
286 fPythia8->ReadString("100443:onIfAll = 13 13");
288 fPythia8->ReadString("553:onMode = off");
289 fPythia8->ReadString("553:onIfAll = 13 13");
291 fPythia8->ReadString("100553:onMode = off");
292 fPythia8->ReadString("100553:onIfAll = 13 13");
294 fPythia8->ReadString("200553:onMode = off");
295 fPythia8->ReadString("200553:onIfAll = 13 13");
297 case kBSemiElectronic:
299 fPythia8->ReadString("511:onMode = off");
300 fPythia8->ReadString("511:onIfAny = 11");
302 fPythia8->ReadString("521:onMode = off");
303 fPythia8->ReadString("521:onIfAny = 11");
305 fPythia8->ReadString("531:onMode = off");
306 fPythia8->ReadString("531:onIfAny = 11");
308 fPythia8->ReadString("5122:onMode = off");
309 fPythia8->ReadString("5122:onIfAny = 11");
311 fPythia8->ReadString("5132:onMode = off");
312 fPythia8->ReadString("5132:onIfAny = 11");
314 fPythia8->ReadString("5232:onMode = off");
315 fPythia8->ReadString("5232:onIfAny = 11");
317 fPythia8->ReadString("5332:onMode = off");
318 fPythia8->ReadString("5332:onIfAny = 11");
320 case kSemiElectronic:
322 fPythia8->ReadString("411:onMode = off");
323 fPythia8->ReadString("411:onIfAll = 11");
325 fPythia8->ReadString("421:onMode = off");
326 fPythia8->ReadString("421:onIfAll = 11");
328 fPythia8->ReadString("431:onMode = off");
329 fPythia8->ReadString("431:onIfAll = 11");
331 fPythia8->ReadString("4122:onMode = off");
332 fPythia8->ReadString("4122:onIfAll = 11");
334 fPythia8->ReadString("4132:onMode = off");
335 fPythia8->ReadString("4132:onIfAll = 11");
337 fPythia8->ReadString("4232:onMode = off");
338 fPythia8->ReadString("4232:onIfAll = 11");
340 fPythia8->ReadString("4332:onMode = off");
341 fPythia8->ReadString("4332:onIfAll = 11");
343 fPythia8->ReadString("511:onMode = off");
344 fPythia8->ReadString("511:onIfAny = 11");
346 fPythia8->ReadString("521:onMode = off");
347 fPythia8->ReadString("521:onIfAny = 11");
349 fPythia8->ReadString("531:onMode = off");
350 fPythia8->ReadString("531:onIfAny = 11");
352 fPythia8->ReadString("5122:onMode = off");
353 fPythia8->ReadString("5122:onIfAny = 11");
355 fPythia8->ReadString("5132:onMode = off");
356 fPythia8->ReadString("5132:onIfAny = 11");
358 fPythia8->ReadString("5232:onMode = off");
359 fPythia8->ReadString("5232:onIfAny = 11");
361 fPythia8->ReadString("5332:onMode = off");
362 fPythia8->ReadString("5332:onIfAny = 11");
366 fPythia8->ReadString("113:onMode = off");
367 fPythia8->ReadString("113:onIfAll = 11 11");
369 fPythia8->ReadString("221:onMode = off");
370 fPythia8->ReadString("221:onIfAll = 11 11");
372 fPythia8->ReadString("223:onMode = off");
373 fPythia8->ReadString("223:onIfAll = 11 11");
375 fPythia8->ReadString("333:onMode = off");
376 fPythia8->ReadString("333:onIfAll = 11 11");
378 fPythia8->ReadString("443:onMode = off");
379 fPythia8->ReadString("443:onIfAll = 11 11");
381 fPythia8->ReadString("100443:onMode = off");
382 fPythia8->ReadString("100443:onIfAll = 11 11");
384 fPythia8->ReadString("553:onMode = off");
385 fPythia8->ReadString("553:onIfAll = 11 11");
387 fPythia8->ReadString("100553:onMode = off");
388 fPythia8->ReadString("100553:onIfAll = 11 11");
390 fPythia8->ReadString("200553:onMode = off");
391 fPythia8->ReadString("200553:onIfAll = 11 11");
394 // B0 -> J/Psi (Psi') X
395 fPythia8->ReadString("511:onMode = off");
396 fPythia8->ReadString("511:onIfAny = 443 100443");
397 // B+/- -> J/Psi (Psi') X
398 fPythia8->ReadString("521:onMode = off");
399 fPythia8->ReadString("521:onIfAny = 443 100443");
400 // B_s -> J/Psi (Psi') X
401 fPythia8->ReadString("531:onMode = off");
402 fPythia8->ReadString("531:onIfAny = 443 100443");
403 // Lambda_b -> J/Psi (Psi') X
404 fPythia8->ReadString("5122:onMode = off");
405 fPythia8->ReadString("5122:onIfAny = 443 100443");
408 fPythia8->ReadString("443:onMode = off");
409 fPythia8->ReadString("443:onIfAll = 13 13");
411 fPythia8->ReadString("100443:onMode = off");
412 fPythia8->ReadString("100443:onIfAll = 13 13");
414 case kBPsiPrimeDiMuon:
416 fPythia8->ReadString("511:onMode = off");
417 fPythia8->ReadString("511:onIfAny = 100443");
419 fPythia8->ReadString("521:onMode = off");
420 fPythia8->ReadString("521:onIfAny = 100443");
422 fPythia8->ReadString("531:onMode = off");
423 fPythia8->ReadString("531:onIfAny = 100443");
424 // Lambda_b -> Psi' X
425 fPythia8->ReadString("5122:onMode = off");
426 fPythia8->ReadString("5122:onIfAny = 100443");
429 fPythia8->ReadString("100443:onMode = off");
430 fPythia8->ReadString("100443:onIfAll = 13 13");
432 case kBJpsiDiElectron:
434 fPythia8->ReadString("511:onMode = off");
435 fPythia8->ReadString("511:onIfAny = 443");
437 fPythia8->ReadString("521:onMode = off");
438 fPythia8->ReadString("521:onIfAny = 443");
440 fPythia8->ReadString("531:onMode = off");
441 fPythia8->ReadString("531:onIfAny = 443");
443 fPythia8->ReadString("5122:onMode = off");
444 fPythia8->ReadString("5122:onIfAny = 443");
447 fPythia8->ReadString("443:onMode = off");
448 fPythia8->ReadString("443:onIfAll = 11 11");
453 fPythia8->ReadString("511:onMode = off");
454 fPythia8->ReadString("511:onIfAny = 443");
456 fPythia8->ReadString("521:onMode = off");
457 fPythia8->ReadString("521:onIfAny = 443");
459 fPythia8->ReadString("531:onMode = off");
460 fPythia8->ReadString("531:onIfAny = 443");
462 fPythia8->ReadString("5122:onMode = off");
463 fPythia8->ReadString("5122:onIfAny = 443");
465 case kBPsiPrimeDiElectron:
467 fPythia8->ReadString("511:onMode = off");
468 fPythia8->ReadString("511:onIfAny = 100443");
470 fPythia8->ReadString("521:onMode = off");
471 fPythia8->ReadString("521:onIfAny = 100443");
473 fPythia8->ReadString("531:onMode = off");
474 fPythia8->ReadString("531:onIfAny = 100443");
475 // Lambda_b -> Psi' X
476 fPythia8->ReadString("5122:onMode = off");
477 fPythia8->ReadString("5122:onIfAny = 100443");
480 fPythia8->ReadString("100443:onMode = off");
481 fPythia8->ReadString("100443:onIfAll = 11 11");
485 fPythia8->ReadString("211:onMode = off");
486 fPythia8->ReadString("211:onIfAny = 13");
490 fPythia8->ReadString("321:onMode = off");
491 fPythia8->ReadString("321:onIfAny = 13");
495 fPythia8->ReadString("211:onMode = off");
496 fPythia8->ReadString("211:onIfAny = 13");
497 fPythia8->ReadString("321:onMode = off");
498 fPythia8->ReadString("321:onIfAny = 13");
502 fPythia8->ReadString("24:onMode = off");
503 fPythia8->ReadString("24:onIfAny = 13");
507 fPythia8->ReadString("24:onMode = off");
508 fPythia8->ReadString("24:onIfAny = 4");
510 case kWToCharmToMuon:
512 fPythia8->ReadString("24:onMode = off");
513 fPythia8->ReadString("24:onIfAny = 4");
515 fPythia8->ReadString("411:onMode = off");
516 fPythia8->ReadString("411:onIfAll = 13");
518 fPythia8->ReadString("421:onMode = off");
519 fPythia8->ReadString("421:onIfAll = 13");
521 fPythia8->ReadString("431:onMode = off");
522 fPythia8->ReadString("431:onIfAll = 13");
524 fPythia8->ReadString("4122:onMode = off");
525 fPythia8->ReadString("4122:onIfAll = 13");
527 fPythia8->ReadString("4132:onMode = off");
528 fPythia8->ReadString("4132:onIfAll = 13");
530 fPythia8->ReadString("4232:onMode = off");
531 fPythia8->ReadString("4232:onIfAll = 13");
533 fPythia8->ReadString("4332:onMode = off");
534 fPythia8->ReadString("4332:onIfAll = 13");
538 fPythia8->ReadString("23:onMode = off");
539 fPythia8->ReadString("23:onIfAll = 13 13");
543 fPythia8->ReadString("23:onMode = off");
544 fPythia8->ReadString("23:onIfAll = 11 11");
549 case kHadronicDWithout4Bodies:
554 fPythia8->ReadString("333:onMode = off");
555 fPythia8->ReadString("333:onIfAll = 321 321");
559 fPythia8->ReadString("3334:onMode = off");
560 fPythia8->ReadString("3334:onIfAll = 3122 321 ");
564 fPythia8->ReadString("3122:onMode = off");
565 fPythia8->ReadString("3122:onIfAll = 2212 211 ");
568 fPythia8->ReadString("5122:onMode = off");
569 fPythia8->ReadString("4122:onMode = off");
570 fPythia8->ReadString("5122:onIfAll = 4122");
571 fPythia8->ReadString("4122:onIfAll = 3122");
576 fPythia8->ReadString("HadronLevel:Decay = off");
581 case kPsiPrimeJpsiDiElectron:
589 Float_t AliDecayerPythia8::GetPartialBranchingRatio(Int_t ipart)
591 // Get the partial branching ration for the forced decay channels
593 Pythia8::Pythia* thePythia = fPythia8->Pythia8();
594 Pythia8::ParticleData & table = thePythia->particleData;
595 Pythia8::ParticleDataEntry* pd = table.particleDataEntryPtr(ipart);
597 Int_t nc = pd->sizeChannels();
600 // Loop over decay channels
601 for (Int_t ic = 0; ic < nc; ic++) {
602 Pythia8::DecayChannel& decCh = pd->channel(ic);
603 for (Int_t i = 0; i < decCh.multiplicity(); i++) {
604 br += decCh.bRatio();
611 Float_t AliDecayerPythia8::GetLifetime(Int_t kf)
613 // Return lifetime of particle
614 Pythia8::Pythia* thePythia = fPythia8->Pythia8();
615 Pythia8::ParticleData& table = thePythia->particleData;
616 Float_t tau = table.tau0(kf);
620 void AliDecayerPythia8::SwitchOffHeavyFlavour()
622 // Switch off heavy flavour production
624 // Maximum number of quark flavours used in pdf
625 fPythia8->ReadString("PDFinProcess:nQuarkIn = 3");
626 // Maximum number of flavors that can be used in showers
627 fPythia8->ReadString("SpaceShower:nQuarkIn = 3");
628 fPythia8->ReadString("TimeShower:nGammaToQuark = 3");
629 fPythia8->ReadString("TimeShower:nGluonToQuark = 3");
633 void AliDecayerPythia8::ForceHadronicD(Int_t optUse4Bodies)
636 // Force golden D decay modes
639 fPythia8->ReadString("313:onMode = off");
640 fPythia8->ReadString("313:onIfAll = 321 211");
642 fPythia8->ReadString("333:onMode = off");
643 fPythia8->ReadString("333:onIfAll = 321 321");
644 // for D0 -> rho0 pi+ k-
645 fPythia8->ReadString("113:onMode = off");
646 fPythia8->ReadString("113:onIfAll = 211 211");
647 // for Lambda_c -> Delta++ K-
648 fPythia8->ReadString("2224:onMode = off");
649 fPythia8->ReadString("2224:onIfAll = 2212 211");
650 // for Lambda_c -> Lambda(1520) K-
651 fPythia8->ReadString("3124:onMode = off");
652 fPythia8->ReadString("3124:onIfAll = 2212 321");
655 fPythia8->ReadString("411:onMode = off");
656 fPythia8->ReadString("421:onMode = off");
657 fPythia8->ReadString("431:onMode = off");
658 fPythia8->ReadString("4112:onMode = off");
659 fPythia8->ReadString("4122:onMode = off");
662 fPythia8->ReadString("411:onIfMatch = 321 211 211");
664 fPythia8->ReadString("411:onIfMatch = 313 211");
666 fPythia8->ReadString("421:onIfMatch = 321 211");
670 fPythia8->ReadString("421:onIfMatch = 321 211 211 211");
672 fPythia8->ReadString("421:onIfMatch = 321 211 113");
674 fPythia8->ReadString("421:onIfMatch = 313 211 211");
678 fPythia8->ReadString("431:onIfMatch = 321 313");
680 fPythia8->ReadString("431:onIfMatch = 333 211");
683 fPythia8->ReadString("4122:onIfMatch = 2212 313");
684 // Lambda_c -> Delta K
685 fPythia8->ReadString("4122:onIfMatch = 2224 321");
686 // Lambda_c -> Lambda(1520) pi
687 fPythia8->ReadString("4122:onIfMatch = 3124 211");
688 // Lambda_c -> p K pi
689 fPythia8->ReadString("4122:onIfMatch = 2212 321 211");
690 // Lambda_c -> Lambda pi
691 fPythia8->ReadString("4122:onIfMatch = 3122 211");
695 //___________________________________________________________________________
696 void AliDecayerPythia8::ReadDecayTable()
698 //to read a decay table (not yet implemented)
702 //___________________________________________________________________________
703 void AliDecayerPythia8::AppendParticle(Int_t pdg, TLorentzVector* p)
705 // Append a particle to the stack
706 fPythia8->Pythia8()->event.append(pdg, 11, 0, 0, p->Px(), p->Py(), p->Pz(), p->E(), p->M());
710 //___________________________________________________________________________
711 void AliDecayerPythia8::ClearEvent()
713 // Clear the event stack
714 fPythia8->Pythia8()->event.clear();