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