2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
22 #include "AliPythia8.h"
25 #include "AliPythiaRndm.h"
30 // Particles produced in string fragmentation point directly to either of the two endpoints
31 // of the string (depending in the side they were generated from).
32 // SetMSTU(16,2); // ????
33 // String drawing almost completely minimizes string length.
34 // Probability that an additional interaction gives two gluons
35 // ... with color connection to nearest neighbours
37 // ... as closed gluon loop
41 // Baryon production model
43 // String fragmentation
45 // sea quarks can be used for baryon formation
47 // choice of max. virtuality for ISR
49 // regularisation scheme of ISR
51 // all resonance decays switched on
53 AliPythia8* AliPythia8::fgAliPythia8=NULL;
55 AliPythia8::AliPythia8():
69 // Default Constructor
72 if (!AliPythiaRndm::GetPythiaRandom())
73 AliPythiaRndm::SetPythiaRandom(GetRandom());
76 AliPythia8::AliPythia8(const AliPythia8& pythia):
94 void AliPythia8::ProcInit(Process_t process, Float_t energy, StrucFunc_t strucfunc)
96 // Initialise the process to generate
97 if (!AliPythiaRndm::GetPythiaRandom())
98 AliPythiaRndm::SetPythiaRandom(GetRandom());
102 fStrucFunc = strucfunc;
103 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
104 ReadString("111:mayDecay = off");
105 ReadString("310:mayDecay = off");
106 ReadString("3122:mayDecay = off");
107 ReadString("3112:mayDecay = off");
108 ReadString("3212:mayDecay = off");
109 ReadString("3222:mayDecay = off");
110 ReadString("3312:mayDecay = off");
111 ReadString("3322:mayDecay = off");
112 ReadString("3334:mayDecay = off");
113 // Select structure function
115 ReadString("PDF:useLHAPDF = on");
116 ReadString(Form("PDF:LHAPDFset = %s", AliStructFuncType::PDFsetName(fStrucFunc).Data()));
118 // Particles produced in string fragmentation point directly to either of the two endpoints
119 // of the string (depending in the side they were generated from).
121 // SetMSTU(16,2); // ????
124 // Pythia initialisation for selected processes//
128 case kPyOldUEQ2ordered: //Old underlying events with Q2 ordered QCD processes
129 // Multiple interactions on.
130 ReadString("PartonLevel:MI = on");
131 // Double Gaussian matter distribution.
132 ReadString("MultipleInteractions:bProfile = 2");
133 ReadString("MultipleInteractions:coreFraction = 0.5");
134 ReadString("MultipleInteractions:coreRadius = 0.4");
136 ReadString("MultipleInteractions:pTmin = 2.0");
137 // Reference energy for pT0 and energy rescaling pace.
138 ReadString("MultipleInteractions:ecmRef = 1800.");
139 ReadString("MultipleInteractions:ecmPow = 0.25");
140 // String drawing almost completely minimizes string length.
143 // ISR and FSR activity.
144 // Q^2 scale of the hard scattering
145 ReadString("SigmaProcess:factorMultFac = 4.");
147 // SetPARJ(81, 0.29);
149 case kPyOldUEQ2ordered2:
150 // Old underlying events with Q2 ordered QCD processes
151 // Multiple interactions on.
152 ReadString("PartonLevel:MI = on");
153 // Double Gaussian matter distribution.
154 ReadString("MultipleInteractions:bProfile = 2");
155 ReadString("MultipleInteractions:coreFraction = 0.5");
156 ReadString("MultipleInteractions:coreRadius = 0.4");
158 ReadString("MultipleInteractions:pTmin = 2.0");
159 // Reference energy for pT0 and energy rescaling pace.
160 ReadString("MultipleInteractions:ecmRef = 1800.");
161 ReadString("MultipleInteractions:ecmPow = 0.16");
162 // String drawing almost completely minimizes string length.
165 // ISR and FSR activity.
166 ReadString("SigmaProcess:factorMultFac = 4.");
171 // Old production mechanism: Old Popcorn
172 ReadString("HardQCD:all = on");
174 // (D=2) Like MSTJ(12)=2 but added prod ofthe 1er rank baryon
176 // (D=1)see can be used to form baryons (BARYON JUNCTION)
181 ReadString("HardQCD:gg2ccbar = on");
182 ReadString("HardQCD:qqbar2ccbar = on");
183 // heavy quark masses
184 ReadString("ParticleData:mcRun = 1.2");
187 ReadString("BeamRemnants:primordialKT = on");
188 ReadString("BeamRemnants:primordialKTsoft = 0.");
189 ReadString("BeamRemnants:primordialKThard = 1.");
190 ReadString("BeamRemnants:halfScaleForKT = 0.");
191 ReadString("BeamRemnants:halfMassForKT = 0.");
194 ReadString("HardQCD:gg2bbbar = on");
195 ReadString("HardQCD:qqbar2bbbar = on");
196 ReadString("ParticleData:mbRun = 4.75");
200 ReadString("Charmonium:gg2QQbar[3S1(1)]g = on");
203 ReadString("Charmonium:all = on");
205 case kPyCharmUnforced:
207 ReadString("HardQCD:gq2qg = on");
209 ReadString("HardQCD:gg2qq = on");
211 ReadString("HardQCD:gg2gg = on");
213 case kPyBeautyUnforced:
215 ReadString("HardQCD:gq2qg = on");
217 ReadString("HardQCD:gg2qq = on");
219 ReadString("HardQCD:gg2gg = on");
222 // Minimum Bias pp-Collisions
225 // select Pythia min. bias model
226 // single diffraction AB-->XB
227 ReadString("SoftQCD:minBias = on");
228 ReadString("SoftQCD:singleDiffractive = on");
229 ReadString("SoftQCD:doubleDiffractive = on");
233 // Minimum Bias pp-Collisions
236 // select Pythia min. bias model
237 ReadString("SoftQCD:minBias = on");
238 ReadString("SoftQCD:singleDiffractive = on");
239 ReadString("SoftQCD:doubleDiffractive = on");
242 // Les Houches Working Group 05 Minimum Bias pp-Collisions: hep-ph/0604120
243 // -> Pythia 6.3 or above is needed
245 ReadString("SoftQCD:minBias = on");
246 ReadString("SoftQCD:singleDiffractive = on");
247 ReadString("SoftQCD:doubleDiffractive = on");
248 ReadString(Form("PDF:LHAPDFset = %s", AliStructFuncType::PDFsetName(kCTEQ6ll).Data()));
252 // ReadString("PartonLevel:MI = on");
253 // Double Gaussian matter distribution.
254 ReadString("MultipleInteractions:bProfile = 2");
255 ReadString("MultipleInteractions:coreFraction = 0.5");
256 ReadString("MultipleInteractions:coreRadius = 0.5");
257 ReadString("MultipleInteractions:expPow = 0.16");
258 ReadString("MultipleInteractions:pTmin = 2.3");
260 // SetPARP(85,0.9); // Regulates gluon prod. mechanism
263 // Minimum Bias pp-Collisions
266 // select Pythia min. bias model
267 ReadString("SoftQCD:minBias = on");
273 ReadString("BeamRemnants:primordialKT = on");
274 ReadString("BeamRemnants:primordialKTsoft = 0.");
275 ReadString("BeamRemnants:primordialKThard = 1.");
276 ReadString("BeamRemnants:halfScaleForKT = 0.");
277 ReadString("BeamRemnants:halfMassForKT = 0.");
279 ReadString("ParticleData:mcRun = 1.20");
280 ReadString("ParticleData:mbRun = 4.78");
288 ReadString("HardQCD:all = on");
290 // Pythia Tune A (CDF)
292 ReadString("PartonLevel:MI = on");
293 ReadString("MultipleInteractions:pTmin = 2.0");
294 ReadString("MultipleInteractions:pT0Ref = 2.8");
295 ReadString("MultipleInteractions:ecmRef = 1800.");
296 ReadString("MultipleInteractions:expPow = 0.25");
297 ReadString("MultipleInteractions:bProfile = 2");
298 ReadString("MultipleInteractions:coreFraction = 0.16");
299 ReadString("MultipleInteractions:coreRadius = 0.4");
300 ReadString("SigmaProcess:factorMultFac = 2.5");
301 // SetPARP(85,0.90) ; // Regulates gluon prod. mechanism
302 // SetPARP(86,0.95); // Regulates gluon prod. mechanism
305 ReadString("PromptPhoton:all = on");
307 case kPyCharmPbPbMNR:
309 case kPyDPlusPbPbMNR:
310 case kPyDPlusStrangePbPbMNR:
311 // Tuning of Pythia parameters aimed to get a resonable agreement
312 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
313 // c-cbar single inclusive and double differential distributions.
314 // This parameter settings are meant to work with Pb-Pb collisions
315 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
316 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
317 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
320 ReadString("BeamRemnants:primordialKT = on");
321 ReadString("BeamRemnants:primordialKTsoft = 0.");
322 ReadString("BeamRemnants:primordialKThard = 1.304");
323 ReadString("BeamRemnants:halfScaleForKT = 0.");
324 ReadString("BeamRemnants:halfMassForKT = 0.");
326 ReadString("ParticleData:mcRun = 1.20");
331 case kPyDPlusStrangepPbMNR:
332 // Tuning of Pythia parameters aimed to get a resonable agreement
333 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
334 // c-cbar single inclusive and double differential distributions.
335 // This parameter settings are meant to work with p-Pb collisions
336 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
337 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
338 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
341 ReadString("BeamRemnants:primordialKT = on");
342 ReadString("BeamRemnants:primordialKTsoft = 0.");
343 ReadString("BeamRemnants:primordialKThard = 1.16");
344 ReadString("BeamRemnants:halfScaleForKT = 0.");
345 ReadString("BeamRemnants:halfMassForKT = 0.");
347 ReadString("ParticleData:mcRun = 1.20");
352 case kPyDPlusStrangeppMNR:
353 case kPyLambdacppMNR:
354 // Tuning of Pythia parameters aimed to get a resonable agreement
355 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
356 // c-cbar single inclusive and double differential distributions.
357 // This parameter settings are meant to work with pp collisions
358 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
359 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
360 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
363 ReadString("BeamRemnants:primordialKT = on");
364 ReadString("BeamRemnants:primordialKTsoft = 0.");
365 ReadString("BeamRemnants:primordialKThard = 1.");
366 ReadString("BeamRemnants:halfScaleForKT = 0.");
367 ReadString("BeamRemnants:halfMassForKT = 0.");
369 ReadString("ParticleData:mcRun = 1.20");
371 case kPyCharmppMNRwmi:
372 // Tuning of Pythia parameters aimed to get a resonable agreement
373 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
374 // c-cbar single inclusive and double differential distributions.
375 // This parameter settings are meant to work with pp collisions
376 // and with kCTEQ5L PDFs.
377 // Added multiple interactions according to ATLAS tune settings.
378 // To get a "reasonable" agreement with MNR results, events have to be
379 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
381 // To get a "perfect" agreement with MNR results, events have to be
382 // generated in four ptHard bins with the following relative
390 ReadString("BeamRemnants:primordialKT = on");
391 ReadString("BeamRemnants:primordialKTsoft = 0.");
392 ReadString("BeamRemnants:primordialKThard = 1.");
393 ReadString("BeamRemnants:halfScaleForKT = 0.");
394 ReadString("BeamRemnants:halfMassForKT = 0.");
396 ReadString("ParticleData:mcRun = 1.20");
399 case kPyBeautyPbPbMNR:
400 // Tuning of Pythia parameters aimed to get a resonable agreement
401 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
402 // b-bbar single inclusive and double differential distributions.
403 // This parameter settings are meant to work with Pb-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.75GeV. Example in ConfigBeautyPPR.C.
409 ReadString("SigmaProcess:factorMultFac = 1.");
411 ReadString("BeamRemnants:primordialKT = on");
412 ReadString("BeamRemnants:primordialKTsoft = 0.");
413 ReadString("BeamRemnants:primordialKThard = 2.035");
414 ReadString("BeamRemnants:halfScaleForKT = 0.");
415 ReadString("BeamRemnants:halfMassForKT = 0.");
417 ReadString("ParticleData:mbRun = 4.75");
419 case kPyBeautypPbMNR:
420 // Tuning of Pythia parameters aimed to get a resonable agreement
421 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
422 // b-bbar single inclusive and double differential distributions.
423 // This parameter settings are meant to work with p-Pb 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.75GeV. Example in ConfigBeautyPPR.C.
429 ReadString("SigmaProcess:factorMultFac = 1.");
431 ReadString("BeamRemnants:primordialKT = on");
432 ReadString("BeamRemnants:primordialKTsoft = 0.");
433 ReadString("BeamRemnants:primordialKThard = 1.6");
434 ReadString("BeamRemnants:halfScaleForKT = 0.");
435 ReadString("BeamRemnants:halfMassForKT = 0.");
437 ReadString("ParticleData:mbRun = 4.75");
440 // Tuning of Pythia parameters aimed to get a resonable agreement
441 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
442 // b-bbar single inclusive and double differential distributions.
443 // This parameter settings are meant to work with pp collisions
444 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
445 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
446 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
449 ReadString("SigmaProcess:factorMultFac = 1.");
451 ReadString("BeamRemnants:primordialKT = on");
452 ReadString("BeamRemnants:primordialKTsoft = 0.");
453 ReadString("BeamRemnants:primordialKThard = 1.0");
454 ReadString("BeamRemnants:halfScaleForKT = 0.");
455 ReadString("BeamRemnants:halfMassForKT = 0.");
457 ReadString("ParticleData:mbRun = 4.75");
459 case kPyBeautyppMNRwmi:
460 // Tuning of Pythia parameters aimed to get a resonable agreement
461 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
462 // b-bbar single inclusive and double differential distributions.
463 // This parameter settings are meant to work with pp collisions
464 // and with kCTEQ5L PDFs.
465 // Added multiple interactions according to ATLAS tune settings.
466 // To get a "reasonable" agreement with MNR results, events have to be
467 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
469 // To get a "perfect" agreement with MNR results, events have to be
470 // generated in four ptHard bins with the following relative
478 ReadString("SigmaProcess:factorMultFac = 1.");
480 ReadString("BeamRemnants:primordialKT = on");
481 ReadString("BeamRemnants:primordialKTsoft = 0.");
482 ReadString("BeamRemnants:primordialKThard = 1.0");
483 ReadString("BeamRemnants:halfScaleForKT = 0.");
484 ReadString("BeamRemnants:halfMassForKT = 0.");
486 ReadString("ParticleData:mbRun = 4.75");
490 //Inclusive production of W+/-
492 ReadString("WeakSingleBoson:ffbar2W = on");
493 // Initial/final parton shower on (Pythia default)
494 // With parton showers on we are generating "W inclusive process"
495 ReadString("PartonLevel:ISR = on");
496 ReadString("PartonLevel:FSR = on");
499 //Inclusive production of Z
501 ReadString("WeakSingleBoson:ffbar2gmZ = on");
502 //only Z included, not gamma
503 ReadString("WeakZ0:gmZmode = 2");
504 // Initial/final parton shower on (Pythia default)
505 // With parton showers on we are generating "Z inclusive process"
506 ReadString("PartonLevel:ISR = on");
507 ReadString("PartonLevel:FSR = on");
512 // SetMSTP(41,1); // all resonance decays switched on
513 Initialize(2212, 2212, fEcms);
516 void AliPythia8::SetNuclei(Int_t /*a1*/, Int_t /*a2*/)
518 // Treat protons as inside nuclei with mass numbers a1 and a2
519 // The MSTP array in the PYPARS common block is used to enable and
520 // select the nuclear structure functions.
521 // MSTP(52) : (D=1) choice of proton and nuclear structure-function library
522 // =1: internal PYTHIA acording to MSTP(51)
523 // =2: PDFLIB proton s.f., with MSTP(51) = 1000xNGROUP+NSET
524 // If the following mass number both not equal zero, nuclear corrections of the stf are used.
525 // MSTP(192) : Mass number of nucleus side 1
526 // MSTP(193) : Mass number of nucleus side 2
533 AliPythia8* AliPythia8::Instance()
535 // Set random number generator
539 fgAliPythia8 = new AliPythia8();
544 void AliPythia8::PrintParticles()
546 // Print list of particl properties
547 ReadString("Main:showAllParticleData");
550 void AliPythia8::ResetDecayTable()
552 // Set default values for pythia decay switches
554 // for (i = 1; i < 501; i++) SetMDCY(i,1,fDefMDCY[i]);
555 // for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
558 void AliPythia8::SetDecayTable()
560 // Set default values for pythia decay switches
563 // for (i = 1; i < 501; i++) fDefMDCY[i] = GetMDCY(i,1);
564 // for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
567 void AliPythia8::Pyclus(Int_t& njet)
569 // Call Pythia clustering algorithm
571 Bool_t ok = fClusterJet.analyze(Pythia8()->event, fYScale, fPtScale, fNJetMin, fNJetMax);
573 if (ok) njet = fClusterJet.size();
576 void AliPythia8::Pycell(Int_t& njet)
578 // Call Pythia jet reconstruction algorithm
580 Bool_t ok = fCellJet.analyze(Pythia8()->event, fMinEtJet, fRJet, fEtSeed);
582 if (ok) njet = fCellJet.size();
585 void AliPythia8::GetJet(Int_t i, Float_t& px, Float_t& py, Float_t& pz, Float_t& e)
588 Float_t et = fCellJet.eT(i);
589 px = et * TMath::Cos(fCellJet.phiWeighted(i));
590 py = et * TMath::Sin(fCellJet.phiWeighted(i));
591 pz = et * TMath::SinH(fCellJet.etaWeighted(i));
592 e = et * TMath::CosH(fCellJet.etaWeighted(i));
595 void AliPythia8::GenerateEvent()
597 // Generate one event
598 AliTPythia8::GenerateEvent();
601 void AliPythia8::GenerateMIEvent()
603 // New multiple interaction scenario
604 AliWarning("Not implemented. No event will be generated");
607 void AliPythia8::PrintStatistics()
609 // End of run statistics
610 AliTPythia8::PrintStatistics();
613 void AliPythia8::EventListing()
615 // End of run statistics
616 AliTPythia8::EventListing();
619 Int_t AliPythia8::ProcessCode()
621 // Returns the subprocess code for the current event
622 printf("Process %5d %5d \n", Pythia8()->info.code(), Pythia8()->info.codeSub());
624 return Pythia8()->info.codeSub();
627 void AliPythia8::ConfigHeavyFlavor()
630 // Default configuration for Heavy Flavor production
634 ReadString("HardQCD:all = on");
636 // No multiple interactions
637 ReadString("PartonLevel:MI = off");
638 ReadString("MultipleInteractions:pTmin = 0.0");
639 ReadString("MultipleInteractions:pT0Ref = 0.0");
641 // Initial/final parton shower on (Pythia default)
642 ReadString("PartonLevel:ISR = on");
643 ReadString("PartonLevel:FSR = on");
646 ReadString("SigmaProcess:alphaSorder = 2");
649 ReadString("SigmaProcess:renormScale2 = 2");
650 ReadString("SigmaProcess:renormMultFac = 1.");
653 void AliPythia8::AtlasTuning()
656 // Configuration for the ATLAS tuning
657 ReadString(Form("PDF:LHAPDFset = %s", AliStructFuncType::PDFsetName(kCTEQ5L).Data()));
658 ReadString("PartonLevel:MI = on");
659 ReadString("MultipleInteractions:pTmin = 1.9");
660 ReadString("MultipleInteractions:pT0Ref = 1.8");
661 ReadString("MultipleInteractions:ecmRef = 1000.");
662 ReadString("MultipleInteractions:expPow = 0.16");
663 ReadString("MultipleInteractions:bProfile = 2");
664 ReadString("MultipleInteractions:coreFraction = 0.16");
665 ReadString("MultipleInteractions:coreRadius = 0.5");
666 // SetPARP(85,0.33); // Regulates gluon prod. mechanism
667 // SetPARP(86,0.66); // Regulates gluon prod. mechanism
668 ReadString("SigmaProcess:factorMultFac = 1.");
671 void AliPythia8::SetPtHardRange(Float_t ptmin, Float_t ptmax)
673 // Set the pt hard range
674 ReadString(Form("PhaseSpace:pTHatMin = %13.3f", ptmin));
675 ReadString(Form("PhaseSpace:pTHatMax = %13.3f", ptmax));
678 void AliPythia8::SetYHardRange(Float_t /*ymin*/, Float_t /*ymax*/)
680 // Set the y hard range
681 printf("YHardRange not implemented in Pythia8 !!!\n");
686 void AliPythia8::SetFragmentation(Int_t flag)
688 // Switch fragmentation on/off
690 ReadString("HadronLevel:Hadronize = on");
692 ReadString("HadronLevel:Hadronize = off");
696 void AliPythia8::SetInitialAndFinalStateRadiation(Int_t flag1, Int_t flag2)
698 // initial state radiation
700 ReadString("PartonLevel:ISR = on");
702 ReadString("PartonLevel:ISR = off");
704 // final state radiation
706 ReadString("PartonLevel:FSR = on");
708 ReadString("PartonLevel:FSR = off");
712 void AliPythia8::SetIntrinsicKt(Float_t kt)
714 ReadString("BeamRemnants:primordialKT = on");
715 ReadString("BeamRemnants:primordialKTsoft = 0.");
716 ReadString(Form("BeamRemnants:primordialKThard = %13.3f", kt));
717 ReadString("BeamRemnants:halfScaleForKT = 0.");
718 ReadString("BeamRemnants:halfMassForKT = 0.");
721 void AliPythia8::SwitchHFOff()
723 // Switch off heavy flavor
724 // Maximum number of quark flavours used in pdf
725 ReadString("PDFinProcess:nQuarkIn = 3");
726 // Maximum number of flavors that can be used in showers
727 ReadString("TimeShower:nGluonToQuark = 3");
728 ReadString("SpaceShower:nQuarkIn = 3");
733 void AliPythia8::SetPycellParameters(Float_t etaMax, Int_t nEta, Int_t nPhi,
734 Float_t thresh, Float_t etseed, Float_t minet, Float_t r)
736 // Set pycell parameters
737 fCellJet = Pythia8::CellJet( etaMax, nEta, nPhi, 2, 0, 0., 0., thresh);
743 void AliPythia8::ModifiedSplitting()
746 // We have to see how to implement this in Pythia8 !!!
748 // Modified splitting probability as a model for quenching
749 // SetPARJ(200, 0.8);
750 // SetMSTJ(41, 1); // QCD radiation only
751 // SetMSTJ(42, 2); // angular ordering
752 // SetMSTJ(44, 2); // option to run alpha_s
753 // SetMSTJ(47, 0); // No correction back to hard scattering element
754 // SetMSTJ(50, 0); // No coherence in first branching
755 // SetPARJ(82, 1.); // Cut off for parton showers
759 void AliPythia8::InitQuenching(Float_t /*cMin*/, Float_t /*cMax*/, Float_t /*k*/, Int_t /*iECMethod*/, Float_t /*zmax*/, Int_t /*ngmax*/)
763 AliWarning("Not implemented !");
766 void AliPythia8::SwitchHadronisationOff()
768 // Switch off hadronisation
769 ReadString("HadronLevel:Hadronize = off");
772 void AliPythia8::SwitchHadronisationOn()
774 // Switch on hadronisarion
775 ReadString("HadronLevel:Hadronize = on");
779 void AliPythia8::GetXandQ(Float_t& x1, Float_t& x2, Float_t& q)
781 // Get x1, x2 and Q for this event
783 q = Pythia8()->info.QFac();
784 x1 = Pythia8()->info.x1();
785 x2 = Pythia8()->info.x2();
789 Float_t AliPythia8::GetXSection()
791 // Get the total cross-section
792 return Pythia8()->info.sigmaGen();
795 Float_t AliPythia8::GetPtHard()
797 // Get the pT hard for this event
798 return Pythia8()->info.pTHat();
804 AliPythia8& AliPythia8::operator=(const AliPythia8& rhs)
806 // Assignment operator
811 void AliPythia8::Copy(TObject&) const
816 Fatal("Copy","Not implemented!\n");
822 void AliPythia8::SetNumberOfParticles(Int_t /*i*/)
824 AliWarning("Not implemented");
827 void AliPythia8::EditEventList(Int_t /*i*/)
829 AliWarning("Not implemented");
832 void AliPythia8::Pyquen(Double_t /*a*/, Int_t /*b*/, Double_t /*c*/)
834 AliWarning("Cannot be used with Pythia8");
837 void AliPythia8::HadronizeEvent()
839 // Needs access to HadronLevel ?
840 AliWarning("Not yet implemented");
843 void AliPythia8::GetQuenchingParameters(Double_t& /*xp*/, Double_t& /*yp*/, Double_t* /*z[4]*/)
845 AliWarning("Not yet implemented");
848 void AliPythia8::LoadEvent(AliStack* /*stack*/, Int_t /*flag*/, Int_t /*reHadr*/)
850 AliWarning("Not yet implemented");