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