]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PYTHIA8/AliDecayerPythia8.cxx
Small update for the status code (Raphaelle)
[u/mrichter/AliRoot.git] / PYTHIA8 / AliDecayerPythia8.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // Implementation of AliDecayer using Pythia8
19 // Author: andreas.morsch@cern.ch
20 #include <TMath.h>
21 #include <TPDGCode.h>
22 #include <TLorentzVector.h>
23 #include "AliTPythia8.h"
24 #include "AliDecayerPythia8.h"
25 #include "ParticleData.h"
26
27 ClassImp(AliDecayerPythia8)
28
29 Bool_t AliDecayerPythia8::fgInit = kFALSE;
30
31 AliDecayerPythia8::AliDecayerPythia8():
32   TVirtualMCDecayer(),
33   fPythia8(new AliTPythia8()),
34   fDebug(0),
35   fDecay(kAll),
36   fHeavyFlavour(kTRUE)
37 {
38     // Constructor
39    fPythia8->Pythia8()->readString("SoftQCD:elastic = on");
40    fPythia8->Pythia8()->init();
41 }
42
43 //___________________________________________________________________________
44 void AliDecayerPythia8::Decay(Int_t pdg, TLorentzVector* p)
45 {
46    // Decay a single particle
47    ClearEvent();
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();
53 }
54
55 //___________________________________________________________________________
56 Int_t AliDecayerPythia8::ImportParticles(TClonesArray *particles)
57 {
58    //import the decay products into particles array
59    return (fPythia8->ImportParticles(particles, "All"));
60 }
61
62
63 void AliDecayerPythia8::Init()
64 {
65 // Initialisation
66 //
67     if (!fgInit) {
68         fgInit = kTRUE;
69         // fPythia->SetDecayTable();
70     }
71
72 // Switch on heavy flavor decays
73     
74     Int_t j;
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             AliTPythia8::Instance()->ReadString(Form("%d:onMode = off", heavy[j]));
80         } else {
81             AliTPythia8::Instance()->ReadString(Form("%d:onMode = on", heavy[j]));
82         }
83     }
84     
85
86 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
87     
88     if (fDecay != kNeutralPion) {
89         AliTPythia8::Instance()->ReadString("111:onMode = off");
90     } else {
91         AliTPythia8::Instance()->ReadString("111:onMode = on");
92     }
93
94     AliTPythia8::Instance()->ReadString("310:onMode = off");
95     AliTPythia8::Instance()->ReadString("3112:onMode = off");
96     AliTPythia8::Instance()->ReadString("3212:onMode = off");
97     AliTPythia8::Instance()->ReadString("3222:onMode = off");
98     AliTPythia8::Instance()->ReadString("3312:onMode = off");
99     AliTPythia8::Instance()->ReadString("3322:onMode = off");
100     AliTPythia8::Instance()->ReadString("3334:onMode = off");
101 // .. Force decay channels
102     ForceDecay();
103 }
104
105 void AliDecayerPythia8::ForceDecay()
106 {
107 // 
108 // Force a particle decay mode
109 // Switch heavy flavour production off if requested
110     if (!fHeavyFlavour) SwitchOffHeavyFlavour();
111 //
112     Decay_t decay = fDecay;
113     AliTPythia8::Instance()->ReadString("HadronLevel:Decay = on");
114     
115     if (decay == kNoDecayHeavy) return;
116
117 //
118 // select mode    
119     switch (decay) 
120     {
121     case kHardMuons:
122 //      B0 -> mu X 
123         AliTPythia8::Instance()->ReadString("511:onMode = off");
124         AliTPythia8::Instance()->ReadString("511:onIfAny = 13 443 100443");
125 //      B+/- -> mu X 
126         AliTPythia8::Instance()->ReadString("521:onMode = off");
127         AliTPythia8::Instance()->ReadString("521:onIfAny = 13 443 100443");
128 //      Bs -> mu X 
129         AliTPythia8::Instance()->ReadString("531:onMode = off");
130         AliTPythia8::Instance()->ReadString("531:onIfAny = 13 443 100443");
131 //      Lambda_b -> mu X 
132         AliTPythia8::Instance()->ReadString("5122:onMode = off");
133         AliTPythia8::Instance()->ReadString("5122:onIfAny = 13 443 100443");
134 //      Sigma_b- -> mu X 
135         AliTPythia8::Instance()->ReadString("5132:onMode = off");
136         AliTPythia8::Instance()->ReadString("5132:onIfAny = 13 443 100443");
137 //      Sigma_b0 -> mu X
138         AliTPythia8::Instance()->ReadString("5232:onMode = off");
139         AliTPythia8::Instance()->ReadString("5232:onIfAny = 13 443 100443");
140 //      Omega_b  -> mu X
141         AliTPythia8::Instance()->ReadString("5332:onMode = off");
142         AliTPythia8::Instance()->ReadString("5332:onIfAny = 13 443 100443");
143 //      Psi' -> mu+ mu-
144         AliTPythia8::Instance()->ReadString("100443:onMode = off");
145         AliTPythia8::Instance()->ReadString("100443:onIfAny = 443");
146 //      Psi  -> mu+ mu-
147         AliTPythia8::Instance()->ReadString("443:onMode = off");
148         AliTPythia8::Instance()->ReadString("443:onIfAll = 13 13");
149 //      D+/- -> mu X
150         AliTPythia8::Instance()->ReadString("411:onMode = off");
151         AliTPythia8::Instance()->ReadString("411:onIfAll = 13");
152 //      D0   -> mu X
153         AliTPythia8::Instance()->ReadString("421:onMode = off");
154         AliTPythia8::Instance()->ReadString("421:onIfAll = 13");
155 //      D_s  -> mu X
156         AliTPythia8::Instance()->ReadString("431:onMode = off");
157         AliTPythia8::Instance()->ReadString("431:onIfAll = 13");
158 //      Lambda_c -> mu X
159         AliTPythia8::Instance()->ReadString("4122:onMode = off");
160         AliTPythia8::Instance()->ReadString("4122:onIfAll = 13");
161 //      Sigma_c  -> mu X
162         AliTPythia8::Instance()->ReadString("4132:onMode = off");
163         AliTPythia8::Instance()->ReadString("4132:onIfAll = 13");
164 //      Sigma_c+ -> mu X
165         AliTPythia8::Instance()->ReadString("4232:onMode = off");
166         AliTPythia8::Instance()->ReadString("4232:onIfAll = 13");
167 //      Omega_c  -> mu X
168         AliTPythia8::Instance()->ReadString("4332:onMode = off");
169         AliTPythia8::Instance()->ReadString("4332:onIfAll = 13");
170
171         break;
172    case kChiToJpsiGammaToMuonMuon:
173 // Chi_1c  -> J/Psi  Gamma
174         AliTPythia8::Instance()->ReadString("20443:onMode = off");
175         AliTPythia8::Instance()->ReadString("20443:onIfAll = 443 22");
176 // Chi_2c  -> J/Psi  Gamma
177         AliTPythia8::Instance()->ReadString("445:onMode = off");
178         AliTPythia8::Instance()->ReadString("445:onIfAll = 443 22");
179 // J/Psi -> mu+ mu-
180         AliTPythia8::Instance()->ReadString("443:onMode = off");
181         AliTPythia8::Instance()->ReadString("443:onIfAll = 13 13");
182         break;
183     case kChiToJpsiGammaToElectronElectron:
184 // Chi_1c  -> J/Psi  Gamma
185         AliTPythia8::Instance()->ReadString("20443:onMode = off");
186         AliTPythia8::Instance()->ReadString("20443:onIfAll = 443 22");
187 // Chi_2c  -> J/Psi  Gamma
188         AliTPythia8::Instance()->ReadString("445:onMode = off");
189         AliTPythia8::Instance()->ReadString("445:onIfAll = 443 22");
190 // J/Psi -> e+ e-
191         AliTPythia8::Instance()->ReadString("443:onMode = off");
192         AliTPythia8::Instance()->ReadString("443:onIfAll = 11 11");
193         break;
194
195     case kBSemiMuonic:
196 //      B0 -> mu X 
197         AliTPythia8::Instance()->ReadString("511:onMode = off");
198         AliTPythia8::Instance()->ReadString("511:onIfAny = 13");
199 //      B+/- -> mu X 
200         AliTPythia8::Instance()->ReadString("521:onMode = off");
201         AliTPythia8::Instance()->ReadString("521:onIfAny = 13");
202 //      B_s -> mu X 
203         AliTPythia8::Instance()->ReadString("531:onMode = off");
204         AliTPythia8::Instance()->ReadString("531:onIfAny = 13");
205 //      Lambda_b -> mu X 
206         AliTPythia8::Instance()->ReadString("5122:onMode = off");
207         AliTPythia8::Instance()->ReadString("5122:onIfAny = 13");
208 //      Sigma_b -> mu X 
209         AliTPythia8::Instance()->ReadString("5132:onMode = off");
210         AliTPythia8::Instance()->ReadString("5132:onIfAny = 13");
211 //      Sigma_b0 -> mu X 
212         AliTPythia8::Instance()->ReadString("5232:onMode = off");
213         AliTPythia8::Instance()->ReadString("5232:onIfAny = 13");
214 //      Omega_b  -> mu X 
215         AliTPythia8::Instance()->ReadString("5332:onMode = off");
216         AliTPythia8::Instance()->ReadString("5332:onIfAny = 13");
217         break;
218     case kSemiMuonic:
219 //      D+- -> mu X
220         AliTPythia8::Instance()->ReadString("411:onMode = off");
221         AliTPythia8::Instance()->ReadString("411:onIfAll = 13");
222 //      D0  -> mu X
223         AliTPythia8::Instance()->ReadString("421:onMode = off");
224         AliTPythia8::Instance()->ReadString("421:onIfAll = 13");
225 //      D_s  -> mu X
226         AliTPythia8::Instance()->ReadString("431:onMode = off");
227         AliTPythia8::Instance()->ReadString("431:onIfAll = 13");
228 //      Lambda_c -> mu X
229         AliTPythia8::Instance()->ReadString("4122:onMode = off");
230         AliTPythia8::Instance()->ReadString("4122:onIfAll = 13");
231 //      Sigma_c  -> mu X
232         AliTPythia8::Instance()->ReadString("4132:onMode = off");
233         AliTPythia8::Instance()->ReadString("4132:onIfAll = 13");
234 //      Sigma  -> mu X
235         AliTPythia8::Instance()->ReadString("4232:onMode = off");
236         AliTPythia8::Instance()->ReadString("4232:onIfAll = 13");
237 //      Omega_c  -> mu X
238         AliTPythia8::Instance()->ReadString("4332:onMode = off");
239         AliTPythia8::Instance()->ReadString("4332:onIfAll = 13");
240 //      B0       -> mu X
241         AliTPythia8::Instance()->ReadString("511:onMode = off");
242         AliTPythia8::Instance()->ReadString("511:onIfAny = 13");
243 //      B+/-     -> mu X
244         AliTPythia8::Instance()->ReadString("521:onMode = off");
245         AliTPythia8::Instance()->ReadString("521:onIfAny = 13");
246 //      B_s      -> mu X
247         AliTPythia8::Instance()->ReadString("531:onMode = off");
248         AliTPythia8::Instance()->ReadString("531:onIfAny = 13");
249 //      Lambda_c -> mu X
250         AliTPythia8::Instance()->ReadString("5122:onMode = off");
251         AliTPythia8::Instance()->ReadString("5122:onIfAny = 13");
252 //      Sigma_c  -> mu X
253         AliTPythia8::Instance()->ReadString("5132:onMode = off");
254         AliTPythia8::Instance()->ReadString("5132:onIfAny = 13");
255 //      Sigma_c  -> mu X
256         AliTPythia8::Instance()->ReadString("5232:onMode = off");
257         AliTPythia8::Instance()->ReadString("5232:onIfAny = 13");
258 //      Omega_c  -> mu X
259         AliTPythia8::Instance()->ReadString("5332:onMode = off");
260         AliTPythia8::Instance()->ReadString("5332:onIfAny = 13");
261
262         break;
263     case kJpsiDiMuon:
264 //      J/Psi-> mu+ mu-
265         AliTPythia8::Instance()->ReadString("443:onMode = off");
266         AliTPythia8::Instance()->ReadString("443:onIfAll = 13 13");
267         break;
268     case kDiMuon:
269 //      Rho -> mu+ mu-
270         AliTPythia8::Instance()->ReadString("113:onMode = off");
271         AliTPythia8::Instance()->ReadString("113:onIfAll = 13 13");
272 //      Eta-> mu+ mu-
273         AliTPythia8::Instance()->ReadString("221:onMode = off");
274         AliTPythia8::Instance()->ReadString("221:onIfAll = 13 13");
275 //      omega-> mu+ mu-
276         AliTPythia8::Instance()->ReadString("223:onMode = off");
277         AliTPythia8::Instance()->ReadString("223:onIfAll = 13 13");
278 //      phi-> mu+ mu-
279         AliTPythia8::Instance()->ReadString("333:onMode = off");
280         AliTPythia8::Instance()->ReadString("333:onIfAll = 13 13");
281 //      J/Psi-> mu+ mu-
282         AliTPythia8::Instance()->ReadString("443:onMode = off");
283         AliTPythia8::Instance()->ReadString("443:onIfAll = 13 13");
284 //      Psi'-> mu+ mu-
285         AliTPythia8::Instance()->ReadString("100443:onMode = off");
286         AliTPythia8::Instance()->ReadString("100443:onIfAll = 13 13");
287 //      Ups-> mu+ mu-
288         AliTPythia8::Instance()->ReadString("553:onMode = off");
289         AliTPythia8::Instance()->ReadString("553:onIfAll = 13 13");
290 //      Ups'-> mu+ mu-
291         AliTPythia8::Instance()->ReadString("100553:onMode = off");
292         AliTPythia8::Instance()->ReadString("100553:onIfAll = 13 13"); 
293 //      Ups''-> mu+ mu-
294         AliTPythia8::Instance()->ReadString("200553:onMode = off");
295         AliTPythia8::Instance()->ReadString("200553:onIfAll = 13 13");
296         break;
297     case kBSemiElectronic:
298 //      B0 - > e+ e-
299         AliTPythia8::Instance()->ReadString("511:onMode = off");
300         AliTPythia8::Instance()->ReadString("511:onIfAny = 11");
301 //      B+- -> e+ e-
302         AliTPythia8::Instance()->ReadString("521:onMode = off");
303         AliTPythia8::Instance()->ReadString("521:onIfAny = 11");
304 //      B_s -> e+ e-
305         AliTPythia8::Instance()->ReadString("531:onMode = off");
306         AliTPythia8::Instance()->ReadString("531:onIfAny = 11");
307 //      Lambda_b -> e+ e-
308         AliTPythia8::Instance()->ReadString("5122:onMode = off");
309         AliTPythia8::Instance()->ReadString("5122:onIfAny = 11");
310 //      Sigma_b -> e+ e-
311         AliTPythia8::Instance()->ReadString("5132:onMode = off");
312         AliTPythia8::Instance()->ReadString("5132:onIfAny = 11");
313 //      Sigma_b -> e+ e-
314         AliTPythia8::Instance()->ReadString("5232:onMode = off");
315         AliTPythia8::Instance()->ReadString("5232:onIfAny = 11");
316 //      Omega_b ->e+ e-
317         AliTPythia8::Instance()->ReadString("5332:onMode = off");
318         AliTPythia8::Instance()->ReadString("5332:onIfAny = 11");
319         break;
320     case kSemiElectronic:
321 //      D+/- -> e X
322         AliTPythia8::Instance()->ReadString("411:onMode = off");
323         AliTPythia8::Instance()->ReadString("411:onIfAll = 11");
324 //      D0 -> e X
325         AliTPythia8::Instance()->ReadString("421:onMode = off");
326         AliTPythia8::Instance()->ReadString("421:onIfAll = 11");
327 //      D_s ->e X
328         AliTPythia8::Instance()->ReadString("431:onMode = off");
329         AliTPythia8::Instance()->ReadString("431:onIfAll = 11");
330 //      Lambda_c -> e X
331         AliTPythia8::Instance()->ReadString("4122:onMode = off");
332         AliTPythia8::Instance()->ReadString("4122:onIfAll = 11");
333 //      Sigma_c -> e X
334         AliTPythia8::Instance()->ReadString("4132:onMode = off");
335         AliTPythia8::Instance()->ReadString("4132:onIfAll = 11");
336 //      Sigma_c -> e X
337         AliTPythia8::Instance()->ReadString("4232:onMode = off");
338         AliTPythia8::Instance()->ReadString("4232:onIfAll = 11");
339 //      Omega_c -> e X
340         AliTPythia8::Instance()->ReadString("4332:onMode = off");
341         AliTPythia8::Instance()->ReadString("4332:onIfAll = 11");
342 //      B0 -> e X
343         AliTPythia8::Instance()->ReadString("511:onMode = off");
344         AliTPythia8::Instance()->ReadString("511:onIfAny = 11");
345 //      B+/- -> e X
346         AliTPythia8::Instance()->ReadString("521:onMode = off");
347         AliTPythia8::Instance()->ReadString("521:onIfAny = 11");
348 //      B_s -> e X
349         AliTPythia8::Instance()->ReadString("531:onMode = off");
350         AliTPythia8::Instance()->ReadString("531:onIfAny = 11");
351 //      Lambda_b -> e X
352         AliTPythia8::Instance()->ReadString("5122:onMode = off");
353         AliTPythia8::Instance()->ReadString("5122:onIfAny = 11");
354 //      Sigma_b -> e X
355         AliTPythia8::Instance()->ReadString("5132:onMode = off");
356         AliTPythia8::Instance()->ReadString("5132:onIfAny = 11");
357 //      Sigma_b -> e X
358         AliTPythia8::Instance()->ReadString("5232:onMode = off");
359         AliTPythia8::Instance()->ReadString("5232:onIfAny = 11");
360 //      Omega_b -> e X
361         AliTPythia8::Instance()->ReadString("5332:onMode = off");
362         AliTPythia8::Instance()->ReadString("5332:onIfAny = 11");
363         break;
364     case kDiElectron:
365 //      Rho -> e+e-
366         AliTPythia8::Instance()->ReadString("113:onMode = off");
367         AliTPythia8::Instance()->ReadString("113:onIfAll = 11 11");
368 //      Eta -> e+e-
369         AliTPythia8::Instance()->ReadString("221:onMode = off");
370         AliTPythia8::Instance()->ReadString("221:onIfAll = 11 11");
371 //      omega -> e+e-
372         AliTPythia8::Instance()->ReadString("223:onMode = off");
373         AliTPythia8::Instance()->ReadString("223:onIfAll = 11 11");
374 //      phi -> e+e-
375         AliTPythia8::Instance()->ReadString("333:onMode = off");
376         AliTPythia8::Instance()->ReadString("333:onIfAll = 11 11");
377 //      J/Psi -> e+e-
378         AliTPythia8::Instance()->ReadString("443:onMode = off");
379         AliTPythia8::Instance()->ReadString("443:onIfAll = 11 11");
380 //      Psi' -> e+e-
381         AliTPythia8::Instance()->ReadString("100443:onMode = off");
382         AliTPythia8::Instance()->ReadString("100443:onIfAll = 11 11");
383 //      Ups -> e+e-
384         AliTPythia8::Instance()->ReadString("553:onMode = off");
385         AliTPythia8::Instance()->ReadString("553:onIfAll = 11 11");
386 //      Ups' -> e+e-
387         AliTPythia8::Instance()->ReadString("100553:onMode = off");
388         AliTPythia8::Instance()->ReadString("100553:onIfAll = 11 11");
389 //      Ups'' -> e+e-
390         AliTPythia8::Instance()->ReadString("200553:onMode = off");
391         AliTPythia8::Instance()->ReadString("200553:onIfAll = 11 11");
392         break;
393     case kBJpsiDiMuon:
394 //      B0   -> J/Psi (Psi') X   
395         AliTPythia8::Instance()->ReadString("511:onMode = off");
396         AliTPythia8::Instance()->ReadString("511:onIfAny = 443 100443");
397 //      B+/-   -> J/Psi (Psi') X   
398         AliTPythia8::Instance()->ReadString("521:onMode = off");
399         AliTPythia8::Instance()->ReadString("521:onIfAny = 443 100443");
400 //      B_s   -> J/Psi (Psi') X   
401         AliTPythia8::Instance()->ReadString("531:onMode = off");
402         AliTPythia8::Instance()->ReadString("531:onIfAny = 443 100443");
403 //      Lambda_b -> J/Psi (Psi') X   
404         AliTPythia8::Instance()->ReadString("5122:onMode = off");
405         AliTPythia8::Instance()->ReadString("5122:onIfAny = 443 100443");
406 //
407 //      J/Psi -> mu+ mu-
408         AliTPythia8::Instance()->ReadString("443:onMode = off");
409         AliTPythia8::Instance()->ReadString("443:onIfAll = 13 13");
410 //      Psi' -> mu+ mu-
411         AliTPythia8::Instance()->ReadString("100443:onMode = off");
412         AliTPythia8::Instance()->ReadString("100443:onIfAll = 13 13");
413         break;
414     case kBPsiPrimeDiMuon:
415 //      B0   -> Psi' X   
416         AliTPythia8::Instance()->ReadString("511:onMode = off");
417         AliTPythia8::Instance()->ReadString("511:onIfAny = 100443");
418 //      B+/-   -> Psi' X   
419         AliTPythia8::Instance()->ReadString("521:onMode = off");
420         AliTPythia8::Instance()->ReadString("521:onIfAny = 100443");
421 //      B_s   -> Psi'  X   
422         AliTPythia8::Instance()->ReadString("531:onMode = off");
423         AliTPythia8::Instance()->ReadString("531:onIfAny = 100443");
424 //      Lambda_b -> Psi' X   
425         AliTPythia8::Instance()->ReadString("5122:onMode = off");
426         AliTPythia8::Instance()->ReadString("5122:onIfAny = 100443");
427 //
428 //      Psi' -> mu+ mu-
429         AliTPythia8::Instance()->ReadString("100443:onMode = off");
430         AliTPythia8::Instance()->ReadString("100443:onIfAll = 13 13");
431         break;
432     case kBJpsiDiElectron:
433 //      B0   -> Psi X   
434         AliTPythia8::Instance()->ReadString("511:onMode = off");
435         AliTPythia8::Instance()->ReadString("511:onIfAny = 443");
436 //      B+/-   -> Psi X   
437         AliTPythia8::Instance()->ReadString("521:onMode = off");
438         AliTPythia8::Instance()->ReadString("521:onIfAny = 443");
439 //      B_s   -> Psi  X   
440         AliTPythia8::Instance()->ReadString("531:onMode = off");
441         AliTPythia8::Instance()->ReadString("531:onIfAny = 443");
442 //      Lambda_b -> Psi X   
443         AliTPythia8::Instance()->ReadString("5122:onMode = off");
444         AliTPythia8::Instance()->ReadString("5122:onIfAny = 443");
445 //
446 //      Psi -> mu+ mu-
447         AliTPythia8::Instance()->ReadString("443:onMode = off");
448         AliTPythia8::Instance()->ReadString("443:onIfAll = 11 11");
449
450         break;
451     case kBJpsi:
452 //      B0   -> Psi X   
453         AliTPythia8::Instance()->ReadString("511:onMode = off");
454         AliTPythia8::Instance()->ReadString("511:onIfAny = 443");
455 //      B+/-   -> Psi X   
456         AliTPythia8::Instance()->ReadString("521:onMode = off");
457         AliTPythia8::Instance()->ReadString("521:onIfAny = 443");
458 //      B_s   -> Psi  X   
459         AliTPythia8::Instance()->ReadString("531:onMode = off");
460         AliTPythia8::Instance()->ReadString("531:onIfAny = 443");
461 //      Lambda_b -> Psi X   
462         AliTPythia8::Instance()->ReadString("5122:onMode = off");
463         AliTPythia8::Instance()->ReadString("5122:onIfAny = 443");
464         break;
465     case kBPsiPrimeDiElectron:
466 //      B0   -> Psi' X   
467         AliTPythia8::Instance()->ReadString("511:onMode = off");
468         AliTPythia8::Instance()->ReadString("511:onIfAny = 100443");
469 //      B+/-   -> Psi' X   
470         AliTPythia8::Instance()->ReadString("521:onMode = off");
471         AliTPythia8::Instance()->ReadString("521:onIfAny = 100443");
472 //      B_s   -> Psi'  X   
473         AliTPythia8::Instance()->ReadString("531:onMode = off");
474         AliTPythia8::Instance()->ReadString("531:onIfAny = 100443");
475 //      Lambda_b -> Psi' X   
476         AliTPythia8::Instance()->ReadString("5122:onMode = off");
477         AliTPythia8::Instance()->ReadString("5122:onIfAny = 100443");
478 //
479 //      Psi' -> mu+ mu-
480         AliTPythia8::Instance()->ReadString("100443:onMode = off");
481         AliTPythia8::Instance()->ReadString("100443:onIfAll = 11 11");
482         break;
483     case kPiToMu:
484 //      pi -> mu nu
485         AliTPythia8::Instance()->ReadString("211:onMode = off");
486         AliTPythia8::Instance()->ReadString("211:onIfAny = 13");
487         break;
488     case kKaToMu:
489 //      K -> mu nu
490         AliTPythia8::Instance()->ReadString("321:onMode = off");
491         AliTPythia8::Instance()->ReadString("321:onIfAny = 13");
492         break;
493     case kAllMuonic:
494 //      pi/K -> mu
495         AliTPythia8::Instance()->ReadString("211:onMode = off");
496         AliTPythia8::Instance()->ReadString("211:onIfAny = 13");
497         AliTPythia8::Instance()->ReadString("321:onMode = off");
498         AliTPythia8::Instance()->ReadString("321:onIfAny = 13");
499         break;
500     case kWToMuon:
501 //      W -> mu X
502         AliTPythia8::Instance()->ReadString("24:onMode = off");
503         AliTPythia8::Instance()->ReadString("24:onIfAny = 13");
504         break;
505     case kWToCharm:
506 //      W -> c X
507         AliTPythia8::Instance()->ReadString("24:onMode = off");
508         AliTPythia8::Instance()->ReadString("24:onIfAny = 4");
509         break;
510     case kWToCharmToMuon:
511 //      W -> c X
512         AliTPythia8::Instance()->ReadString("24:onMode = off");
513         AliTPythia8::Instance()->ReadString("24:onIfAny = 4");
514 //      D+- -> mu X
515         AliTPythia8::Instance()->ReadString("411:onMode = off");
516         AliTPythia8::Instance()->ReadString("411:onIfAll = 13");
517 //      D0 -> mu X
518         AliTPythia8::Instance()->ReadString("421:onMode = off");
519         AliTPythia8::Instance()->ReadString("421:onIfAll = 13");
520 //      D_s -> mu X
521         AliTPythia8::Instance()->ReadString("431:onMode = off");
522         AliTPythia8::Instance()->ReadString("431:onIfAll = 13");
523 //      Lambda_c -> mu X
524         AliTPythia8::Instance()->ReadString("4122:onMode = off");
525         AliTPythia8::Instance()->ReadString("4122:onIfAll = 13");
526 //      Sigma_c -> mu X
527         AliTPythia8::Instance()->ReadString("4132:onMode = off");
528         AliTPythia8::Instance()->ReadString("4132:onIfAll = 13");
529 //      Sigma_c -> mu X
530         AliTPythia8::Instance()->ReadString("4232:onMode = off");
531         AliTPythia8::Instance()->ReadString("4232:onIfAll = 13");
532 //      Omega_c -> mu X
533         AliTPythia8::Instance()->ReadString("4332:onMode = off");
534         AliTPythia8::Instance()->ReadString("4332:onIfAll = 13");
535         break;
536     case kZDiMuon:
537 //      Z -> mu+ mu-
538         AliTPythia8::Instance()->ReadString("23:onMode = off");
539         AliTPythia8::Instance()->ReadString("23:onIfAll = 13 13");
540         break;
541     case kZDiElectron:
542 //      Z -> e+ e-
543         AliTPythia8::Instance()->ReadString("23:onMode = off");
544         AliTPythia8::Instance()->ReadString("23:onIfAll = 11 11");
545         break;
546     case kHadronicD:
547         ForceHadronicD(1);
548         break;
549     case kHadronicDWithout4Bodies:
550         ForceHadronicD(0);
551         break;
552     case kPhiKK:
553         // Phi-> K+ K-
554         AliTPythia8::Instance()->ReadString("333:onMode = off");
555         AliTPythia8::Instance()->ReadString("333:onIfAll = 321 321");
556         break;
557     case kOmega:
558         // Omega -> Lambda K
559         AliTPythia8::Instance()->ReadString("3334:onMode = off");
560         AliTPythia8::Instance()->ReadString("3334:onIfAll = 3122 321 ");
561         break;
562     case kLambda:
563         // Lambda -> p pi-
564         AliTPythia8::Instance()->ReadString("3122:onMode = off");
565         AliTPythia8::Instance()->ReadString("3122:onIfAll = 2212 211 ");
566         break;
567     case kAll:
568         break;
569     case kNoDecay:
570         AliTPythia8::Instance()->ReadString("HadronLevel:Decay = off");
571         break;
572     case kNoDecayHeavy:
573     case kNoDecayBeauty:
574     case kNeutralPion:
575     case kPsiPrimeJpsiDiElectron: 
576     case kElectronEM:
577     case kGammaEM:
578     case kDiElectronEM:
579         break;
580     }
581 }
582
583 Float_t AliDecayerPythia8::GetPartialBranchingRatio(Int_t ipart)
584 {
585     // Get the partial branching ration for the forced decay channels
586     
587     Pythia8::Pythia* thePythia       = AliTPythia8::Instance()->Pythia8();
588     Pythia8::ParticleData & table    = thePythia->particleData;
589     Pythia8::ParticleDataEntry* pd   = table.particleDataEntryPtr(ipart);
590
591     Int_t nc = pd->sizeChannels();
592     Float_t br = 0.;
593 //
594 //  Loop over decay channels
595     for (Int_t ic = 0; ic < nc; ic++) {
596       Pythia8::DecayChannel& decCh = pd->channel(ic);
597         for (Int_t i = 0; i < decCh.multiplicity(); i++) {
598             br += decCh.bRatio();
599         }
600     }
601     return (br);
602 }
603
604
605 Float_t AliDecayerPythia8::GetLifetime(Int_t kf)
606 {
607     // Return lifetime of particle
608     Pythia8::Pythia* thePythia       = AliTPythia8::Instance()->Pythia8();
609     Pythia8::ParticleData& table     = thePythia->particleData;
610     Float_t tau = table.tau0(kf);
611     return (tau);
612 }
613
614 void  AliDecayerPythia8::SwitchOffHeavyFlavour()
615 {
616     // Switch off heavy flavour production
617     //
618 // Maximum number of quark flavours used in pdf 
619     AliTPythia8::Instance()->ReadString("PDFinProcess:nQuarkIn = 3");
620 // Maximum number of flavors that can be used in showers
621     AliTPythia8::Instance()->ReadString("SpaceShower:nQuarkIn = 3");
622     AliTPythia8::Instance()->ReadString("TimeShower:nGammaToQuark = 3");
623     AliTPythia8::Instance()->ReadString("TimeShower:nGluonToQuark = 3");
624 }
625
626
627 void AliDecayerPythia8::ForceHadronicD(Int_t optUse4Bodies)
628 {
629 //
630 // Force golden D decay modes
631 //
632     // K* -> K pi
633     AliTPythia8::Instance()->ReadString("313:onMode = off");
634     AliTPythia8::Instance()->ReadString("313:onIfAll = 321 211");
635     // for Ds -> Phi pi+
636     AliTPythia8::Instance()->ReadString("333:onMode = off");
637     AliTPythia8::Instance()->ReadString("333:onIfAll = 321 321");
638     // for D0 -> rho0 pi+ k-
639     AliTPythia8::Instance()->ReadString("113:onMode = off");
640     AliTPythia8::Instance()->ReadString("113:onIfAll = 211 211");
641     // for Lambda_c -> Delta++ K-
642     AliTPythia8::Instance()->ReadString("2224:onMode = off");
643     AliTPythia8::Instance()->ReadString("2224:onIfAll = 2212 211");
644     // for Lambda_c -> Lambda(1520) K-
645     AliTPythia8::Instance()->ReadString("3124:onMode = off");
646     AliTPythia8::Instance()->ReadString("3124:onIfAll = 2212 321");
647
648
649     AliTPythia8::Instance()->ReadString("411:onMode = off");
650     AliTPythia8::Instance()->ReadString("421:onMode = off");
651     AliTPythia8::Instance()->ReadString("431:onMode = off");
652     AliTPythia8::Instance()->ReadString("4112:onMode = off");
653     AliTPythia8::Instance()->ReadString("4122:onMode = off");
654
655     // D+/- -> K pi pi 
656     AliTPythia8::Instance()->ReadString("411:onIfMatch = 321 211 211");
657     // D+/- -> K* pi
658     AliTPythia8::Instance()->ReadString("411:onIfMatch = 313 211");
659     // D0 -> K pi
660     AliTPythia8::Instance()->ReadString("421:onIfMatch = 321 211");
661
662     if (optUse4Bodies) {
663         // D0 -> K pi pi pi
664         AliTPythia8::Instance()->ReadString("421:onIfMatch = 321 211 211 211");
665         // D0 -> K pi rho
666         AliTPythia8::Instance()->ReadString("421:onIfMatch = 321 211 113");
667         // D0 -> K*0 pi pi
668         AliTPythia8::Instance()->ReadString("421:onIfMatch = 313 211 211");
669     }
670     
671     // D_s -> K K*
672     AliTPythia8::Instance()->ReadString("431:onIfMatch = 321 313");
673     // D_s -> Phi pi
674     AliTPythia8::Instance()->ReadString("431:onIfMatch = 333 211");
675
676     // Lambda_c -> p K*
677     AliTPythia8::Instance()->ReadString("4122:onIfMatch = 2212 313");
678     // Lambda_c -> Delta K
679     AliTPythia8::Instance()->ReadString("4122:onIfMatch = 2224 321");
680     // Lambda_c -> Lambda(1520) pi
681     AliTPythia8::Instance()->ReadString("4122:onIfMatch = 3124 211");
682     // Lambda_c -> p K pi
683     AliTPythia8::Instance()->ReadString("4122:onIfMatch = 2212 321 211");
684     // Lambda_c -> Lambda pi 
685     AliTPythia8::Instance()->ReadString("4122:onIfMatch = 3122 211");
686
687 }
688
689 //___________________________________________________________________________
690 void    AliDecayerPythia8::ReadDecayTable()
691 {
692    //to read a decay table (not yet implemented)
693 }
694
695
696 //___________________________________________________________________________
697 void AliDecayerPythia8::AppendParticle(Int_t pdg, TLorentzVector* p)
698 {
699    // Append a particle to the stack
700    fPythia8->Pythia8()->event.append(pdg, 11, 0, 0, p->Px(), p->Py(), p->Pz(), p->E(), p->M());   
701 }
702
703
704 //___________________________________________________________________________
705 void AliDecayerPythia8::ClearEvent()
706 {
707    // Clear the event stack
708    fPythia8->Pythia8()->event.clear();
709 }
710