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