10 #ifndef Pythia8_ProcessContainer_H 11 #define Pythia8_ProcessContainer_H 13 #include "Pythia8/Basics.h" 14 #include "Pythia8/BeamParticle.h" 15 #include "Pythia8/Event.h" 16 #include "Pythia8/Info.h" 17 #include "Pythia8/ParticleData.h" 18 #include "Pythia8/PartonDistributions.h" 19 #include "Pythia8/PhaseSpace.h" 20 #include "Pythia8/PhysicsBase.h" 21 #include "Pythia8/PythiaStdlib.h" 22 #include "Pythia8/ResonanceDecays.h" 23 #include "Pythia8/Settings.h" 24 #include "Pythia8/SigmaProcess.h" 25 #include "Pythia8/SigmaTotal.h" 26 #include "Pythia8/SigmaOnia.h" 27 #include "Pythia8/StandardModel.h" 28 #include "Pythia8/SusyCouplings.h" 29 #include "Pythia8/SLHAinterface.h" 30 #include "Pythia8/UserHooks.h" 45 PhaseSpacePtr phaseSpacePtrIn = 0) :
46 sigmaProcessPtr(sigmaProcessPtrIn),
47 phaseSpacePtr(phaseSpacePtrIn), resDecaysPtr(), gammaKinPtr(),
48 matchInOut(), idRenameBeams(), setLifetime(), setQuarkMass(),
49 setLeptonMass(), idNewM(), smearHadronMass(), idSmearHadIn(),
50 idSmearHadrons(), mRecalculate(), mNewM(), isLHA(), isNonDiff(),
51 isResolved(), isDiffA(), isDiffB(), isDiffC(), isQCD3body(),
52 allowNegSig(), isSameSave(), increaseMaximum(), canVetoResDecay(),
53 lhaStrat(), lhaStratAbs(), processCode(), useStrictLHEFscales(),
54 isAsymLHA(), betazLHA(), newSigmaMx(), nTry(), nSel(), nAcc(),
55 nTryStat(), sigmaMx(), sigmaSgn(), sigmaSum(), sigma2Sum(), sigmaNeg(),
56 sigmaAvg(), sigmaFin(), deltaFin(), weightNow(), wtAccSum(),
57 beamAhasResGamma(), beamBhasResGamma(), beamHasResGamma(),
58 beamHasGamma(), beamAgammaMode(), beamBgammaMode(), gammaModeEvent(),
59 approximatedGammaFlux(), doMerging(), nTryRequested(), nSelRequested(),
60 nAccRequested(), sigmaTemp(), sigma2Temp(), normVar3() {}
69 {lhaUpPtr = lhaUpPtrIn; setLifetime = 0;
70 if (settingsPtrIn && rndmPtrIn) {
72 setLifetime = settingsPtrIn->mode(
"LesHouches:setLifetime");
75 if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
76 if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
82 void newECM(
double eCM) {phaseSpacePtr->newECM(eCM);}
107 void setBeamModes(
bool setVMD =
false,
bool isSampled =
true);
110 string name()
const {
return sigmaProcessPtr->name();}
111 int code()
const {
return sigmaProcessPtr->code();}
112 int nFinal()
const {
return sigmaProcessPtr->nFinal();}
113 bool isSUSY()
const {
return sigmaProcessPtr->isSUSY();}
114 bool isNonDiffractive()
const {
return isNonDiff;}
115 bool isSoftQCD()
const {
return (code() > 100 && code() < 107);}
116 bool isElastic()
const {
return (code() == 102);}
120 double sigmaMax()
const {
return sigmaMx;}
121 long nTried()
const {
return nTry;}
122 long nSelected()
const {
return nSel;}
123 long nAccepted()
const {
return nAcc;}
124 double weightSum()
const {
return wtAccSum;}
125 double sigmaSelMC(
bool doAccumulate =
true)
126 {
if (nTry > nTryStat && doAccumulate) sigmaDelta();
return sigmaAvg;}
127 double sigmaMC(
bool doAccumulate =
true)
128 {
if (nTry > nTryStat && doAccumulate) sigmaDelta();
return sigmaFin;}
129 double deltaMC(
bool doAccumulate =
true)
130 {
if (nTry > nTryStat && doAccumulate) sigmaDelta();
return deltaFin;}
137 int id1()
const {
return sigmaProcessPtr->id(1);}
138 int id2()
const {
return sigmaProcessPtr->id(2);}
139 double x1()
const {
return phaseSpacePtr->x1();}
140 double x2()
const {
return phaseSpacePtr->x2();}
141 double Q2Fac()
const {
return sigmaProcessPtr->Q2Fac();}
142 double mHat()
const {
return sqrtpos(phaseSpacePtr->sHat());}
143 double pTHat()
const {
return phaseSpacePtr->pTHat();}
147 int lhaStrategy()
const {
return lhaStrat;}
151 int subCodeLHA(
int i)
const {
return codeLHA[i];}
152 long nTriedLHA(
int i)
const {
return nTryLHA[i];}
153 long nSelectedLHA(
int i)
const {
return nSelLHA[i];}
154 long nAcceptedLHA(
int i)
const {
return nAccLHA[i];}
157 void isSame(
bool isSameIn) { isSameSave = isSameIn;}
158 bool isSame()
const {
return isSameSave;}
163 static const int N12SAMPLE, N3SAMPLE;
166 SigmaProcessPtr sigmaProcessPtr;
169 PhaseSpacePtr phaseSpacePtr;
182 int idRenameBeams, setLifetime, setQuarkMass, setLeptonMass, idNewM[9],
184 vector<int> idSmearHadIn, idSmearHadrons;
185 double mRecalculate, mNewM[9];
188 bool isLHA, isNonDiff, isResolved, isDiffA, isDiffB, isDiffC, isQCD3body,
189 allowNegSig, isSameSave, increaseMaximum, canVetoResDecay;
190 int lhaStrat, lhaStratAbs, processCode;
191 bool useStrictLHEFscales;
199 long nTry, nSel, nAcc, nTryStat;
200 double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
201 sigmaFin, deltaFin, weightNow, wtAccSum;
204 bool beamAhasResGamma, beamBhasResGamma, beamHasResGamma, beamHasGamma;
205 int beamAgammaMode, beamBgammaMode, gammaModeEvent;
208 bool approximatedGammaFlux;
215 vector<long> nTryLHA, nSelLHA, nAccLHA;
218 long nTryRequested, nSelRequested, nAccRequested;
221 double sigmaTemp, sigma2Temp, normVar3;
244 bool init2(vector<ProcessContainer*>& container2Ptrs,
Info* infoPtr);
249 void setupIdVecs(
Settings& settings);
250 bool allowIdVals(
int idCheck1,
int idCheck2);
253 vector<int> idVecA, idVecB;
ProcessContainer(SigmaProcessPtr sigmaProcessPtrIn=0, PhaseSpacePtr phaseSpacePtrIn=0)
Constructor.
Definition: ProcessContainer.h:44
void accumulate()
Accumulate statistics after user veto.
Definition: ProcessContainer.cc:450
double sigmaMaxSwitch()
New value for switched beam identity or energy (for SoftQCD processes).
Definition: ProcessContainer.h:133
Definition: PhysicsBase.h:27
void newECM(double eCM)
Update the CM energy of the event.
Definition: ProcessContainer.h:82
The Event class holds all info on the generated event.
Definition: Event.h:408
bool trialProcess()
Generate a trial event; accepted or not.
Definition: ProcessContainer.cc:268
Definition: ProcessContainer.h:39
Definition: SLHAinterface.h:27
bool constructProcess(Event &process, bool isHardest=true)
Give the hard subprocess (with option for a second hard subprocess).
Definition: ProcessContainer.cc:533
void isSame(bool isSameIn)
When two hard processes set or get info whether process is matched.
Definition: ProcessContainer.h:157
bool init(bool isFirst, ResonanceDecays *resDecaysPtrIn, SLHAinterface *slhaInterfacePtr, GammaKinematics *gammaKinPtrIn)
Initialize phase space and counters.
Definition: ProcessContainer.cc:48
Rndm * rndmPtr
Pointer to the random number generator.
Definition: PhysicsBase.h:93
SetupContainers()
Constructor.
Definition: ProcessContainer.h:238
A helper class used to setup the onia processes.
Definition: SigmaOnia.h:63
int id1() const
Some kinematics quantities.
Definition: ProcessContainer.h:137
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: PhysicsBase.h:84
bool newSigmaMax() const
Member functions for info on generation process.
Definition: ProcessContainer.h:119
int codeLHASize() const
Info on Les Houches events.
Definition: ProcessContainer.h:150
void setLHAPtr(LHAupPtr lhaUpPtrIn, ParticleData *particleDataPtrIn=0, Settings *settingsPtrIn=0, Rndm *rndmPtrIn=0)
Store or replace Les Houches pointer.
Definition: ProcessContainer.h:67
void reset()
Reset statistics on events generated so far.
Definition: ProcessContainer.cc:1385
Class to sample the virtuality and transverse momentum of emitted photons.
Definition: GammaKinematics.h:23
string name() const
Process name and code, and the number of final-state particles.
Definition: ProcessContainer.h:110
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
void setBeamModes(bool setVMD=false, bool isSampled=true)
Set the photon modes according to the process.
Definition: ProcessContainer.cc:496
Info * infoPtr
Definition: PhysicsBase.h:78
bool decayResonances(Event &process)
Do resonance decays.
Definition: ProcessContainer.cc:1325
Definition: ProcessContainer.h:233
Definition: ResonanceDecays.h:28
bool constructDecays(Event &process)
Give the Les Houches decay chain for external resonance.
Definition: ProcessContainer.cc:1201
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
void updateBeamIDs()
Switch to new beam particle identities; for similar hadrons only.
Definition: ProcessContainer.h:79
bool constructState()
Pick flavours and colour flow of process.
Definition: ProcessContainer.cc:476
Definition: Settings.h:196
bool isLHAContainer() const
Tell whether container is for Les Houches events.
Definition: ProcessContainer.h:146
double sqrtpos(const double &x)
Avoid problem with negative square root argument (from roundoff).
Definition: PythiaStdlib.h:191