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