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