]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PYTHIA6/AliPythia6.cxx
Ptotections added. (M. van Leeuwen)
[u/mrichter/AliRoot.git] / PYTHIA6 / AliPythia6.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id: AliPythia.cxx,v 1.40 2007/10/09 08:43:24 morsch Exp $ */
18
19 #include "AliPythia6.h"
20 #include "AliStack.h"
21 #include "AliPythiaRndm.h"
22 #include "AliFastGlauber.h"
23 #include "AliQuenchingWeights.h"
24
25 #include "TVector3.h"
26 #include "TParticle.h"
27 #include "PyquenCommon.h"
28
29 ClassImp(AliPythia6)
30
31 #ifndef WIN32
32 # define pyclus pyclus_
33 # define pycell pycell_
34 # define pyshow pyshow_
35 # define pyrobo pyrobo_
36 # define pyquen pyquen_
37 # define pyevnw pyevnw_
38 # define type_of_call
39 #else
40 # define pyclus PYCLUS
41 # define pycell PYCELL
42 # define pyrobo PYROBO
43 # define pyquen PYQUEN
44 # define pyevnw PYEVNW
45 # define type_of_call _stdcall
46 #endif
47
48 extern "C" void type_of_call pyclus(Int_t & );
49 extern "C" void type_of_call pycell(Int_t & );
50 extern "C" void type_of_call pyshow(Int_t &, Int_t &, Double_t &);
51 extern "C" void type_of_call pyrobo(Int_t &, Int_t &, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &);
52 extern "C" void type_of_call pyquen(Double_t &, Int_t &, Double_t &);
53 extern "C" void type_of_call pyevnw();
54
55
56 //_____________________________________________________________________________
57
58 AliPythia6* AliPythia6::fgAliPythia=NULL;
59
60 AliPythia6::AliPythia6():
61     TPythia6(),
62     AliPythiaBase(),
63     fProcess(kPyMb),
64     fEcms(0.),
65     fStrucFunc(kCTEQ5L),
66     fXJet(0.),
67     fYJet(0.),
68     fNGmax(30),
69     fZmax(0.97),
70     fGlauber(0),
71     fQuenchingWeights(0)
72 {
73 // Default Constructor
74 //
75 //  Set random number
76     if (!AliPythiaRndm::GetPythiaRandom()) 
77       AliPythiaRndm::SetPythiaRandom(GetRandom());
78     fGlauber          = 0;
79     fQuenchingWeights = 0;
80 }
81
82 AliPythia6::AliPythia6(const AliPythia6& pythia):
83     TPythia6(),
84     AliPythiaBase(),
85     fProcess(kPyMb),
86     fEcms(0.),
87     fStrucFunc(kCTEQ5L),
88     fXJet(0.),
89     fYJet(0.),
90     fNGmax(30),
91     fZmax(0.97),
92     fGlauber(0),
93     fQuenchingWeights(0)
94 {
95     // Copy Constructor
96     pythia.Copy(*this);
97 }
98
99 void AliPythia6::ProcInit(Process_t process, Float_t energy, StrucFunc_t strucfunc)
100 {
101 // Initialise the process to generate 
102     if (!AliPythiaRndm::GetPythiaRandom()) 
103       AliPythiaRndm::SetPythiaRandom(GetRandom());
104     
105     fProcess = process;
106     fEcms = energy;
107     fStrucFunc = strucfunc;
108 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
109     SetMDCY(Pycomp(111) ,1,0);
110     SetMDCY(Pycomp(310) ,1,0);
111     SetMDCY(Pycomp(3122),1,0);
112     SetMDCY(Pycomp(3112),1,0);
113     SetMDCY(Pycomp(3212),1,0);
114     SetMDCY(Pycomp(3222),1,0);
115     SetMDCY(Pycomp(3312),1,0);
116     SetMDCY(Pycomp(3322),1,0);
117     SetMDCY(Pycomp(3334),1,0);
118     // Select structure function 
119     SetMSTP(52,2);
120     SetMSTP(51,AliStructFuncType::PDFsetIndex(strucfunc));
121     // Particles produced in string fragmentation point directly to either of the two endpoints
122     // of the string (depending in the side they were generated from).
123     SetMSTU(16,2);
124
125 //
126 // Pythia initialisation for selected processes//
127 //
128 // Make MSEL clean
129 //
130     for (Int_t i=1; i<= 200; i++) {
131         SetMSUB(i,0);
132     }
133 //  select charm production
134     switch (process) 
135     {
136     case kPyOldUEQ2ordered:  //Old underlying events with Q2 ordered QCD processes
137 //        Multiple interactions on.
138         SetMSTP(81,1);
139 // Double Gaussian matter distribution.
140         SetMSTP(82,4);
141         SetPARP(83,0.5);
142         SetPARP(84,0.4);
143 //  pT0.
144         SetPARP(82,2.0);
145 //  Reference energy for pT0 and energy rescaling pace.
146         SetPARP(89,1800);
147         SetPARP(90,0.25);
148 //  String drawing almost completely minimizes string length.
149         SetPARP(85,0.9);
150         SetPARP(86,0.95);
151 // ISR and FSR activity.
152         SetPARP(67,4);
153         SetPARP(71,4);
154 // Lambda_FSR scale.
155         SetPARJ(81,0.29);
156         break;
157     case kPyOldUEQ2ordered2:   
158 // Old underlying events with Q2 ordered QCD processes
159 // Multiple interactions on.
160         SetMSTP(81,1);
161 // Double Gaussian matter distribution.
162         SetMSTP(82,4);
163         SetPARP(83,0.5);
164         SetPARP(84,0.4);
165 // pT0.
166         SetPARP(82,2.0);
167 // Reference energy for pT0 and energy rescaling pace.
168         SetPARP(89,1800);
169         SetPARP(90,0.16);  // here is the difference with  kPyOldUEQ2ordered
170 // String drawing almost completely minimizes string length.
171         SetPARP(85,0.9);
172         SetPARP(86,0.95);
173 // ISR and FSR activity.
174         SetPARP(67,4);
175         SetPARP(71,4);
176 // Lambda_FSR scale.
177         SetPARJ(81,0.29);       
178         break;
179     case kPyOldPopcorn:  
180 // Old production mechanism: Old Popcorn
181         SetMSEL(1);
182         SetMSTJ(12,3); 
183 // (D=2) Like MSTJ(12)=2 but added prod ofthe 1er rank baryon
184         SetMSTP(88,2); 
185 // (D=1)see can be used to form  baryons (BARYON JUNCTION)
186         SetMSTJ(1,1);  
187         AtlasTuning();
188         break;
189     case kPyCharm:
190         SetMSEL(4);
191 //  heavy quark masses
192
193         SetPMAS(4,1,1.2);
194 //
195 //    primordial pT
196         SetMSTP(91,1);
197         SetPARP(91,1.);
198         SetPARP(93,5.);
199 //
200         break;
201     case kPyBeauty:
202         SetMSEL(5);
203         SetPMAS(5,1,4.75);
204         break;
205     case kPyJpsi:
206         SetMSEL(0);
207 // gg->J/Psi g
208         SetMSUB(86,1);
209         break;
210     case kPyJpsiChi:
211         SetMSEL(0);
212 // gg->J/Psi g
213         SetMSUB(86,1);
214 // gg-> chi_0c g
215         SetMSUB(87,1);
216 // gg-> chi_1c g
217         SetMSUB(88,1);
218 // gg-> chi_2c g
219         SetMSUB(89,1);  
220         break;
221     case kPyCharmUnforced:
222         SetMSEL(0);
223 // gq->qg   
224         SetMSUB(28,1);
225 // gg->qq
226         SetMSUB(53,1);
227 // gg->gg
228         SetMSUB(68,1);
229         break;
230     case kPyBeautyUnforced:
231         SetMSEL(0);
232 // gq->qg   
233         SetMSUB(28,1);
234 // gg->qq
235         SetMSUB(53,1);
236 // gg->gg
237         SetMSUB(68,1);
238         break;
239     case kPyMb:
240 // Minimum Bias pp-Collisions
241 //
242 //   
243 //      select Pythia min. bias model
244         SetMSEL(0);
245         SetMSUB(92,1);             // single diffraction AB-->XB
246         SetMSUB(93,1);             // single diffraction AB-->AX
247         SetMSUB(94,1);             // double diffraction
248         SetMSUB(95,1);             // low pt production
249
250         AtlasTuning();
251         break;
252     case kPyMbDefault:
253 // Minimum Bias pp-Collisions
254 //
255 //   
256 //      select Pythia min. bias model
257         SetMSEL(0);
258         SetMSUB(92,1);             // single diffraction AB-->XB
259         SetMSUB(93,1);             // single diffraction AB-->AX
260         SetMSUB(94,1);             // double diffraction
261         SetMSUB(95,1);             // low pt production
262
263         break;
264     case kPyLhwgMb:
265 // Les Houches Working Group 05 Minimum Bias pp-Collisions: hep-ph/0604120
266 //  -> Pythia 6.3 or above is needed
267 //   
268         SetMSEL(0);
269         SetMSUB(92,1);             // single diffraction AB-->XB
270         SetMSUB(93,1);             // single diffraction AB-->AX
271         SetMSUB(94,1);             // double diffraction
272         SetMSUB(95,1);             // low pt production
273         SetMSTP(51,AliStructFuncType::PDFsetIndex(kCTEQ6ll));
274         SetMSTP(52,2);
275         SetMSTP(68,1);
276         SetMSTP(70,2);
277         SetMSTP(81,1);             // Multiple Interactions ON
278         SetMSTP(82,4);             // Double Gaussian Model
279         SetMSTP(88,1);
280
281         SetPARP(82,2.3);           // [GeV]    PT_min at Ref. energy
282         SetPARP(83,0.5);           // Core density in proton matter distribution (def.value)
283         SetPARP(84,0.5);           // Core radius
284         SetPARP(85,0.9);           // Regulates gluon prod. mechanism
285         SetPARP(90,0.2);           // 2*epsilon (exponent in power law)
286
287         break;
288     case kPyMbNonDiffr:
289 // Minimum Bias pp-Collisions
290 //
291 //   
292 //      select Pythia min. bias model
293         SetMSEL(0);
294         SetMSUB(95,1);             // low pt production
295
296         AtlasTuning();
297         break;
298     case kPyMbMSEL1:
299         ConfigHeavyFlavor();
300 // Intrinsic <kT^2>
301         SetMSTP(91,1);// Width (1=gaussian) primordial kT dist. inside hadrons
302         SetPARP(91,1.);     // <kT^2> = PARP(91,1.)^2
303         SetPARP(93,5.);     // Upper cut-off
304 // Set Q-quark mass
305         SetPMAS(4,1,1.2);   // Charm quark mass
306         SetPMAS(5,1,4.78);  // Beauty quark mass
307         SetPARP(71,4.);     // Defaut value
308 // Atlas Tuning
309         AtlasTuning();
310         break;
311     case kPyJets:
312 //
313 //  QCD Jets
314 //
315         SetMSEL(1);
316  // Pythia Tune A (CDF)
317  //
318        SetPARP(67,2.5);           // Regulates Initial State Radiation (value from best fit to D0 dijet analysis)
319        SetMSTP(82,4);             // Double Gaussian Model
320        SetPARP(82,2.0);           // [GeV]    PT_min at Ref. energy
321        SetPARP(84,0.4);           // Core radius
322        SetPARP(85,0.90) ;         // Regulates gluon prod. mechanism
323        SetPARP(86,0.95);          // Regulates gluon prod. mechanism
324        SetPARP(89,1800.);         // [GeV]   Ref. energy
325        SetPARP(90,0.25);          // 2*epsilon (exponent in power law)
326        break;
327     case kPyDirectGamma:
328         SetMSEL(10);
329         break;
330     case kPyCharmPbPbMNR:
331     case kPyD0PbPbMNR:
332     case kPyDPlusPbPbMNR:
333     case kPyDPlusStrangePbPbMNR:
334       // Tuning of Pythia parameters aimed to get a resonable agreement
335       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
336       // c-cbar single inclusive and double differential distributions.
337       // This parameter settings are meant to work with Pb-Pb collisions
338       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
339       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
340       // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
341         ConfigHeavyFlavor();
342       // Intrinsic <kT>
343       SetMSTP(91,1);
344       SetPARP(91,1.304);
345       SetPARP(93,6.52);
346       // Set c-quark mass
347       SetPMAS(4,1,1.2);
348       break;
349     case kPyCharmpPbMNR:
350     case kPyD0pPbMNR:
351     case kPyDPluspPbMNR:
352     case kPyDPlusStrangepPbMNR:
353       // Tuning of Pythia parameters aimed to get a resonable agreement
354       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
355       // c-cbar single inclusive and double differential distributions.
356       // This parameter settings are meant to work with p-Pb collisions
357       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
358       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
359       // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
360         ConfigHeavyFlavor();
361       // Intrinsic <kT>
362         SetMSTP(91,1);
363         SetPARP(91,1.16);
364         SetPARP(93,5.8);
365         
366       // Set c-quark mass
367         SetPMAS(4,1,1.2);
368       break;
369     case kPyCharmppMNR:
370     case kPyD0ppMNR:
371     case kPyDPlusppMNR:
372     case kPyDPlusStrangeppMNR:
373       // Tuning of Pythia parameters aimed to get a resonable agreement
374       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
375       // c-cbar single inclusive and double differential distributions.
376       // This parameter settings are meant to work with pp collisions
377       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
378       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
379       // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
380         ConfigHeavyFlavor();
381       // Intrinsic <kT^2>
382         SetMSTP(91,1);
383         SetPARP(91,1.);
384         SetPARP(93,5.);
385         
386       // Set c-quark mass
387         SetPMAS(4,1,1.2);
388       break;
389     case kPyCharmppMNRwmi:
390       // Tuning of Pythia parameters aimed to get a resonable agreement
391       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
392       // c-cbar single inclusive and double differential distributions.
393       // This parameter settings are meant to work with pp collisions
394       // and with kCTEQ5L PDFs.
395       // Added multiple interactions according to ATLAS tune settings.
396       // To get a "reasonable" agreement with MNR results, events have to be 
397       // generated with the minimum ptHard (AliGenPythia::SetPtHard)
398       // set to 2.76 GeV.
399       // To get a "perfect" agreement with MNR results, events have to be 
400       // generated in four ptHard bins with the following relative 
401       // normalizations:
402       // 2.76-3 GeV: 25%
403       //    3-4 GeV: 40%
404       //    4-8 GeV: 29%
405       //     >8 GeV:  6%
406         ConfigHeavyFlavor();
407       // Intrinsic <kT^2>
408         SetMSTP(91,1);
409         SetPARP(91,1.);
410         SetPARP(93,5.);
411
412       // Set c-quark mass
413         SetPMAS(4,1,1.2);
414         AtlasTuning();
415         break;
416     case kPyBeautyPbPbMNR:
417       // Tuning of Pythia parameters aimed to get a resonable agreement
418       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
419       // b-bbar single inclusive and double differential distributions.
420       // This parameter settings are meant to work with Pb-Pb collisions
421       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
422       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
423       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
424         ConfigHeavyFlavor();
425       // QCD scales
426         SetPARP(67,1.0);
427         SetPARP(71,1.0);
428       // Intrinsic <kT>
429         SetMSTP(91,1);
430         SetPARP(91,2.035);
431         SetPARP(93,10.17);
432       // Set b-quark mass
433         SetPMAS(5,1,4.75);
434       break;
435     case kPyBeautypPbMNR:
436       // Tuning of Pythia parameters aimed to get a resonable agreement
437       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
438       // b-bbar single inclusive and double differential distributions.
439       // This parameter settings are meant to work with p-Pb collisions
440       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
441       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
442       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
443         ConfigHeavyFlavor();
444       // QCD scales
445         SetPARP(67,1.0);
446         SetPARP(71,1.0);
447       // Intrinsic <kT>
448         SetMSTP(91,1);
449         SetPARP(91,1.60);
450         SetPARP(93,8.00);
451       // Set b-quark mass
452         SetPMAS(5,1,4.75);
453       break;
454     case kPyBeautyppMNR:
455       // Tuning of Pythia parameters aimed to get a resonable agreement
456       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
457       // b-bbar single inclusive and double differential distributions.
458       // This parameter settings are meant to work with pp collisions
459       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
460       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
461       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
462         ConfigHeavyFlavor();
463       // QCD scales
464         SetPARP(67,1.0);
465         SetPARP(71,1.0);
466         
467         // Intrinsic <kT>
468         SetMSTP(91,1);
469         SetPARP(91,1.);
470         SetPARP(93,5.);
471         
472         // Set b-quark mass
473         SetPMAS(5,1,4.75);
474       break;
475      case kPyBeautyppMNRwmi:
476       // Tuning of Pythia parameters aimed to get a resonable agreement
477       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
478       // b-bbar single inclusive and double differential distributions.
479       // This parameter settings are meant to work with pp collisions
480       // and with kCTEQ5L PDFs.
481       // Added multiple interactions according to ATLAS tune settings.
482       // To get a "reasonable" agreement with MNR results, events have to be 
483       // generated with the minimum ptHard (AliGenPythia::SetPtHard)
484       // set to 2.76 GeV.
485       // To get a "perfect" agreement with MNR results, events have to be 
486       // generated in four ptHard bins with the following relative 
487       // normalizations:
488       // 2.76-4 GeV:  5% 
489       //    4-6 GeV: 31%
490       //    6-8 GeV: 28%
491       //     >8 GeV: 36%
492          ConfigHeavyFlavor();
493       // QCD scales
494          SetPARP(67,1.0);
495          SetPARP(71,1.0);
496          
497          // Intrinsic <kT>
498          SetMSTP(91,1);
499          SetPARP(91,1.);
500          SetPARP(93,5.);
501
502       // Set b-quark mass
503          SetPMAS(5,1,4.75);
504
505          AtlasTuning();
506          break; 
507     case kPyW:
508
509       //Inclusive production of W+/-
510       SetMSEL(0);
511       //f fbar -> W+ 
512       SetMSUB(2,1);
513       //        //f fbar -> g W+
514       //        SetMSUB(16,1);
515       //        //f fbar -> gamma W+
516       //        SetMSUB(20,1);
517       //        //f g -> f W+  
518       //        SetMSUB(31,1);
519       //        //f gamma -> f W+
520       //        SetMSUB(36,1);
521       
522       // Initial/final parton shower on (Pythia default)
523       // With parton showers on we are generating "W inclusive process"
524       SetMSTP(61,1); //Initial QCD & QED showers on
525       SetMSTP(71,1); //Final QCD & QED showers on
526       
527       break;  
528
529     case kPyZ:
530
531       //Inclusive production of Z
532       SetMSEL(0);
533       //f fbar -> Z/gamma
534       SetMSUB(1,1);
535       
536       //       // f fbar -> g Z/gamma
537       //       SetMSUB(15,1);
538       //       // f fbar -> gamma Z/gamma
539       //       SetMSUB(19,1);
540       //       // f g -> f Z/gamma
541       //       SetMSUB(30,1);
542       //       // f gamma -> f Z/gamma
543       //       SetMSUB(35,1);
544       
545       //only Z included, not gamma
546       SetMSTP(43,2);
547       
548       // Initial/final parton shower on (Pythia default)
549       // With parton showers on we are generating "Z inclusive process"
550       SetMSTP(61,1); //Initial QCD & QED showers on
551       SetMSTP(71,1); //Final QCD & QED showers on
552       
553       break;  
554
555     }
556 //
557 //  Initialize PYTHIA
558     SetMSTP(41,1);   // all resonance decays switched on
559     Initialize("CMS","p","p",fEcms);
560     
561 }
562
563 Int_t AliPythia6::CheckedLuComp(Int_t kf)
564 {
565 // Check Lund particle code (for debugging)
566     Int_t kc=Pycomp(kf);
567     return kc;
568 }
569
570 void AliPythia6::SetNuclei(Int_t a1, Int_t a2)
571 {
572 // Treat protons as inside nuclei with mass numbers a1 and a2  
573 //    The MSTP array in the PYPARS common block is used to enable and 
574 //    select the nuclear structure functions. 
575 //    MSTP(52)  : (D=1) choice of proton and nuclear structure-function library
576 //            =1: internal PYTHIA acording to MSTP(51) 
577 //            =2: PDFLIB proton  s.f., with MSTP(51)  = 1000xNGROUP+NSET
578 //    If the following mass number both not equal zero, nuclear corrections of the stf are used.
579 //    MSTP(192) : Mass number of nucleus side 1
580 //    MSTP(193) : Mass number of nucleus side 2
581     SetMSTP(52,2);
582     SetMSTP(192, a1);
583     SetMSTP(193, a2);  
584 }
585         
586
587 AliPythia6* AliPythia6::Instance()
588
589 // Set random number generator 
590     if (fgAliPythia) {
591         return fgAliPythia;
592     } else {
593         fgAliPythia = new AliPythia6();
594         return fgAliPythia;
595     }
596 }
597
598 void AliPythia6::PrintParticles()
599
600 // Print list of particl properties
601     Int_t np = 0;
602     char*   name = new char[16];    
603     for (Int_t kf=0; kf<1000000; kf++) {
604         for (Int_t c = 1;  c > -2; c-=2) {
605             Int_t kc = Pycomp(c*kf);
606             if (kc) {
607                 Float_t mass  = GetPMAS(kc,1);
608                 Float_t width = GetPMAS(kc,2);  
609                 Float_t tau   = GetPMAS(kc,4);
610
611                 Pyname(kf,name);
612         
613                 np++;
614                 
615                 printf("\n mass, width, tau: %6d %s %10.3f %10.3e %10.3e", 
616                        c*kf, name, mass, width, tau);
617             }
618         }
619     }
620     printf("\n Number of particles %d \n \n", np);
621 }
622
623 void  AliPythia6::ResetDecayTable()
624 {
625 //  Set default values for pythia decay switches
626     Int_t i;
627     for (i = 1; i <  501; i++) SetMDCY(i,1,fDefMDCY[i]);
628     for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
629 }
630
631 void  AliPythia6::SetDecayTable()
632 {
633 //  Set default values for pythia decay switches
634 //
635     Int_t i;
636     for (i = 1; i <  501; i++) fDefMDCY[i] = GetMDCY(i,1);
637     for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
638 }
639
640 void  AliPythia6::Pyclus(Int_t& njet)
641 {
642 //  Call Pythia clustering algorithm
643 //
644     pyclus(njet);
645 }
646
647 void  AliPythia6::Pycell(Int_t& njet)
648 {
649 //  Call Pythia jet reconstruction algorithm
650 //
651     pycell(njet);
652 }
653
654 void AliPythia6::GetJet(Int_t i, Float_t& px, Float_t& py, Float_t& pz, Float_t& e)
655 {
656     // Get jet number i
657     Int_t n = GetN();
658     px    = GetPyjets()->P[0][n+i];
659     py    = GetPyjets()->P[1][n+i];
660     pz    = GetPyjets()->P[2][n+i];
661     e     = GetPyjets()->P[3][n+i];
662 }
663
664 void  AliPythia6::Pyshow(Int_t ip1, Int_t ip2, Double_t qmax)
665 {
666 //  Call Pythia showering
667 //
668     pyshow(ip1, ip2, qmax);
669 }
670
671 void AliPythia6::Pyrobo(Int_t imi, Int_t ima, Double_t the, Double_t phi, Double_t bex, Double_t bey, Double_t bez)
672 {
673     pyrobo(imi, ima, the, phi, bex, bey, bez);
674 }
675
676
677
678 void AliPythia6::InitQuenching(Float_t cMin, Float_t cMax, Float_t k, Int_t iECMethod, Float_t zmax, Int_t ngmax)
679 {
680 // Initializes 
681 // (1) The quenching model using quenching weights according to C. Salgado and U. Wiedemann
682 // (2) The nuclear geometry using the Glauber Model
683 //     
684     
685     fGlauber = new AliFastGlauber();
686     fGlauber->Init(2);
687     fGlauber->SetCentralityClass(cMin, cMax); 
688
689     fQuenchingWeights = new AliQuenchingWeights();
690     fQuenchingWeights->InitMult();
691     fQuenchingWeights->SetK(k);
692     fQuenchingWeights->SetECMethod(AliQuenchingWeights::kECMethod(iECMethod));
693     fNGmax = ngmax;
694     fZmax  = zmax;
695     
696 }
697
698
699 void  AliPythia6::Quench()
700 {
701 //
702 //
703 //  Simple Jet Quenching routine:
704 //  =============================
705 //  The jet formed by all final state partons radiated by the parton created 
706 //  in the hard collisions is quenched by a factor (1-z) using light cone variables in 
707 //  the initial parton reference frame:
708 //  (E + p_z)new = (1-z) (E + p_z)old
709 //
710 //
711 //
712 //
713 //   The lost momentum is first balanced by one gluon with virtuality > 0.   
714 //   Subsequently the gluon splits to yield two gluons with E = p.
715 //
716 //
717 // 
718     static Float_t eMean = 0.;
719     static Int_t   icall = 0;
720     
721     Double_t p0[4][5];
722     Double_t p1[4][5];
723     Double_t p2[4][5];
724     Int_t   klast[4] = {-1, -1, -1, -1};
725
726     Int_t numpart   = fPyjets->N;
727     Double_t px = 0., py = 0., pz = 0., e = 0., m = 0., p = 0., pt = 0., theta = 0., phi = 0.;
728     Double_t pxq[4], pyq[4], pzq[4], eq[4], yq[4], mq[4], pq[4], phiq[4], thetaq[4], ptq[4];
729     Bool_t  quenched[4];
730     Double_t wjtKick[4];
731     Int_t nGluon[4];
732     Int_t qPdg[4];
733     Int_t   imo, kst, pdg;
734     
735 //
736 //  Sore information about Primary partons
737 //
738 //  j =
739 //  0, 1 partons from hard scattering
740 //  2, 3 partons from initial state radiation
741 // 
742     for (Int_t i = 2; i <= 7; i++) {
743         Int_t j = 0;
744         // Skip gluons that participate in hard scattering
745         if (i == 4 || i == 5) continue;
746         // Gluons from hard Scattering
747         if (i == 6 || i == 7) {
748             j = i - 6;
749             pxq[j]    = fPyjets->P[0][i];
750             pyq[j]    = fPyjets->P[1][i];
751             pzq[j]    = fPyjets->P[2][i];
752             eq[j]     = fPyjets->P[3][i];
753             mq[j]     = fPyjets->P[4][i];
754         } else {
755             // Gluons from initial state radiation
756             //
757             // Obtain 4-momentum vector from difference between original parton and parton after gluon 
758             // radiation. Energy is calculated independently because initial state radition does not 
759             // conserve strictly momentum and energy for each partonic system independently.
760             //
761             // Not very clean. Should be improved !
762             //
763             //
764             j = i;
765             pxq[j]    = fPyjets->P[0][i] - fPyjets->P[0][i+2];
766             pyq[j]    = fPyjets->P[1][i] - fPyjets->P[1][i+2];
767             pzq[j]    = fPyjets->P[2][i] - fPyjets->P[2][i+2];
768             mq[j]     = fPyjets->P[4][i];
769             eq[j]     = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j] + mq[j] * mq[j]);
770         }
771 //
772 //  Calculate some kinematic variables
773 //
774         yq[j]     = 0.5 * TMath::Log((eq[j] + pzq[j] + 1.e-14) / (eq[j] - pzq[j] + 1.e-14));
775         pq[j]     = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j]);
776         phiq[j]   = TMath::Pi()+TMath::ATan2(-pyq[j], -pxq[j]);
777         ptq[j]    = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j]);
778         thetaq[j] = TMath::ATan2(ptq[j], pzq[j]);
779         qPdg[j]   =  fPyjets->K[1][i];
780     }
781   
782     Double_t int0[4];
783     Double_t int1[4];
784     
785     fGlauber->GetI0I1ForPythiaAndXY(4, phiq, int0, int1, fXJet, fYJet, 15.);
786
787     for (Int_t j = 0; j < 4; j++) {
788         //
789         // Quench only central jets and with E > 10.
790         //
791
792
793         Int_t itype = (qPdg[j] == 21) ? 2 : 1;
794         Double_t eloss = fQuenchingWeights->GetELossRandomKFast(itype, int0[j], int1[j], eq[j]);
795
796         if (TMath::Abs(yq[j]) > 2.5 || eq[j] < 10.) {
797             fZQuench[j] = 0.;
798         } else {
799             if (eq[j] > 40. && TMath::Abs(yq[j]) < 0.5) {
800                 icall ++;
801                 eMean += eloss;
802             }
803             //
804             // Extra pt
805             Double_t l =   fQuenchingWeights->CalcLk(int0[j], int1[j]);     
806             wjtKick[j] = TMath::Sqrt(l *  fQuenchingWeights->CalcQk(int0[j], int1[j]));
807             //
808             // Fractional energy loss
809             fZQuench[j] = eloss / eq[j];
810             //
811             // Avoid complete loss
812             //
813             if (fZQuench[j] > fZmax) fZQuench[j] = fZmax;
814             //
815             // Some debug printing
816
817             
818 //          printf("Initial parton # %3d, Type %3d Energy %10.3f Phi %10.3f Length %10.3f Loss %10.3f Kick %10.3f Mean: %10.3f %10.3f\n", 
819 //                 j, itype, eq[j], phiq[j], l, eloss, wjtKick[j], eMean / Float_t(icall+1), yq[j]);
820             
821 //          fZQuench[j] = 0.8;
822 //          while (fZQuench[j] >= 0.95)  fZQuench[j] = gRandom->Exp(0.2);
823         }
824         
825         quenched[j] = (fZQuench[j] > 0.01);
826     } // primary partons
827     
828     
829
830     Double_t pNew[1000][4];
831     Int_t    kNew[1000];
832     Int_t icount = 0;
833     Double_t zquench[4];
834     
835 //
836 //  System Loop    
837     for (Int_t isys = 0; isys < 4; isys++) {
838 //      Skip to next system if not quenched.
839         if (!quenched[isys]) continue;
840         
841         nGluon[isys]   = 1 + Int_t(fZQuench[isys] / (1. - fZQuench[isys]));
842         if (nGluon[isys] > fNGmax) nGluon[isys] = fNGmax;
843         zquench[isys] = 1. - TMath::Power(1. - fZQuench[isys], 1./Double_t(nGluon[isys]));
844         wjtKick[isys]  = wjtKick[isys] / TMath::Sqrt(Double_t(nGluon[isys]));
845
846
847         
848         Int_t igMin = -1;
849         Int_t igMax = -1;
850         Double_t pg[4] = {0., 0., 0., 0.};
851         
852 //
853 // Loop on radiation events
854
855         for (Int_t iglu = 0; iglu < nGluon[isys]; iglu++) {
856             while (1) {
857                 icount = 0;
858                 for (Int_t k = 0; k < 4; k++)
859                 {
860                     p0[isys][k] = 0.;
861                     p1[isys][k] = 0.;
862                     p2[isys][k] = 0.;
863                 }
864 //      Loop over partons
865                 for (Int_t i = 0; i < numpart; i++)
866                 {
867                     imo =  fPyjets->K[2][i];
868                     kst =  fPyjets->K[0][i];
869                     pdg =  fPyjets->K[1][i];
870                     
871                 
872                 
873 //      Quarks and gluons only
874                     if (pdg != 21 && TMath::Abs(pdg) > 6) continue;
875 //      Particles from hard scattering only
876                     
877                     if (imo > 8 && imo < 1000) imo = fPyjets->K[2][imo - 1];
878                     Int_t imom = imo % 1000;
879                     if ((isys == 0 || isys == 1) && ((imom != (isys + 7)))) continue;
880                     if ((isys == 2 || isys == 3) && ((imom != (isys + 1)))) continue;               
881                     
882                     
883 //      Skip comment lines
884                     if (kst != 1 && kst != 2) continue;
885 //
886 //      Parton kinematic
887                     px    = fPyjets->P[0][i];
888                     py    = fPyjets->P[1][i];
889                     pz    = fPyjets->P[2][i];
890                     e     = fPyjets->P[3][i];
891                     m     = fPyjets->P[4][i];
892                     pt    = TMath::Sqrt(px * px + py * py);
893                     p     = TMath::Sqrt(px * px + py * py + pz * pz); 
894                     phi   = TMath::Pi() + TMath::ATan2(-py, -px);
895                     theta = TMath::ATan2(pt, pz);
896                 
897 //
898 //      Save 4-momentum sum for balancing
899                     Int_t index = isys;
900                     
901                     p0[index][0] += px;
902                     p0[index][1] += py;
903                     p0[index][2] += pz;
904                     p0[index][3] += e;
905                 
906                     klast[index] = i;
907                     
908 //
909 //      Fractional energy loss
910                     Double_t z = zquench[index];
911                     
912                     
913 //      Don't fully quench radiated gluons
914 //
915                     if (imo > 1000) {
916 //      This small factor makes sure that the gluons are not too close in phase space to avoid recombination
917 //
918
919                         z = 0.02;
920                     }
921 //                  printf("z: %d %f\n", imo, z);
922                     
923
924 //
925                     
926                     //
927                     //
928                     //      Transform into frame in which initial parton is along z-axis
929                     //
930                     TVector3 v(px, py, pz);
931                     v.RotateZ(-phiq[index]);  v.RotateY(-thetaq[index]);
932                     Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pl  = v.Z();
933
934                     Double_t jt  = TMath::Sqrt(pxs * pxs + pys * pys);
935                     Double_t mt2 = jt * jt + m * m;
936                     Double_t zmax = 1.;     
937                     //
938                     // Kinematic limit on z
939                     //
940                     if (m > 0.) zmax = 1. - m / TMath::Sqrt(m * m + jt * jt);
941                     //
942                     // Change light-cone kinematics rel. to initial parton
943                     //  
944                     Double_t eppzOld = e + pl;
945                     Double_t empzOld = e - pl;
946                     
947                     Double_t eppzNew = (1. - z) * eppzOld;
948                     Double_t empzNew = empzOld - mt2 * z / eppzOld;
949                     Double_t eNew    = 0.5 * (eppzNew + empzNew);
950                     Double_t plNew   = 0.5 * (eppzNew - empzNew);
951                     
952                     Double_t jtNew;
953                     //
954                     // if mt very small (or sometimes even < 0 for numerical reasons) set it to 0
955                     Double_t mt2New = eppzNew * empzNew;
956                     if (mt2New < 1.e-8) mt2New = 0.;
957                     if (z < zmax) {
958                         if (m * m > mt2New) {
959                             //
960                             // This should not happen 
961                             //
962                             Fatal("Quench()", "This should never happen %e %e %e!", m, eppzNew, empzNew);
963                             jtNew = 0;
964                         } else {
965                             jtNew    = TMath::Sqrt(mt2New - m * m);
966                         }
967                     } else {
968                         // If pT is to small (probably a leading massive particle) we scale only the energy
969                         // This can cause negative masses of the radiated gluon
970                         // Let's hope for the best ...
971                         jtNew = jt;
972                         eNew  = TMath::Sqrt(plNew * plNew + mt2);
973                         
974                     }
975                     //
976                     //     Calculate new px, py
977                     //
978                     Double_t pxNew = 0;
979                     Double_t pyNew = 0;
980
981                     if (jt > 0.) {
982                         pxNew   = jtNew / jt * pxs;
983                         pyNew   = jtNew / jt * pys;     
984                     }
985                     
986 //                  Double_t dpx = pxs - pxNew;
987 //                  Double_t dpy = pys - pyNew;
988 //                  Double_t dpz = pl  - plNew;
989 //                  Double_t de  = e   - eNew;
990 //                  Double_t dmass2 = de * de  - dpx * dpx - dpy * dpy - dpz * dpz;
991 //                  printf("New mass (1) %e %e %e %e %e %e %e \n", dmass2, jt, jtNew, pl, plNew, e, eNew);
992 //                  printf("New mass (2) %e %e \n", pxNew, pyNew);
993                     //
994                     //      Rotate back
995                     //  
996                     TVector3 w(pxNew, pyNew, plNew);
997                     w.RotateY(thetaq[index]); w.RotateZ(phiq[index]);
998                     pxNew = w.X(); pyNew = w.Y(); plNew = w.Z();
999                 
1000                     p1[index][0] += pxNew;
1001                     p1[index][1] += pyNew;
1002                     p1[index][2] += plNew;
1003                     p1[index][3] += eNew;       
1004                     //
1005                     // Updated 4-momentum vectors
1006                     //
1007                     pNew[icount][0]  = pxNew;
1008                     pNew[icount][1]  = pyNew;
1009                     pNew[icount][2]  = plNew;
1010                     pNew[icount][3]  = eNew;
1011                     kNew[icount]     = i;
1012                     icount++;
1013                 } // parton loop
1014                 //
1015                 // Check if there was phase-space for quenching
1016                 //
1017
1018                 if (icount == 0) quenched[isys] = kFALSE;
1019                 if (!quenched[isys]) break;
1020                 
1021                 for (Int_t j = 0; j < 4; j++) 
1022                 {
1023                     p2[isys][j] = p0[isys][j] - p1[isys][j];
1024                 }
1025                 p2[isys][4] = p2[isys][3] * p2[isys][3] - p2[isys][0] * p2[isys][0] - p2[isys][1] * p2[isys][1] - p2[isys][2] * p2[isys][2];
1026                 if (p2[isys][4] > 0.) {
1027                     p2[isys][4] = TMath::Sqrt(p2[isys][4]);
1028                     break;
1029                 } else {
1030                     printf("Warning negative mass squared in system %d %f ! \n", isys, zquench[isys]);
1031                     printf("4-Momentum: %10.3e %10.3e %10.3e %10.3e %10.3e \n", p2[isys][0], p2[isys][1], p2[isys][2], p2[isys][3], p2[isys][4]);
1032                     if (p2[isys][4] < -0.01) {
1033                         printf("Negative mass squared !\n");
1034                         // Here we have to put the gluon back to mass shell
1035                         // This will lead to a small energy imbalance
1036                         p2[isys][4]  = 0.;
1037                         p2[isys][3]  = TMath::Sqrt(p2[isys][0] * p2[isys][0] + p2[isys][1] * p2[isys][1] + p2[isys][2] * p2[isys][2]);
1038                         break;
1039                     } else {
1040                         p2[isys][4] = 0.;
1041                         break;
1042                     }
1043                 }
1044                 /*
1045                 zHeavy *= 0.98;
1046                 printf("zHeavy lowered to %f\n", zHeavy);
1047                 if (zHeavy < 0.01) {
1048                     printf("No success ! \n");
1049                     icount = 0;
1050                     quenched[isys] = kFALSE;
1051                     break;
1052                 }
1053                 */
1054             } // iteration on z (while)
1055             
1056 //          Update  event record
1057             for (Int_t k = 0; k < icount; k++) {
1058 //              printf("%6d %6d %10.3e %10.3e %10.3e %10.3e\n", k, kNew[k], pNew[k][0],pNew[k][1], pNew[k][2], pNew[k][3] );
1059                 fPyjets->P[0][kNew[k]] = pNew[k][0];
1060                 fPyjets->P[1][kNew[k]] = pNew[k][1];
1061                 fPyjets->P[2][kNew[k]] = pNew[k][2];
1062                 fPyjets->P[3][kNew[k]] = pNew[k][3];
1063             }
1064             //
1065             // Add the gluons
1066             //
1067             Int_t ish = 0;    
1068             Int_t iGlu;
1069             if (!quenched[isys]) continue;
1070 //
1071 //      Last parton from shower i
1072             Int_t in = klast[isys];
1073 //
1074 //      Continue if no parton in shower i selected
1075             if (in == -1) continue;
1076 //  
1077 //      If this is the second initial parton and it is behind the first move pointer by previous ish
1078             if (isys == 1 && klast[1] > klast[0]) in += ish;
1079 //
1080 //      Starting index
1081             
1082 //          jmin = in - 1;
1083 // How many additional gluons will be generated
1084             ish  = 1;
1085             if (p2[isys][4] > 0.05) ish = 2;
1086 //
1087 //      Position of gluons
1088             iGlu = numpart;
1089             if (iglu == 0) igMin = iGlu;
1090             igMax = iGlu;
1091             numpart += ish;
1092             (fPyjets->N) += ish;
1093             
1094             if (ish == 1) {
1095                 fPyjets->P[0][iGlu] = p2[isys][0];
1096                 fPyjets->P[1][iGlu] = p2[isys][1];
1097                 fPyjets->P[2][iGlu] = p2[isys][2];
1098                 fPyjets->P[3][iGlu] = p2[isys][3];
1099                 fPyjets->P[4][iGlu] = p2[isys][4];
1100                 
1101                 fPyjets->K[0][iGlu] = 1;
1102                 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu] = 1;
1103                 fPyjets->K[1][iGlu] = 21;       
1104                 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1105                 fPyjets->K[3][iGlu] = -1;       
1106                 fPyjets->K[4][iGlu] = -1;
1107                 
1108                 pg[0] += p2[isys][0];
1109                 pg[1] += p2[isys][1];
1110                 pg[2] += p2[isys][2];
1111                 pg[3] += p2[isys][3];
1112             } else {
1113                 //
1114                 // Split gluon in rest frame.
1115                 //
1116                 Double_t bx   =  p2[isys][0] / p2[isys][3];
1117                 Double_t by   =  p2[isys][1] / p2[isys][3];
1118                 Double_t bz   =  p2[isys][2] / p2[isys][3];
1119                 Double_t pst  =  p2[isys][4] / 2.;
1120                 //
1121                 // Isotropic decay ????
1122                 Double_t cost = 2. * gRandom->Rndm() - 1.;
1123                 Double_t sint = TMath::Sqrt(1. - cost * cost);
1124                 Double_t phis =  2. * TMath::Pi() * gRandom->Rndm();
1125                 
1126                 Double_t pz1 =   pst * cost;
1127                 Double_t pz2 =  -pst * cost;
1128                 Double_t pt1 =   pst * sint;
1129                 Double_t pt2 =  -pst * sint;
1130                 Double_t px1 =   pt1 * TMath::Cos(phis);
1131                 Double_t py1 =   pt1 * TMath::Sin(phis);            
1132                 Double_t px2 =   pt2 * TMath::Cos(phis);
1133                 Double_t py2 =   pt2 * TMath::Sin(phis);            
1134                 
1135                 fPyjets->P[0][iGlu] = px1;
1136                 fPyjets->P[1][iGlu] = py1;
1137                 fPyjets->P[2][iGlu] = pz1;
1138                 fPyjets->P[3][iGlu] = pst;
1139                 fPyjets->P[4][iGlu] = 0.;
1140                 
1141                 fPyjets->K[0][iGlu] = 1 ;
1142                 fPyjets->K[1][iGlu] = 21;       
1143                 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1144                 fPyjets->K[3][iGlu] = -1;       
1145                 fPyjets->K[4][iGlu] = -1;
1146                 
1147                 fPyjets->P[0][iGlu+1] = px2;
1148                 fPyjets->P[1][iGlu+1] = py2;
1149                 fPyjets->P[2][iGlu+1] = pz2;
1150                 fPyjets->P[3][iGlu+1] = pst;
1151                 fPyjets->P[4][iGlu+1] = 0.;
1152                 
1153                 fPyjets->K[0][iGlu+1] = 1;
1154                 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu+1] = 1;
1155                 fPyjets->K[1][iGlu+1] = 21;     
1156                 fPyjets->K[2][iGlu+1] = fPyjets->K[2][in] + 1000;
1157                 fPyjets->K[3][iGlu+1] = -1;     
1158                 fPyjets->K[4][iGlu+1] = -1;
1159                 SetMSTU(1,0);
1160                 SetMSTU(2,0);
1161                 //
1162                 // Boost back
1163                 //
1164                 Pyrobo(iGlu + 1, iGlu + 2, 0., 0., bx, by, bz);
1165             }
1166 /*    
1167             for (Int_t ig = iGlu; ig < iGlu+ish; ig++) {
1168                 Double_t px, py, pz;
1169                 px = fPyjets->P[0][ig]; 
1170                 py = fPyjets->P[1][ig]; 
1171                 pz = fPyjets->P[2][ig]; 
1172                 TVector3 v(px, py, pz);
1173                 v.RotateZ(-phiq[isys]);
1174                 v.RotateY(-thetaq[isys]);
1175                 Double_t pxs     = v.X(); Double_t pys = v.Y(); Double_t pzs  = v.Z();     
1176                 Double_t r       = AliPythiaRndm::GetPythiaRandom()->Rndm();
1177                 Double_t jtKick  = 0.3 * TMath::Sqrt(-TMath::Log(r));
1178                 if (ish == 2)   jtKick  = wjtKick[i] * TMath::Sqrt(-TMath::Log(r)) / TMath::Sqrt(2.);
1179                 Double_t phiKick = 2. * TMath::Pi() * AliPythiaRndm::GetPythiaRandom()->Rndm();
1180                 pxs += jtKick * TMath::Cos(phiKick);
1181                 pys += jtKick * TMath::Sin(phiKick);
1182                 TVector3 w(pxs, pys, pzs);
1183                 w.RotateY(thetaq[isys]);
1184                 w.RotateZ(phiq[isys]);
1185                 fPyjets->P[0][ig] = w.X(); 
1186                 fPyjets->P[1][ig] = w.Y(); 
1187                 fPyjets->P[2][ig] = w.Z(); 
1188                 fPyjets->P[2][ig] = w.Mag();
1189             }
1190 */
1191         } // kGluon         
1192         
1193         
1194     // Check energy conservation
1195         Double_t pxs = 0.;
1196         Double_t pys = 0.;
1197         Double_t pzs = 0.;      
1198         Double_t es  = 14000.;
1199         
1200         for (Int_t i = 0; i < numpart; i++)
1201         {
1202             kst =  fPyjets->K[0][i];
1203             if (kst != 1 && kst != 2) continue;
1204             pxs += fPyjets->P[0][i];
1205             pys += fPyjets->P[1][i];
1206             pzs += fPyjets->P[2][i];        
1207             es  -= fPyjets->P[3][i];        
1208         }
1209         if (TMath::Abs(pxs) > 1.e-2 ||
1210             TMath::Abs(pys) > 1.e-2 ||
1211             TMath::Abs(pzs) > 1.e-1) {
1212             printf("%e %e %e %e\n", pxs, pys, pzs, es);
1213 //              Fatal("Quench()", "4-Momentum non-conservation");
1214         }
1215         
1216     } // end quenching loop (systems)
1217 // Clean-up
1218     for (Int_t i = 0; i < numpart; i++)
1219     {
1220         imo =  fPyjets->K[2][i];
1221         if (imo > 1000) {
1222             fPyjets->K[2][i] = fPyjets->K[2][i] % 1000;
1223         }
1224     }
1225 //      this->Pylist(1);
1226 } // end quench
1227
1228
1229 void AliPythia6::Pyquen(Double_t a, Int_t ibf, Double_t b)
1230 {
1231     // Igor Lokthine's quenching routine
1232     // http://lokhtin.web.cern.ch/lokhtin/pyquen/pyquen.txt
1233
1234     pyquen(a, ibf, b);
1235 }
1236
1237 void AliPythia6::SetPyquenParameters(Double_t t0, Double_t tau0, Int_t nf, Int_t iengl, Int_t iangl)
1238 {
1239     // Set the parameters for the PYQUEN package.
1240     // See comments in PyquenCommon.h
1241     
1242     
1243     PYQPAR.t0    = t0;
1244     PYQPAR.tau0  = tau0;
1245     PYQPAR.nf    = nf;
1246     PYQPAR.iengl = iengl;
1247     PYQPAR.iangl = iangl;
1248 }
1249
1250 void  AliPythia6::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
1251 {
1252 //
1253 // Load event into Pythia Common Block
1254 //
1255
1256     Int_t npart = stack -> GetNprimary();
1257     Int_t n0 = 0;
1258     
1259     if (!flag) {
1260         GetPyjets()->N = npart;
1261     } else {
1262         n0 = GetPyjets()->N;
1263         GetPyjets()->N = n0 + npart;
1264     }
1265     
1266     
1267     for (Int_t part = 0; part < npart; part++) {
1268         TParticle *mPart = stack->Particle(part);
1269         
1270         Int_t kf     =  mPart->GetPdgCode();
1271         Int_t ks     =  mPart->GetStatusCode();
1272         Int_t idf    =  mPart->GetFirstDaughter();
1273         Int_t idl    =  mPart->GetLastDaughter();
1274         
1275         if (reHadr) {
1276             if (ks == 11 || ks == 12) {
1277                 ks  -= 10;
1278                 idf  = -1;
1279                 idl  = -1;
1280             }
1281         }
1282         
1283         Float_t px = mPart->Px();
1284         Float_t py = mPart->Py();
1285         Float_t pz = mPart->Pz();
1286         Float_t e  = mPart->Energy();
1287         Float_t m  = mPart->GetCalcMass();
1288         
1289         
1290         (GetPyjets())->P[0][part+n0] = px;
1291         (GetPyjets())->P[1][part+n0] = py;
1292         (GetPyjets())->P[2][part+n0] = pz;
1293         (GetPyjets())->P[3][part+n0] = e;
1294         (GetPyjets())->P[4][part+n0] = m;
1295         
1296         (GetPyjets())->K[1][part+n0] = kf;
1297         (GetPyjets())->K[0][part+n0] = ks;
1298         (GetPyjets())->K[3][part+n0] = idf + 1;
1299         (GetPyjets())->K[4][part+n0] = idl + 1;
1300         (GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1;
1301     }
1302 }
1303
1304
1305 void AliPythia6::Pyevnw()
1306 {
1307     // New multiple interaction scenario
1308     pyevnw();
1309 }
1310
1311 void AliPythia6::GetQuenchingParameters(Double_t& xp, Double_t& yp, Double_t z[4])
1312 {
1313     // Return event specific quenching parameters
1314     xp = fXJet;
1315     yp = fYJet;
1316     for (Int_t i = 0; i < 4; i++) z[i] = fZQuench[i];
1317
1318 }
1319
1320 void AliPythia6::ConfigHeavyFlavor()
1321 {
1322     //
1323     // Default configuration for Heavy Flavor production
1324     //
1325     // All QCD processes
1326     //
1327     SetMSEL(1);
1328     
1329     // No multiple interactions
1330     SetMSTP(81,0);
1331     SetPARP(81, 0.);
1332     SetPARP(82, 0.);    
1333     // Initial/final parton shower on (Pythia default)
1334     SetMSTP(61,1);
1335     SetMSTP(71,1);
1336     
1337     // 2nd order alpha_s
1338     SetMSTP(2,2);
1339     
1340     // QCD scales
1341     SetMSTP(32,2);
1342     SetPARP(34,1.0);
1343 }
1344
1345 void AliPythia6::AtlasTuning()
1346 {
1347     //
1348     // Configuration for the ATLAS tuning
1349         SetMSTP(51,AliStructFuncType::PDFsetIndex(kCTEQ5L));
1350         SetMSTP(81,1);             // Multiple Interactions ON
1351         SetMSTP(82,4);             // Double Gaussian Model
1352         SetPARP(81,1.9);           // Min. pt for multiple interactions (default in 6.2-14) 
1353         SetPARP(82,1.8);           // [GeV]    PT_min at Ref. energy
1354         SetPARP(89,1000.);         // [GeV]   Ref. energy
1355         SetPARP(90,0.16);          // 2*epsilon (exponent in power law)
1356         SetPARP(83,0.5);           // Core density in proton matter distribution (def.value)
1357         SetPARP(84,0.5);           // Core radius
1358         SetPARP(85,0.33);          // Regulates gluon prod. mechanism
1359         SetPARP(86,0.66);          // Regulates gluon prod. mechanism
1360         SetPARP(67,1);             // Regulates Initial State Radiation
1361 }
1362
1363 void AliPythia6::SetPtHardRange(Float_t ptmin, Float_t ptmax)
1364 {
1365     // Set the pt hard range
1366     SetCKIN(3, ptmin);
1367     SetCKIN(4, ptmax);
1368 }
1369
1370 void AliPythia6::SetYHardRange(Float_t ymin, Float_t ymax)
1371 {
1372     // Set the y hard range
1373     SetCKIN(7, ymin);
1374     SetCKIN(8, ymax);
1375 }
1376
1377
1378 void AliPythia6::SetFragmentation(Int_t flag)
1379 {
1380     // Switch fragmentation on/off
1381     SetMSTP(111, flag);
1382 }
1383
1384 void AliPythia6::SetInitialAndFinalStateRadiation(Int_t flag1, Int_t flag2)
1385 {
1386 //  initial state radiation    
1387     SetMSTP(61, flag1);
1388 //  final state radiation
1389     SetMSTP(71, flag2);
1390 }
1391
1392 void AliPythia6::SetIntrinsicKt(Float_t kt)
1393 {
1394     // Set the inreinsic kt
1395     if (kt > 0.) {
1396         SetMSTP(91,1);
1397         SetPARP(91,kt); 
1398         SetPARP(93, 4. * kt);
1399     } else {
1400         SetMSTP(91,0);
1401     }
1402 }
1403
1404 void AliPythia6::SwitchHFOff()
1405 {
1406     // Switch off heavy flavor
1407     // Maximum number of quark flavours used in pdf 
1408     SetMSTP(58, 3);
1409     // Maximum number of flavors that can be used in showers
1410     SetMSTJ(45, 3);     
1411 }
1412
1413 void AliPythia6::SetPycellParameters(Float_t etamax, Int_t neta, Int_t nphi,
1414                                        Float_t thresh, Float_t etseed, Float_t minet, Float_t r)
1415 {
1416 // Set pycell parameters
1417     SetPARU(51, etamax);
1418     SetMSTU(51, neta);
1419     SetMSTU(52, nphi);
1420     SetPARU(58, thresh);
1421     SetPARU(52, etseed);
1422     SetPARU(53, minet);
1423     SetPARU(54, r);
1424     SetMSTU(54,  2);
1425 }
1426
1427 void AliPythia6::ModifiedSplitting()
1428 {
1429     // Modified splitting probability as a model for quenching
1430     SetPARJ(200, 0.8);
1431     SetMSTJ(41, 1);  // QCD radiation only
1432     SetMSTJ(42, 2);  // angular ordering
1433     SetMSTJ(44, 2);  // option to run alpha_s
1434     SetMSTJ(47, 0);  // No correction back to hard scattering element
1435     SetMSTJ(50, 0);  // No coherence in first branching
1436     SetPARJ(82, 1.); // Cut off for parton showers
1437 }
1438
1439 void AliPythia6::SwitchHadronisationOff()
1440 {
1441     // Switch off hadronisarion
1442     SetMSTJ(1, 0);
1443 }
1444
1445 void AliPythia6::SwitchHadronisationOn()
1446 {
1447     // Switch on hadronisarion
1448     SetMSTJ(1, 1);
1449 }
1450
1451
1452 void AliPythia6::GetXandQ(Float_t& x1, Float_t& x2, Float_t& q)
1453 {
1454     // Get x1, x2 and Q for this event
1455     
1456     q  = GetVINT(51);
1457     x1 = GetVINT(41);
1458     x2 = GetVINT(42);
1459 }
1460
1461 Float_t AliPythia6::GetXSection()
1462 {
1463     // Get the total cross-section
1464     return (GetPARI(1));
1465 }
1466
1467 Float_t AliPythia6::GetPtHard()
1468 {
1469     // Get the pT hard for this event
1470     return GetVINT(47);
1471 }
1472
1473 Int_t AliPythia6::ProcessCode()
1474 {
1475     // Get the subprocess code
1476     return GetMSTI(1);
1477 }
1478
1479 void AliPythia6::PrintStatistics()
1480 {
1481     // End of run statistics
1482     Pystat(1);
1483 }
1484
1485 void AliPythia6::EventListing()
1486 {
1487     // End of run statistics
1488     Pylist(2);
1489 }
1490
1491 AliPythia6& AliPythia6::operator=(const  AliPythia6& rhs)
1492 {
1493 // Assignment operator
1494     rhs.Copy(*this);
1495     return *this;
1496 }
1497
1498  void AliPythia6::Copy(TObject&) const
1499 {
1500     //
1501     // Copy 
1502     //
1503     Fatal("Copy","Not implemented!\n");
1504 }