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("Beams:primordialKT = on");
188 ReadString("Beams:primordialKTsoft = 0.");
189 ReadString("Beams:primordialKThard = 1.");
190 ReadString("Beams:halfScaleForKT = 0.");
191 ReadString("Beams: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("Beams:primordialKT = on");
274 ReadString("Beams:primordialKTsoft = 0.");
275 ReadString("Beams:primordialKThard = 1.");
276 ReadString("Beams:halfScaleForKT = 0.");
277 ReadString("Beams: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("Beams:primordialKT = on");
321 ReadString("Beams:primordialKTsoft = 0.");
322 ReadString("Beams:primordialKThard = 1.304");
323 ReadString("Beams:halfScaleForKT = 0.");
324 ReadString("Beams: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("Beams:primordialKT = on");
342 ReadString("Beams:primordialKTsoft = 0.");
343 ReadString("Beams:primordialKThard = 1.16");
344 ReadString("Beams:halfScaleForKT = 0.");
345 ReadString("Beams:halfMassForKT = 0.");
347 ReadString("ParticleData:mcRun = 1.20");
352 case kPyDPlusStrangeppMNR:
353 // Tuning of Pythia parameters aimed to get a resonable agreement
354 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
355 // c-cbar single inclusive and double differential distributions.
356 // This parameter settings are meant to work with pp collisions
357 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
358 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
359 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
362 ReadString("Beams:primordialKT = on");
363 ReadString("Beams:primordialKTsoft = 0.");
364 ReadString("Beams:primordialKThard = 1.");
365 ReadString("Beams:halfScaleForKT = 0.");
366 ReadString("Beams:halfMassForKT = 0.");
368 ReadString("ParticleData:mcRun = 1.20");
370 case kPyCharmppMNRwmi:
371 // Tuning of Pythia parameters aimed to get a resonable agreement
372 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
373 // c-cbar single inclusive and double differential distributions.
374 // This parameter settings are meant to work with pp collisions
375 // and with kCTEQ5L PDFs.
376 // Added multiple interactions according to ATLAS tune settings.
377 // To get a "reasonable" agreement with MNR results, events have to be
378 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
380 // To get a "perfect" agreement with MNR results, events have to be
381 // generated in four ptHard bins with the following relative
389 ReadString("Beams:primordialKT = on");
390 ReadString("Beams:primordialKTsoft = 0.");
391 ReadString("Beams:primordialKThard = 1.");
392 ReadString("Beams:halfScaleForKT = 0.");
393 ReadString("Beams:halfMassForKT = 0.");
395 ReadString("ParticleData:mcRun = 1.20");
398 case kPyBeautyPbPbMNR:
399 // Tuning of Pythia parameters aimed to get a resonable agreement
400 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
401 // b-bbar single inclusive and double differential distributions.
402 // This parameter settings are meant to work with Pb-Pb collisions
403 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
404 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
405 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
408 ReadString("SigmaProcess:factorMultFac = 1.");
410 ReadString("Beams:primordialKT = on");
411 ReadString("Beams:primordialKTsoft = 0.");
412 ReadString("Beams:primordialKThard = 2.035");
413 ReadString("Beams:halfScaleForKT = 0.");
414 ReadString("Beams:halfMassForKT = 0.");
416 ReadString("ParticleData:mbRun = 4.75");
418 case kPyBeautypPbMNR:
419 // Tuning of Pythia parameters aimed to get a resonable agreement
420 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
421 // b-bbar single inclusive and double differential distributions.
422 // This parameter settings are meant to work with p-Pb collisions
423 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
424 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
425 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
428 ReadString("SigmaProcess:factorMultFac = 1.");
430 ReadString("Beams:primordialKT = on");
431 ReadString("Beams:primordialKTsoft = 0.");
432 ReadString("Beams:primordialKThard = 1.6");
433 ReadString("Beams:halfScaleForKT = 0.");
434 ReadString("Beams:halfMassForKT = 0.");
436 ReadString("ParticleData:mbRun = 4.75");
439 // Tuning of Pythia parameters aimed to get a resonable agreement
440 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
441 // b-bbar single inclusive and double differential distributions.
442 // This parameter settings are meant to work with pp collisions
443 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
444 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
445 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
448 ReadString("SigmaProcess:factorMultFac = 1.");
450 ReadString("Beams:primordialKT = on");
451 ReadString("Beams:primordialKTsoft = 0.");
452 ReadString("Beams:primordialKThard = 1.0");
453 ReadString("Beams:halfScaleForKT = 0.");
454 ReadString("Beams:halfMassForKT = 0.");
456 ReadString("ParticleData:mbRun = 4.75");
458 case kPyBeautyppMNRwmi:
459 // Tuning of Pythia parameters aimed to get a resonable agreement
460 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
461 // b-bbar single inclusive and double differential distributions.
462 // This parameter settings are meant to work with pp collisions
463 // and with kCTEQ5L PDFs.
464 // Added multiple interactions according to ATLAS tune settings.
465 // To get a "reasonable" agreement with MNR results, events have to be
466 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
468 // To get a "perfect" agreement with MNR results, events have to be
469 // generated in four ptHard bins with the following relative
477 ReadString("SigmaProcess:factorMultFac = 1.");
479 ReadString("Beams:primordialKT = on");
480 ReadString("Beams:primordialKTsoft = 0.");
481 ReadString("Beams:primordialKThard = 1.0");
482 ReadString("Beams:halfScaleForKT = 0.");
483 ReadString("Beams:halfMassForKT = 0.");
485 ReadString("ParticleData:mbRun = 4.75");
489 //Inclusive production of W+/-
491 ReadString("WeakSingleBoson:ffbar2W = on");
492 // Initial/final parton shower on (Pythia default)
493 // With parton showers on we are generating "W inclusive process"
494 ReadString("PartonLevel:ISR = on");
495 ReadString("PartonLevel:FSR = on");
498 //Inclusive production of Z
500 ReadString("WeakSingleBoson:ffbar2gmZ = on");
501 //only Z included, not gamma
502 ReadString("WeakZ0:gmZmode = 2");
503 // Initial/final parton shower on (Pythia default)
504 // With parton showers on we are generating "Z inclusive process"
505 ReadString("PartonLevel:ISR = on");
506 ReadString("PartonLevel:FSR = on");
511 // SetMSTP(41,1); // all resonance decays switched on
512 Initialize(2212, 2212, fEcms);
515 void AliPythia8::SetNuclei(Int_t /*a1*/, Int_t /*a2*/)
517 // Treat protons as inside nuclei with mass numbers a1 and a2
518 // The MSTP array in the PYPARS common block is used to enable and
519 // select the nuclear structure functions.
520 // MSTP(52) : (D=1) choice of proton and nuclear structure-function library
521 // =1: internal PYTHIA acording to MSTP(51)
522 // =2: PDFLIB proton s.f., with MSTP(51) = 1000xNGROUP+NSET
523 // If the following mass number both not equal zero, nuclear corrections of the stf are used.
524 // MSTP(192) : Mass number of nucleus side 1
525 // MSTP(193) : Mass number of nucleus side 2
532 AliPythia8* AliPythia8::Instance()
534 // Set random number generator
538 fgAliPythia8 = new AliPythia8();
543 void AliPythia8::PrintParticles()
545 // Print list of particl properties
546 ReadString("Main:showAllParticleData");
549 void AliPythia8::ResetDecayTable()
551 // Set default values for pythia decay switches
553 // for (i = 1; i < 501; i++) SetMDCY(i,1,fDefMDCY[i]);
554 // for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
557 void AliPythia8::SetDecayTable()
559 // Set default values for pythia decay switches
562 // for (i = 1; i < 501; i++) fDefMDCY[i] = GetMDCY(i,1);
563 // for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
566 void AliPythia8::Pyclus(Int_t& njet)
568 // Call Pythia clustering algorithm
570 Bool_t ok = fClusterJet.analyze(Pythia8()->event, fYScale, fPtScale, fNJetMin, fNJetMax);
572 if (ok) njet = fClusterJet.size();
575 void AliPythia8::Pycell(Int_t& njet)
577 // Call Pythia jet reconstruction algorithm
579 Bool_t ok = fCellJet.analyze(Pythia8()->event, fMinEtJet, fRJet, fEtSeed);
581 if (ok) njet = fCellJet.size();
584 void AliPythia8::GetJet(Int_t i, Float_t& px, Float_t& py, Float_t& pz, Float_t& e)
587 Float_t et = fCellJet.eT(i);
588 px = et * TMath::Cos(fCellJet.phiWeighted(i));
589 py = et * TMath::Sin(fCellJet.phiWeighted(i));
590 pz = et * TMath::SinH(fCellJet.etaWeighted(i));
591 e = et * TMath::CosH(fCellJet.etaWeighted(i));
594 void AliPythia8::GenerateEvent()
596 // Generate one event
597 TPythia8::GenerateEvent();
600 void AliPythia8::GenerateMIEvent()
602 // New multiple interaction scenario
603 AliWarning("Not implemented. No event will be generated");
606 void AliPythia8::PrintStatistics()
608 // End of run statistics
609 TPythia8::PrintStatistics();
612 void AliPythia8::EventListing()
614 // End of run statistics
615 TPythia8::EventListing();
618 Int_t AliPythia8::ProcessCode()
620 // Returns the subprocess code for the current event
621 return Pythia8()->info.codeSub();
624 void AliPythia8::ConfigHeavyFlavor()
627 // Default configuration for Heavy Flavor production
631 ReadString("HardQCD:all = on");
633 // No multiple interactions
634 ReadString("PartonLevel:MI = off");
635 ReadString("MultipleInteractions:pTmin = 0.0");
636 ReadString("MultipleInteractions:pT0Ref = 0.0");
638 // Initial/final parton shower on (Pythia default)
639 ReadString("PartonLevel:ISR = on");
640 ReadString("PartonLevel:FSR = on");
643 ReadString("SigmaProcess:alphaSorder = 2");
646 ReadString("SigmaProcess:renormScale2 = 2");
647 ReadString("SigmaProcess:renormMultFac = 1.");
650 void AliPythia8::AtlasTuning()
653 // Configuration for the ATLAS tuning
654 ReadString(Form("PDF:LHAPDFset = %s", AliStructFuncType::PDFsetName(kCTEQ5L).Data()));
655 ReadString("PartonLevel:MI = on");
656 ReadString("MultipleInteractions:pTmin = 1.9");
657 ReadString("MultipleInteractions:pT0Ref = 1.8");
658 ReadString("MultipleInteractions:ecmRef = 1000.");
659 ReadString("MultipleInteractions:expPow = 0.16");
660 ReadString("MultipleInteractions:bProfile = 2");
661 ReadString("MultipleInteractions:coreFraction = 0.16");
662 ReadString("MultipleInteractions:coreRadius = 0.5");
663 // SetPARP(85,0.33); // Regulates gluon prod. mechanism
664 // SetPARP(86,0.66); // Regulates gluon prod. mechanism
665 ReadString("SigmaProcess:factorMultFac = 1.");
668 void AliPythia8::SetPtHardRange(Float_t ptmin, Float_t ptmax)
670 // Set the pt hard range
671 ReadString(Form("PhaseSpace:pTHatMin = %13.3f", ptmin));
672 ReadString(Form("PhaseSpace:pTHatMax = %13.3f", ptmax));
675 void AliPythia8::SetYHardRange(Float_t /*ymin*/, Float_t /*ymax*/)
677 // Set the y hard range
678 printf("YHardRange not implemented in Pythia8 !!!\n");
683 void AliPythia8::SetFragmentation(Int_t flag)
685 // Switch fragmentation on/off
687 ReadString("HadronLevel:Hadronize = on");
689 ReadString("HadronLevel:Hadronize = off");
693 void AliPythia8::SetInitialAndFinalStateRadiation(Int_t flag1, Int_t flag2)
695 // initial state radiation
697 ReadString("PartonLevel:ISR = on");
699 ReadString("PartonLevel:ISR = off");
701 // final state radiation
703 ReadString("PartonLevel:FSR = on");
705 ReadString("PartonLevel:FSR = off");
709 void AliPythia8::SetIntrinsicKt(Float_t kt)
711 ReadString("Beams:primordialKT = on");
712 ReadString("Beams:primordialKTsoft = 0.");
713 ReadString(Form("Beams:primordialKThard = %13.3f", kt));
714 ReadString("Beams:halfScaleForKT = 0.");
715 ReadString("Beams:halfMassForKT = 0.");
718 void AliPythia8::SwitchHFOff()
720 // Switch off heavy flavor
721 // Maximum number of quark flavours used in pdf
722 ReadString("PDFinProcess:nQuarkIn = 3");
723 // Maximum number of flavors that can be used in showers
724 ReadString("TimeShower:nGluonToQuark = 3");
725 ReadString("SpaceShower:nQuarkIn = 3");
730 void AliPythia8::SetPycellParameters(Float_t etaMax, Int_t nEta, Int_t nPhi,
731 Float_t thresh, Float_t etseed, Float_t minet, Float_t r)
733 // Set pycell parameters
734 fCellJet = Pythia8::CellJet( etaMax, nEta, nPhi, 2, 0, 0., 0., thresh);
740 void AliPythia8::ModifiedSplitting()
743 // We have to see how to implement this in Pythia8 !!!
745 // Modified splitting probability as a model for quenching
746 // SetPARJ(200, 0.8);
747 // SetMSTJ(41, 1); // QCD radiation only
748 // SetMSTJ(42, 2); // angular ordering
749 // SetMSTJ(44, 2); // option to run alpha_s
750 // SetMSTJ(47, 0); // No correction back to hard scattering element
751 // SetMSTJ(50, 0); // No coherence in first branching
752 // SetPARJ(82, 1.); // Cut off for parton showers
756 void AliPythia8::InitQuenching(Float_t /*cMin*/, Float_t /*cMax*/, Float_t /*k*/, Int_t /*iECMethod*/, Float_t /*zmax*/, Int_t /*ngmax*/)
760 AliWarning("Not implemented !");
763 void AliPythia8::SwitchHadronisationOff()
765 // Switch off hadronisation
766 ReadString("HadronLevel:Hadronize = off");
769 void AliPythia8::SwitchHadronisationOn()
771 // Switch on hadronisarion
772 ReadString("HadronLevel:Hadronize = on");
776 void AliPythia8::GetXandQ(Float_t& x1, Float_t& x2, Float_t& q)
778 // Get x1, x2 and Q for this event
780 q = Pythia8()->info.QFac();
781 x1 = Pythia8()->info.x1();
782 x2 = Pythia8()->info.x2();
786 Float_t AliPythia8::GetXSection()
788 // Get the total cross-section
789 return Pythia8()->info.sigmaGen();
792 Float_t AliPythia8::GetPtHard()
794 // Get the pT hard for this event
795 return Pythia8()->info.pTHat();
801 AliPythia8& AliPythia8::operator=(const AliPythia8& rhs)
803 // Assignment operator
808 void AliPythia8::Copy(TObject&) const
813 Fatal("Copy","Not implemented!\n");
819 void AliPythia8::SetNumberOfParticles(Int_t /*i*/)
821 AliWarning("Not implemented");
824 void AliPythia8::EditEventList(Int_t /*i*/)
826 AliWarning("Not implemented");
829 void AliPythia8::Pyquen(Double_t /*a*/, Int_t /*b*/, Double_t /*c*/)
831 AliWarning("Cannot be used with Pythia8");
834 void AliPythia8::HadronizeEvent()
836 // Needs access to HadronLevel ?
837 AliWarning("Not yet implemented");
840 void AliPythia8::GetQuenchingParameters(Double_t& /*xp*/, Double_t& /*yp*/, Double_t* /*z[4]*/)
842 AliWarning("Not yet implemented");
845 void AliPythia8::LoadEvent(AliStack* /*stack*/, Int_t /*flag*/, Int_t /*reHadr*/)
847 AliWarning("Not yet implemented");