]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PYTHIA6/AliPythia.cxx
Option for patching of the omega dalitz decay.
[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       // Tuning of Pythia parameters aimed to get a resonable agreement
421       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
422       // c-cbar single inclusive and double differential distributions.
423       // This parameter settings are meant to work with pp collisions
424       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
425       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
426       // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
427         ConfigHeavyFlavor();
428       // Intrinsic <kT^2>
429         SetMSTP(91,1);
430         SetPARP(91,1.);
431         SetPARP(93,5.);
432         
433       // Set c-quark mass
434         SetPMAS(4,1,1.2);
435       break;
436     case kPyCharmppMNRwmi:
437       // Tuning of Pythia parameters aimed to get a resonable agreement
438       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
439       // c-cbar single inclusive and double differential distributions.
440       // This parameter settings are meant to work with pp collisions
441       // and with kCTEQ5L PDFs.
442       // Added multiple interactions according to ATLAS tune settings.
443       // To get a "reasonable" agreement with MNR results, events have to be 
444       // generated with the minimum ptHard (AliGenPythia::SetPtHard)
445       // set to 2.76 GeV.
446       // To get a "perfect" agreement with MNR results, events have to be 
447       // generated in four ptHard bins with the following relative 
448       // normalizations:
449       // 2.76-3 GeV: 25%
450       //    3-4 GeV: 40%
451       //    4-8 GeV: 29%
452       //     >8 GeV:  6%
453         ConfigHeavyFlavor();
454       // Intrinsic <kT^2>
455         SetMSTP(91,1);
456         SetPARP(91,1.);
457         SetPARP(93,5.);
458
459       // Set c-quark mass
460         SetPMAS(4,1,1.2);
461         AtlasTuning();
462         break;
463     case kPyBeautyPbPbMNR:
464       // Tuning of Pythia parameters aimed to get a resonable agreement
465       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
466       // b-bbar single inclusive and double differential distributions.
467       // This parameter settings are meant to work with Pb-Pb collisions
468       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
469       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
470       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
471         ConfigHeavyFlavor();
472       // QCD scales
473         SetPARP(67,1.0);
474         SetPARP(71,1.0);
475       // Intrinsic <kT>
476         SetMSTP(91,1);
477         SetPARP(91,2.035);
478         SetPARP(93,10.17);
479       // Set b-quark mass
480         SetPMAS(5,1,4.75);
481       break;
482     case kPyBeautypPbMNR:
483       // Tuning of Pythia parameters aimed to get a resonable agreement
484       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
485       // b-bbar single inclusive and double differential distributions.
486       // This parameter settings are meant to work with p-Pb collisions
487       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
488       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
489       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
490         ConfigHeavyFlavor();
491       // QCD scales
492         SetPARP(67,1.0);
493         SetPARP(71,1.0);
494       // Intrinsic <kT>
495         SetMSTP(91,1);
496         SetPARP(91,1.60);
497         SetPARP(93,8.00);
498       // Set b-quark mass
499         SetPMAS(5,1,4.75);
500       break;
501     case kPyBeautyppMNR:
502       // Tuning of Pythia parameters aimed to get a resonable agreement
503       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
504       // b-bbar single inclusive and double differential distributions.
505       // This parameter settings are meant to work with pp collisions
506       // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
507       // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
508       // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
509         ConfigHeavyFlavor();
510       // QCD scales
511         SetPARP(67,1.0);
512         SetPARP(71,1.0);
513         
514         // Intrinsic <kT>
515         SetMSTP(91,1);
516         SetPARP(91,1.);
517         SetPARP(93,5.);
518         
519         // Set b-quark mass
520         SetPMAS(5,1,4.75);
521       break;
522      case kPyBeautyJets:
523      case kPyBeautyppMNRwmi:
524       // Tuning of Pythia parameters aimed to get a resonable agreement
525       // between with the NLO calculation by Mangano, Nason, Ridolfi for the
526       // b-bbar single inclusive and double differential distributions.
527       // This parameter settings are meant to work with pp collisions
528       // and with kCTEQ5L PDFs.
529       // Added multiple interactions according to ATLAS tune settings.
530       // To get a "reasonable" agreement with MNR results, events have to be 
531       // generated with the minimum ptHard (AliGenPythia::SetPtHard)
532       // set to 2.76 GeV.
533       // To get a "perfect" agreement with MNR results, events have to be 
534       // generated in four ptHard bins with the following relative 
535       // normalizations:
536       // 2.76-4 GeV:  5% 
537       //    4-6 GeV: 31%
538       //    6-8 GeV: 28%
539       //     >8 GeV: 36%
540          ConfigHeavyFlavor();
541       // QCD scales
542          SetPARP(67,1.0);
543          SetPARP(71,1.0);
544          
545          // Intrinsic <kT>
546          SetMSTP(91,1);
547          SetPARP(91,1.);
548          SetPARP(93,5.);
549
550       // Set b-quark mass
551          SetPMAS(5,1,4.75);
552
553          AtlasTuning();
554          break; 
555     case kPyW:
556
557       //Inclusive production of W+/-
558       SetMSEL(0);
559       //f fbar -> W+ 
560       SetMSUB(2,1);
561       //        //f fbar -> g W+
562       //        SetMSUB(16,1);
563       //        //f fbar -> gamma W+
564       //        SetMSUB(20,1);
565       //        //f g -> f W+  
566       //        SetMSUB(31,1);
567       //        //f gamma -> f W+
568       //        SetMSUB(36,1);
569       
570       // Initial/final parton shower on (Pythia default)
571       // With parton showers on we are generating "W 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     case kPyZ:
578
579       //Inclusive production of Z
580       SetMSEL(0);
581       //f fbar -> Z/gamma
582       SetMSUB(1,1);
583       
584       //       // f fbar -> g Z/gamma
585       //       SetMSUB(15,1);
586       //       // f fbar -> gamma Z/gamma
587       //       SetMSUB(19,1);
588       //       // f g -> f Z/gamma
589       //       SetMSUB(30,1);
590       //       // f gamma -> f Z/gamma
591       //       SetMSUB(35,1);
592       
593       //only Z included, not gamma
594       SetMSTP(43,2);
595       
596       // Initial/final parton shower on (Pythia default)
597       // With parton showers on we are generating "Z inclusive process"
598       SetMSTP(61,1); //Initial QCD & QED showers on
599       SetMSTP(71,1); //Final QCD & QED showers on
600       
601       break;  
602
603     }
604 //
605 //  Initialize PYTHIA
606 //
607 //  Select the tune
608     if (itune > -1) Pytune(itune);
609     
610 //  
611     SetMSTP(41,1);   // all resonance decays switched on
612     Initialize("CMS","p","p",fEcms);
613     fOmegaDalitz.Init();
614 }
615
616 Int_t AliPythia::CheckedLuComp(Int_t kf)
617 {
618 // Check Lund particle code (for debugging)
619     Int_t kc=Pycomp(kf);
620     printf("\n Lucomp kf,kc %d %d",kf,kc);
621     return kc;
622 }
623
624 void AliPythia::SetNuclei(Int_t a1, Int_t a2, Int_t pdf)
625 {
626 // Treat protons as inside nuclei with mass numbers a1 and a2  
627 //    The MSTP array in the PYPARS common block is used to enable and 
628 //    select the nuclear structure functions. 
629 //    MSTP(52)  : (D=1) choice of proton and nuclear structure-function library
630 //            =1: internal PYTHIA acording to MSTP(51) 
631 //            =2: PDFLIB proton  s.f., with MSTP(51)  = 1000xNGROUP+NSET
632 //    If the following mass number both not equal zero, nuclear corrections of the stf are used.
633 //    MSTP(192) : Mass number of nucleus side 1
634 //    MSTP(193) : Mass number of nucleus side 2
635 //    MSTP(194) : Nuclear structure function: 0: EKS98 1:EPS08
636     SetMSTP(52,2);
637     SetMSTP(192, a1);
638     SetMSTP(193, a2); 
639     SetMSTP(194, pdf);
640 }
641         
642
643 AliPythia* AliPythia::Instance()
644
645 // Set random number generator 
646     if (fgAliPythia) {
647         return fgAliPythia;
648     } else {
649         fgAliPythia = new AliPythia();
650         return fgAliPythia;
651     }
652 }
653
654 void AliPythia::PrintParticles()
655
656 // Print list of particl properties
657     Int_t np = 0;
658     char*   name = new char[16];    
659     for (Int_t kf=0; kf<1000000; kf++) {
660         for (Int_t c = 1;  c > -2; c-=2) {
661             Int_t kc = Pycomp(c*kf);
662             if (kc) {
663                 Float_t mass  = GetPMAS(kc,1);
664                 Float_t width = GetPMAS(kc,2);  
665                 Float_t tau   = GetPMAS(kc,4);
666
667                 Pyname(kf,name);
668         
669                 np++;
670                 
671                 printf("\n mass, width, tau: %6d %s %10.3f %10.3e %10.3e", 
672                        c*kf, name, mass, width, tau);
673             }
674         }
675     }
676     printf("\n Number of particles %d \n \n", np);
677 }
678
679 void  AliPythia::ResetDecayTable()
680 {
681 //  Set default values for pythia decay switches
682     Int_t i;
683     for (i = 1; i <  501; i++) SetMDCY(i,1,fDefMDCY[i]);
684     for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
685 }
686
687 void  AliPythia::SetDecayTable()
688 {
689 //  Set default values for pythia decay switches
690 //
691     Int_t i;
692     for (i = 1; i <  501; i++) fDefMDCY[i] = GetMDCY(i,1);
693     for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
694 }
695
696 void  AliPythia::Pyclus(Int_t& njet)
697 {
698 //  Call Pythia clustering algorithm
699 //
700     pyclus(njet);
701 }
702
703 void  AliPythia::Pycell(Int_t& njet)
704 {
705 //  Call Pythia jet reconstruction algorithm
706 //
707     pycell(njet);
708 }
709
710 void  AliPythia::Pyshow(Int_t ip1, Int_t ip2, Double_t qmax)
711 {
712 //  Call Pythia jet reconstruction algorithm
713 //
714     pyshow(ip1, ip2, qmax);
715 }
716
717 void AliPythia::Pyrobo(Int_t imi, Int_t ima, Double_t the, Double_t phi, Double_t bex, Double_t bey, Double_t bez)
718 {
719     pyrobo(imi, ima, the, phi, bex, bey, bez);
720 }
721
722 void AliPythia::Pytune(Int_t itune)
723 {
724 /*
725 C
726 C ITUNE    NAME (detailed descriptions below)
727 C     0 Default : No settings changed => linked Pythia version's defaults.
728 C ====== Old UE, Q2-ordered showers ==========================================
729 C   100       A : Rick Field's CDF Tune A 
730 C   101      AW : Rick Field's CDF Tune AW
731 C   102      BW : Rick Field's CDF Tune BW
732 C   103      DW : Rick Field's CDF Tune DW
733 C   104     DWT : Rick Field's CDF Tune DW with slower UE energy scaling
734 C   105      QW : Rick Field's CDF Tune QW (NB: needs CTEQ6.1M pdfs externally)
735 C   106 ATLAS-DC2: Arthur Moraes' (old) ATLAS tune (ATLAS DC2 / Rome)
736 C   107     ACR : Tune A modified with annealing CR
737 C   108      D6 : Rick Field's CDF Tune D6 (NB: needs CTEQ6L pdfs externally)
738 C   109     D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
739 C ====== Intermediate Models =================================================
740 C   200    IM 1 : Intermediate model: new UE, Q2-ordered showers, annealing CR
741 C   201     APT : Tune A modified to use pT-ordered final-state showers
742 C ====== New UE, interleaved pT-ordered showers, annealing CR ================
743 C   300      S0 : Sandhoff-Skands Tune 0 
744 C   301      S1 : Sandhoff-Skands Tune 1
745 C   302      S2 : Sandhoff-Skands Tune 2
746 C   303     S0A : S0 with "Tune A" UE energy scaling
747 C   304    NOCR : New UE "best try" without colour reconnections
748 C   305     Old : New UE, original (primitive) colour reconnections
749 C   306 ATLAS-CSC: Arthur Moraes' (new) ATLAS tune (needs CTEQ6L externally)
750 C ======= The Uppsala models =================================================
751 C   ( NB! must be run with special modified Pythia 6.215 version )
752 C   ( available from http://www.isv.uu.se/thep/MC/scigal/        )
753 C   400   GAL 0 : Generalized area-law model. Old parameters
754 C   401   SCI 0 : Soft-Colour-Interaction model. Old parameters
755 C   402   GAL 1 : Generalized area-law model. Tevatron MB retuned (Skands)
756 */
757     pytune(itune);
758 }
759
760 void AliPythia::Py2ent(Int_t idx, Int_t pdg1, Int_t pdg2, Double_t p){
761   // Inset 2-parton system at line idx
762   py2ent(idx, pdg1, pdg2, p);
763 }
764
765
766 void AliPythia::InitQuenching(Float_t cMin, Float_t cMax, Float_t k, Int_t iECMethod, Float_t zmax, Int_t ngmax)
767 {
768 // Initializes 
769 // (1) The quenching model using quenching weights according to C. Salgado and U. Wiedemann
770 // (2) The nuclear geometry using the Glauber Model
771 //     
772     
773     fGlauber = AliFastGlauber::Instance();
774     fGlauber->Init(2);
775     fGlauber->SetCentralityClass(cMin, cMax); 
776
777     fQuenchingWeights = new AliQuenchingWeights();
778     fQuenchingWeights->InitMult();
779     fQuenchingWeights->SetK(k);
780     fQuenchingWeights->SetECMethod(AliQuenchingWeights::kECMethod(iECMethod));
781     fNGmax = ngmax;
782     fZmax  = zmax;
783     
784 }
785
786
787 void  AliPythia::Quench()
788 {
789 //
790 //
791 //  Simple Jet Quenching routine:
792 //  =============================
793 //  The jet formed by all final state partons radiated by the parton created 
794 //  in the hard collisions is quenched by a factor (1-z) using light cone variables in 
795 //  the initial parton reference frame:
796 //  (E + p_z)new = (1-z) (E + p_z)old
797 //
798 //
799 //
800 //
801 //   The lost momentum is first balanced by one gluon with virtuality > 0.   
802 //   Subsequently the gluon splits to yield two gluons with E = p.
803 //
804 //
805 // 
806     static Float_t eMean = 0.;
807     static Int_t   icall = 0;
808     
809     Double_t p0[4][5];
810     Double_t p1[4][5];
811     Double_t p2[4][5];
812     Int_t   klast[4] = {-1, -1, -1, -1};
813
814     Int_t numpart   = fPyjets->N;
815     Double_t px = 0., py = 0., pz = 0., e = 0., m = 0., p = 0., pt = 0., theta = 0., phi = 0.;
816     Double_t pxq[4], pyq[4], pzq[4], eq[4], yq[4], mq[4], pq[4], phiq[4], thetaq[4], ptq[4];
817     Bool_t  quenched[4];
818     Double_t wjtKick[4];
819     Int_t nGluon[4];
820     Int_t qPdg[4];
821     Int_t   imo, kst, pdg;
822     
823 //
824 //  Sore information about Primary partons
825 //
826 //  j =
827 //  0, 1 partons from hard scattering
828 //  2, 3 partons from initial state radiation
829 // 
830     for (Int_t i = 2; i <= 7; i++) {
831         Int_t j = 0;
832         // Skip gluons that participate in hard scattering
833         if (i == 4 || i == 5) continue;
834         // Gluons from hard Scattering
835         if (i == 6 || i == 7) {
836             j = i - 6;
837             pxq[j]    = fPyjets->P[0][i];
838             pyq[j]    = fPyjets->P[1][i];
839             pzq[j]    = fPyjets->P[2][i];
840             eq[j]     = fPyjets->P[3][i];
841             mq[j]     = fPyjets->P[4][i];
842         } else {
843             // Gluons from initial state radiation
844             //
845             // Obtain 4-momentum vector from difference between original parton and parton after gluon 
846             // radiation. Energy is calculated independently because initial state radition does not 
847             // conserve strictly momentum and energy for each partonic system independently.
848             //
849             // Not very clean. Should be improved !
850             //
851             //
852             j = i;
853             pxq[j]    = fPyjets->P[0][i] - fPyjets->P[0][i+2];
854             pyq[j]    = fPyjets->P[1][i] - fPyjets->P[1][i+2];
855             pzq[j]    = fPyjets->P[2][i] - fPyjets->P[2][i+2];
856             mq[j]     = fPyjets->P[4][i];
857             eq[j]     = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j] + mq[j] * mq[j]);
858         }
859 //
860 //  Calculate some kinematic variables
861 //
862         yq[j]     = 0.5 * TMath::Log((eq[j] + pzq[j] + 1.e-14) / (eq[j] - pzq[j] + 1.e-14));
863         pq[j]     = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j]);
864         phiq[j]   = TMath::Pi()+TMath::ATan2(-pyq[j], -pxq[j]);
865         ptq[j]    = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j]);
866         thetaq[j] = TMath::ATan2(ptq[j], pzq[j]);
867         qPdg[j]   =  fPyjets->K[1][i];
868     }
869   
870     Double_t int0[4];
871     Double_t int1[4];
872     
873     fGlauber->GetI0I1ForPythiaAndXY(4, phiq, int0, int1, fXJet, fYJet, 15.);
874
875     for (Int_t j = 0; j < 4; j++) {
876         //
877         // Quench only central jets and with E > 10.
878         //
879
880
881         Int_t itype = (qPdg[j] == 21) ? 2 : 1;
882         Double_t eloss = fQuenchingWeights->GetELossRandomKFast(itype, int0[j], int1[j], eq[j]);
883
884         if (TMath::Abs(yq[j]) > 2.5 || eq[j] < 10.) {
885             fZQuench[j] = 0.;
886         } else {
887             if (eq[j] > 40. && TMath::Abs(yq[j]) < 0.5) {
888                 icall ++;
889                 eMean += eloss;
890             }
891             //
892             // Extra pt
893             Double_t l =   fQuenchingWeights->CalcLk(int0[j], int1[j]);     
894             wjtKick[j] = TMath::Sqrt(l *  fQuenchingWeights->CalcQk(int0[j], int1[j]));
895             //
896             // Fractional energy loss
897             fZQuench[j] = eloss / eq[j];
898             //
899             // Avoid complete loss
900             //
901             if (fZQuench[j] > fZmax) fZQuench[j] = fZmax;
902             //
903             // Some debug printing
904
905             
906 //          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", 
907 //                 j, itype, eq[j], phiq[j], l, eloss, wjtKick[j], eMean / Float_t(icall+1), yq[j]);
908             
909 //          fZQuench[j] = 0.8;
910 //          while (fZQuench[j] >= 0.95)  fZQuench[j] = gRandom->Exp(0.2);
911         }
912         
913         quenched[j] = (fZQuench[j] > 0.01);
914     } // primary partons
915     
916     
917
918     Double_t pNew[1000][4];
919     Int_t    kNew[1000];
920     Int_t icount = 0;
921     Double_t zquench[4];
922     
923 //
924 //  System Loop    
925     for (Int_t isys = 0; isys < 4; isys++) {
926 //      Skip to next system if not quenched.
927         if (!quenched[isys]) continue;
928         
929         nGluon[isys]   = 1 + Int_t(fZQuench[isys] / (1. - fZQuench[isys]));
930         if (nGluon[isys] > fNGmax) nGluon[isys] = fNGmax;
931         zquench[isys] = 1. - TMath::Power(1. - fZQuench[isys], 1./Double_t(nGluon[isys]));
932         wjtKick[isys]  = wjtKick[isys] / TMath::Sqrt(Double_t(nGluon[isys]));
933
934
935         
936         Int_t igMin = -1;
937         Int_t igMax = -1;
938         Double_t pg[4] = {0., 0., 0., 0.};
939         
940 //
941 // Loop on radiation events
942
943         for (Int_t iglu = 0; iglu < nGluon[isys]; iglu++) {
944             while (1) {
945                 icount = 0;
946                 for (Int_t k = 0; k < 4; k++)
947                 {
948                     p0[isys][k] = 0.;
949                     p1[isys][k] = 0.;
950                     p2[isys][k] = 0.;
951                 }
952 //      Loop over partons
953                 for (Int_t i = 0; i < numpart; i++)
954                 {
955                     imo =  fPyjets->K[2][i];
956                     kst =  fPyjets->K[0][i];
957                     pdg =  fPyjets->K[1][i];
958                     
959                 
960                 
961 //      Quarks and gluons only
962                     if (pdg != 21 && TMath::Abs(pdg) > 6) continue;
963 //      Particles from hard scattering only
964                     
965                     if (imo > 8 && imo < 1000) imo = fPyjets->K[2][imo - 1];
966                     Int_t imom = imo % 1000;
967                     if ((isys == 0 || isys == 1) && ((imom != (isys + 7)))) continue;
968                     if ((isys == 2 || isys == 3) && ((imom != (isys + 1)))) continue;               
969                     
970                     
971 //      Skip comment lines
972                     if (kst != 1 && kst != 2) continue;
973 //
974 //      Parton kinematic
975                     px    = fPyjets->P[0][i];
976                     py    = fPyjets->P[1][i];
977                     pz    = fPyjets->P[2][i];
978                     e     = fPyjets->P[3][i];
979                     m     = fPyjets->P[4][i];
980                     pt    = TMath::Sqrt(px * px + py * py);
981                     p     = TMath::Sqrt(px * px + py * py + pz * pz); 
982                     phi   = TMath::Pi() + TMath::ATan2(-py, -px);
983                     theta = TMath::ATan2(pt, pz);
984                 
985 //
986 //      Save 4-momentum sum for balancing
987                     Int_t index = isys;
988                     
989                     p0[index][0] += px;
990                     p0[index][1] += py;
991                     p0[index][2] += pz;
992                     p0[index][3] += e;
993                 
994                     klast[index] = i;
995                     
996 //
997 //      Fractional energy loss
998                     Double_t z = zquench[index];
999                     
1000                     
1001 //      Don't fully quench radiated gluons
1002 //
1003                     if (imo > 1000) {
1004 //      This small factor makes sure that the gluons are not too close in phase space to avoid recombination
1005 //
1006
1007                         z = 0.02;
1008                     }
1009 //                  printf("z: %d %f\n", imo, z);
1010                     
1011
1012 //
1013                     
1014                     //
1015                     //
1016                     //      Transform into frame in which initial parton is along z-axis
1017                     //
1018                     TVector3 v(px, py, pz);
1019                     v.RotateZ(-phiq[index]);  v.RotateY(-thetaq[index]);
1020                     Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pl  = v.Z();
1021
1022                     Double_t jt  = TMath::Sqrt(pxs * pxs + pys * pys);
1023                     Double_t mt2 = jt * jt + m * m;
1024                     Double_t zmax = 1.;     
1025                     //
1026                     // Kinematic limit on z
1027                     //
1028                     if (m > 0.) zmax = 1. - m / TMath::Sqrt(m * m + jt * jt);
1029                     //
1030                     // Change light-cone kinematics rel. to initial parton
1031                     //  
1032                     Double_t eppzOld = e + pl;
1033                     Double_t empzOld = e - pl;
1034                     
1035                     Double_t eppzNew = (1. - z) * eppzOld;
1036                     Double_t empzNew = empzOld - mt2 * z / eppzOld;
1037                     Double_t eNew    = 0.5 * (eppzNew + empzNew);
1038                     Double_t plNew   = 0.5 * (eppzNew - empzNew);
1039                     
1040                     Double_t jtNew;
1041                     //
1042                     // if mt very small (or sometimes even < 0 for numerical reasons) set it to 0
1043                     Double_t mt2New = eppzNew * empzNew;
1044                     if (mt2New < 1.e-8) mt2New = 0.;
1045                     if (z < zmax) {
1046                         if (m * m > mt2New) {
1047                             //
1048                             // This should not happen 
1049                             //
1050                             Fatal("Quench()", "This should never happen %e %e %e!", m, eppzNew, empzNew);
1051                             jtNew = 0;
1052                         } else {
1053                             jtNew    = TMath::Sqrt(mt2New - m * m);
1054                         }
1055                     } else {
1056                         // If pT is to small (probably a leading massive particle) we scale only the energy
1057                         // This can cause negative masses of the radiated gluon
1058                         // Let's hope for the best ...
1059                         jtNew = jt;
1060                         eNew  = TMath::Sqrt(plNew * plNew + mt2);
1061                         
1062                     }
1063                     //
1064                     //     Calculate new px, py
1065                     //
1066                     Double_t pxNew   = 0;
1067                     Double_t pyNew   = 0;
1068                     
1069                     if (jt>0) {
1070                       pxNew = jtNew / jt * pxs;
1071                       pyNew = jtNew / jt * pys;
1072                     }   
1073 //                  Double_t dpx = pxs - pxNew;
1074 //                  Double_t dpy = pys - pyNew;
1075 //                  Double_t dpz = pl  - plNew;
1076 //                  Double_t de  = e   - eNew;
1077 //                  Double_t dmass2 = de * de  - dpx * dpx - dpy * dpy - dpz * dpz;
1078 //                  printf("New mass (1) %e %e %e %e %e %e %e \n", dmass2, jt, jtNew, pl, plNew, e, eNew);
1079 //                  printf("New mass (2) %e %e \n", pxNew, pyNew);
1080                     //
1081                     //      Rotate back
1082                     //  
1083                     TVector3 w(pxNew, pyNew, plNew);
1084                     w.RotateY(thetaq[index]); w.RotateZ(phiq[index]);
1085                     pxNew = w.X(); pyNew = w.Y(); plNew = w.Z();
1086                 
1087                     p1[index][0] += pxNew;
1088                     p1[index][1] += pyNew;
1089                     p1[index][2] += plNew;
1090                     p1[index][3] += eNew;       
1091                     //
1092                     // Updated 4-momentum vectors
1093                     //
1094                     pNew[icount][0]  = pxNew;
1095                     pNew[icount][1]  = pyNew;
1096                     pNew[icount][2]  = plNew;
1097                     pNew[icount][3]  = eNew;
1098                     kNew[icount]     = i;
1099                     icount++;
1100                 } // parton loop
1101                 //
1102                 // Check if there was phase-space for quenching
1103                 //
1104
1105                 if (icount == 0) quenched[isys] = kFALSE;
1106                 if (!quenched[isys]) break;
1107                 
1108                 for (Int_t j = 0; j < 4; j++) 
1109                 {
1110                     p2[isys][j] = p0[isys][j] - p1[isys][j];
1111                 }
1112                 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];
1113                 if (p2[isys][4] > 0.) {
1114                     p2[isys][4] = TMath::Sqrt(p2[isys][4]);
1115                     break;
1116                 } else {
1117                     printf("Warning negative mass squared in system %d %f ! \n", isys, zquench[isys]);
1118                     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]);
1119                     if (p2[isys][4] < -0.01) {
1120                         printf("Negative mass squared !\n");
1121                         // Here we have to put the gluon back to mass shell
1122                         // This will lead to a small energy imbalance
1123                         p2[isys][4]  = 0.;
1124                         p2[isys][3]  = TMath::Sqrt(p2[isys][0] * p2[isys][0] + p2[isys][1] * p2[isys][1] + p2[isys][2] * p2[isys][2]);
1125                         break;
1126                     } else {
1127                         p2[isys][4] = 0.;
1128                         break;
1129                     }
1130                 }
1131                 /*
1132                 zHeavy *= 0.98;
1133                 printf("zHeavy lowered to %f\n", zHeavy);
1134                 if (zHeavy < 0.01) {
1135                     printf("No success ! \n");
1136                     icount = 0;
1137                     quenched[isys] = kFALSE;
1138                     break;
1139                 }
1140                 */
1141             } // iteration on z (while)
1142             
1143 //          Update  event record
1144             for (Int_t k = 0; k < icount; k++) {
1145 //              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] );
1146                 fPyjets->P[0][kNew[k]] = pNew[k][0];
1147                 fPyjets->P[1][kNew[k]] = pNew[k][1];
1148                 fPyjets->P[2][kNew[k]] = pNew[k][2];
1149                 fPyjets->P[3][kNew[k]] = pNew[k][3];
1150             }
1151             //
1152             // Add the gluons
1153             //
1154             Int_t ish = 0;    
1155             Int_t iGlu;
1156             if (!quenched[isys]) continue;
1157 //
1158 //      Last parton from shower i
1159             Int_t in = klast[isys];
1160 //
1161 //      Continue if no parton in shower i selected
1162             if (in == -1) continue;
1163 //  
1164 //      If this is the second initial parton and it is behind the first move pointer by previous ish
1165             if (isys == 1 && klast[1] > klast[0]) in += ish;
1166 //
1167 //      Starting index
1168             
1169 //          jmin = in - 1;
1170 // How many additional gluons will be generated
1171             ish  = 1;
1172             if (p2[isys][4] > 0.05) ish = 2;
1173 //
1174 //      Position of gluons
1175             iGlu = numpart;
1176             if (iglu == 0) igMin = iGlu;
1177             igMax = iGlu;
1178             numpart += ish;
1179             (fPyjets->N) += ish;
1180             
1181             if (ish == 1) {
1182                 fPyjets->P[0][iGlu] = p2[isys][0];
1183                 fPyjets->P[1][iGlu] = p2[isys][1];
1184                 fPyjets->P[2][iGlu] = p2[isys][2];
1185                 fPyjets->P[3][iGlu] = p2[isys][3];
1186                 fPyjets->P[4][iGlu] = p2[isys][4];
1187                 
1188                 fPyjets->K[0][iGlu] = 1;
1189                 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu] = 1;
1190                 fPyjets->K[1][iGlu] = 21;       
1191                 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1192                 fPyjets->K[3][iGlu] = -1;       
1193                 fPyjets->K[4][iGlu] = -1;
1194                 
1195                 pg[0] += p2[isys][0];
1196                 pg[1] += p2[isys][1];
1197                 pg[2] += p2[isys][2];
1198                 pg[3] += p2[isys][3];
1199             } else {
1200                 //
1201                 // Split gluon in rest frame.
1202                 //
1203                 Double_t bx   =  p2[isys][0] / p2[isys][3];
1204                 Double_t by   =  p2[isys][1] / p2[isys][3];
1205                 Double_t bz   =  p2[isys][2] / p2[isys][3];
1206                 Double_t pst  =  p2[isys][4] / 2.;
1207                 //
1208                 // Isotropic decay ????
1209                 Double_t cost = 2. * gRandom->Rndm() - 1.;
1210                 Double_t sint = TMath::Sqrt((1.-cost)*(1.+cost));
1211                 Double_t phis =  2. * TMath::Pi() * gRandom->Rndm();
1212                 
1213                 Double_t pz1 =   pst * cost;
1214                 Double_t pz2 =  -pst * cost;
1215                 Double_t pt1 =   pst * sint;
1216                 Double_t pt2 =  -pst * sint;
1217                 Double_t px1 =   pt1 * TMath::Cos(phis);
1218                 Double_t py1 =   pt1 * TMath::Sin(phis);            
1219                 Double_t px2 =   pt2 * TMath::Cos(phis);
1220                 Double_t py2 =   pt2 * TMath::Sin(phis);            
1221                 
1222                 fPyjets->P[0][iGlu] = px1;
1223                 fPyjets->P[1][iGlu] = py1;
1224                 fPyjets->P[2][iGlu] = pz1;
1225                 fPyjets->P[3][iGlu] = pst;
1226                 fPyjets->P[4][iGlu] = 0.;
1227                 
1228                 fPyjets->K[0][iGlu] = 1 ;
1229                 fPyjets->K[1][iGlu] = 21;       
1230                 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1231                 fPyjets->K[3][iGlu] = -1;       
1232                 fPyjets->K[4][iGlu] = -1;
1233                 
1234                 fPyjets->P[0][iGlu+1] = px2;
1235                 fPyjets->P[1][iGlu+1] = py2;
1236                 fPyjets->P[2][iGlu+1] = pz2;
1237                 fPyjets->P[3][iGlu+1] = pst;
1238                 fPyjets->P[4][iGlu+1] = 0.;
1239                 
1240                 fPyjets->K[0][iGlu+1] = 1;
1241                 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu+1] = 1;
1242                 fPyjets->K[1][iGlu+1] = 21;     
1243                 fPyjets->K[2][iGlu+1] = fPyjets->K[2][in] + 1000;
1244                 fPyjets->K[3][iGlu+1] = -1;     
1245                 fPyjets->K[4][iGlu+1] = -1;
1246                 SetMSTU(1,0);
1247                 SetMSTU(2,0);
1248                 //
1249                 // Boost back
1250                 //
1251                 Pyrobo(iGlu + 1, iGlu + 2, 0., 0., bx, by, bz);
1252             }
1253 /*    
1254             for (Int_t ig = iGlu; ig < iGlu+ish; ig++) {
1255                 Double_t px, py, pz;
1256                 px = fPyjets->P[0][ig]; 
1257                 py = fPyjets->P[1][ig]; 
1258                 pz = fPyjets->P[2][ig]; 
1259                 TVector3 v(px, py, pz);
1260                 v.RotateZ(-phiq[isys]);
1261                 v.RotateY(-thetaq[isys]);
1262                 Double_t pxs     = v.X(); Double_t pys = v.Y(); Double_t pzs  = v.Z();     
1263                 Double_t r       = AliPythiaRndm::GetPythiaRandom()->Rndm();
1264                 Double_t jtKick  = 0.3 * TMath::Sqrt(-TMath::Log(r));
1265                 if (ish == 2)   jtKick  = wjtKick[i] * TMath::Sqrt(-TMath::Log(r)) / TMath::Sqrt(2.);
1266                 Double_t phiKick = 2. * TMath::Pi() * AliPythiaRndm::GetPythiaRandom()->Rndm();
1267                 pxs += jtKick * TMath::Cos(phiKick);
1268                 pys += jtKick * TMath::Sin(phiKick);
1269                 TVector3 w(pxs, pys, pzs);
1270                 w.RotateY(thetaq[isys]);
1271                 w.RotateZ(phiq[isys]);
1272                 fPyjets->P[0][ig] = w.X(); 
1273                 fPyjets->P[1][ig] = w.Y(); 
1274                 fPyjets->P[2][ig] = w.Z(); 
1275                 fPyjets->P[2][ig] = w.Mag();
1276             }
1277 */
1278         } // kGluon         
1279         
1280         
1281     // Check energy conservation
1282         Double_t pxs = 0.;
1283         Double_t pys = 0.;
1284         Double_t pzs = 0.;      
1285         Double_t es  = 14000.;
1286         
1287         for (Int_t i = 0; i < numpart; i++)
1288         {
1289             kst =  fPyjets->K[0][i];
1290             if (kst != 1 && kst != 2) continue;
1291             pxs += fPyjets->P[0][i];
1292             pys += fPyjets->P[1][i];
1293             pzs += fPyjets->P[2][i];        
1294             es  -= fPyjets->P[3][i];        
1295         }
1296         if (TMath::Abs(pxs) > 1.e-2 ||
1297             TMath::Abs(pys) > 1.e-2 ||
1298             TMath::Abs(pzs) > 1.e-1) {
1299             printf("%e %e %e %e\n", pxs, pys, pzs, es);
1300 //              Fatal("Quench()", "4-Momentum non-conservation");
1301         }
1302         
1303     } // end quenching loop (systems)
1304 // Clean-up
1305     for (Int_t i = 0; i < numpart; i++)
1306     {
1307         imo =  fPyjets->K[2][i];
1308         if (imo > 1000) {
1309             fPyjets->K[2][i] = fPyjets->K[2][i] % 1000;
1310         }
1311     }
1312 //      this->Pylist(1);
1313 } // end quench
1314
1315
1316 void AliPythia::Pyquen(Double_t a, Int_t ibf, Double_t b)
1317 {
1318     // Igor Lokthine's quenching routine
1319     // http://lokhtin.web.cern.ch/lokhtin/pyquen/pyquen.txt
1320
1321     pyquen(a, ibf, b);
1322 }
1323
1324 void AliPythia::SetPyquenParameters(Double_t t0, Double_t tau0, Int_t nf, Int_t iengl, Int_t iangl)
1325 {
1326     // Set the parameters for the PYQUEN package.
1327     // See comments in PyquenCommon.h
1328     
1329     
1330     PYQPAR.t0    = t0;
1331     PYQPAR.tau0  = tau0;
1332     PYQPAR.nf    = nf;
1333     PYQPAR.iengl = iengl;
1334     PYQPAR.iangl = iangl;
1335 }
1336
1337
1338 void AliPythia::Pyevnw()
1339 {
1340     // New multiple interaction scenario
1341     pyevnw();
1342 }
1343
1344 void  AliPythia::Pyshowq(Int_t ip1, Int_t ip2, Double_t qmax)
1345 {
1346     //  Call medium-modified Pythia jet reconstruction algorithm
1347     //
1348     pyshowq(ip1, ip2, qmax);
1349 }
1350  void AliPythia::Qpygin0()                                                                      
1351  {                                                                                              
1352      // New multiple interaction scenario                                                       
1353      qpygin0();                                                                                 
1354  }
1355
1356 void AliPythia::GetQuenchingParameters(Double_t& xp, Double_t& yp, Double_t z[4])
1357 {
1358     // Return event specific quenching parameters
1359     xp = fXJet;
1360     yp = fYJet;
1361     for (Int_t i = 0; i < 4; i++) z[i] = fZQuench[i];
1362
1363 }
1364
1365 void AliPythia::ConfigHeavyFlavor()
1366 {
1367     //
1368     // Default configuration for Heavy Flavor production
1369     //
1370     // All QCD processes
1371     //
1372     SetMSEL(1);
1373     
1374     // No multiple interactions
1375     SetMSTP(81,0);
1376     SetPARP(81, 0.);
1377     SetPARP(82, 0.);    
1378     // Initial/final parton shower on (Pythia default)
1379     SetMSTP(61,1);
1380     SetMSTP(71,1);
1381     
1382     // 2nd order alpha_s
1383     SetMSTP(2,2);
1384     
1385     // QCD scales
1386     SetMSTP(32,2);
1387     SetPARP(34,1.0);
1388 }
1389
1390 void AliPythia::AtlasTuning()
1391 {
1392     //
1393     // Configuration for the ATLAS tuning
1394     if (fItune > -1) return;
1395     printf("ATLAS TUNE \n");
1396     
1397     SetMSTP(51, AliStructFuncType::PDFsetIndex(kCTEQ5L));      // CTEQ5L pdf
1398     SetMSTP(81,1);             // Multiple Interactions ON
1399     SetMSTP(82,4);             // Double Gaussian Model
1400     SetPARP(81,1.9);           // Min. pt for multiple interactions (default in 6.2-14) 
1401     SetPARP(82,1.8);           // [GeV]    PT_min at Ref. energy
1402     SetPARP(89,1000.);         // [GeV]   Ref. energy
1403     SetPARP(90,0.16);          // 2*epsilon (exponent in power law)
1404     SetPARP(83,0.5);           // Core density in proton matter distribution (def.value)
1405     SetPARP(84,0.5);           // Core radius
1406     SetPARP(85,0.33);          // Regulates gluon prod. mechanism
1407     SetPARP(86,0.66);          // Regulates gluon prod. mechanism
1408     SetPARP(67,1);             // Regulates Initial State Radiation
1409 }
1410
1411 void AliPythia::AtlasTuning_MC09()
1412 {
1413     //
1414     // Configuration for the ATLAS tuning
1415     if (fItune > -1) return;
1416     printf("ATLAS New TUNE MC09\n");
1417     SetMSTP(81,21);             // treatment for MI, ISR, FSR and beam remnants: MI on, new model
1418     SetMSTP(82, 4);             // Double Gaussian Model
1419     SetMSTP(52, 2);             // External PDF
1420     SetMSTP(51, 20650);         // MRST LO*
1421   
1422     
1423     SetMSTP(70, 0);             // (was 2: def manual 1, def code 0) virtuality scale for ISR 
1424     SetMSTP(72, 1);             // (was 0: def 1) maximum scale for FSR
1425     SetMSTP(88, 1);             // (was 0: def 1) strategy for qq junction to di-quark or baryon in beam remnant
1426     SetMSTP(90, 0);             // (was 1: def 0) strategy of compensate the primordial kT
1427
1428     SetPARP(78, 0.3);           // the amount of color reconnection in the final state
1429     SetPARP(80, 0.1);           // probability of color partons kicked out from beam remnant
1430     SetPARP(82, 2.3);           // [GeV]    PT_min at Ref. energy    
1431     SetPARP(83, 0.8);           // Core density in proton matter distribution (def.value)    
1432     SetPARP(84, 0.7);           // Core radius
1433     SetPARP(90, 0.25);          //  2*epsilon (exponent in power law)
1434     SetPARJ(81, 0.29);          // (was 0.14: def 0.29) Labmda value in running alpha_s for parton showers
1435
1436     SetMSTP(95, 6);
1437     SetPARJ(41, 0.3);           // a and b parameters of the symmm. Lund FF
1438     SetPARJ(42, 0.58);
1439     SetPARJ(46, 0.75);          // mod. of the Lund FF for heavy end-point quarks
1440     SetPARP(89,1800.);         // [GeV]   Ref. energy
1441 }
1442
1443 AliPythia& AliPythia::operator=(const  AliPythia& rhs)
1444 {
1445 // Assignment operator
1446     rhs.Copy(*this);
1447     return *this;
1448 }
1449
1450  void AliPythia::Copy(TObject&) const
1451 {
1452     //
1453     // Copy 
1454     //
1455     Fatal("Copy","Not implemented!\n");
1456 }
1457
1458 void AliPythia::DalitzDecays()
1459 {
1460
1461   //
1462   // Replace all omega dalitz decays with the correct matrix element decays
1463   //
1464   Int_t nt = fPyjets->N;
1465   for (Int_t i = 0; i < nt; i++) {
1466     if (fPyjets->K[1][i] != 223)                       continue;        
1467     Int_t fd = fPyjets->K[3][i] - 1;
1468     Int_t ld = fPyjets->K[4][i] - 1;
1469     if (fd < 0)                                        continue;
1470     if ((ld - fd) != 2)                                continue;
1471     if ((fPyjets->K[1][fd] != 111) ||
1472         (TMath::Abs(fPyjets->K[1][fd+1]) != 11))       continue;
1473     TLorentzVector omega(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1474     fOmegaDalitz.Decay(223, &omega);
1475     for (Int_t j = 0; j < 3; j++) {
1476       for (Int_t k = 0; k < 4; k++) {
1477         TLorentzVector vec = (fOmegaDalitz.Products())[2-j];
1478         fPyjets->P[k][fd+j] = vec[k];
1479       }
1480     }
1481   }
1482 }