1 // ProcessContainer.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2013 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
6 // This file contains the collected machinery of a process.
7 // ProcessContainer: contains information on a particular process.
8 // SetupContainers: administrates the selection/creation of processes.
10 #ifndef Pythia8_ProcessContainer_H
11 #define Pythia8_ProcessContainer_H
14 #include "BeamParticle.h"
17 #include "ParticleData.h"
18 #include "PartonDistributions.h"
19 #include "PhaseSpace.h"
20 #include "PythiaStdlib.h"
21 #include "ResonanceDecays.h"
23 #include "SigmaProcess.h"
24 #include "SigmaTotal.h"
25 #include "StandardModel.h"
26 #include "SusyCouplings.h"
27 #include "SusyLesHouches.h"
28 #include "UserHooks.h"
32 //==========================================================================
34 // The ProcessContainer class combines pointers to matrix element and
35 // phase space generator with general generation info.
37 class ProcessContainer {
42 ProcessContainer(SigmaProcess* sigmaProcessPtrIn = 0,
43 bool externalPtrIn = false) : sigmaProcessPtr(sigmaProcessPtrIn),
44 externalPtr(externalPtrIn), phaseSpacePtr(0) {}
46 // Destructor. Do not destroy external sigmaProcessPtr.
47 ~ProcessContainer() {delete phaseSpacePtr;
48 if (!externalPtr) delete sigmaProcessPtr;}
50 // Initialize phase space and counters.
51 bool init(bool isFirst, Info* infoPtrIn, Settings& settings,
52 ParticleData* particleDataPtrIn, Rndm* rndmPtrIn, BeamParticle* beamAPtr,
53 BeamParticle* beamBPtr, Couplings* couplings, SigmaTotal* sigmaTotPtr,
54 ResonanceDecays* resDecaysPtrIn, SusyLesHouches* slhaPtr,
55 UserHooks* userHooksPtr);
57 // Store or replace Les Houches pointer.
58 void setLHAPtr( LHAup* lhaUpPtrIn, ParticleData* particleDataPtrIn = 0)
59 {lhaUpPtr = lhaUpPtrIn;
60 if (particleDataPtrIn != 0) particleDataPtr = particleDataPtrIn;
61 if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
62 if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
64 // Update the CM energy of the event.
65 void newECM(double eCM) {phaseSpacePtr->newECM(eCM);}
67 // Generate a trial event; accepted or not.
70 // Pick flavours and colour flow of process.
71 bool constructState();
73 // Give the hard subprocess (with option for a second hard subprocess).
74 bool constructProcess( Event& process, bool isHardest = true);
76 // Give the Les Houches decay chain for external resonance.
77 bool constructDecays( Event& process);
79 // Do resonance decays.
80 bool decayResonances( Event& process);
82 // Accumulate statistics after user veto.
85 // Reset statistics on events generated so far.
88 // Process name and code, and the number of final-state particles.
89 string name() const {return sigmaProcessPtr->name();}
90 int code() const {return sigmaProcessPtr->code();}
91 int nFinal() const {return sigmaProcessPtr->nFinal();}
92 bool isSUSY() const {return sigmaProcessPtr->isSUSY();}
94 // Member functions for info on generation process.
95 bool newSigmaMax() const {return newSigmaMx;}
96 double sigmaMax() const {return sigmaMx;}
97 long nTried() const {return nTry;}
98 long nSelected() const {return nSel;}
99 long nAccepted() const {return nAcc;}
100 double weightSum() const {return wtAccSum;}
101 double sigmaSelMC() {if (nTry > nTryStat) sigmaDelta(); return sigmaAvg;}
102 double sigmaMC() {if (nTry > nTryStat) sigmaDelta(); return sigmaFin;}
103 double deltaMC() {if (nTry > nTryStat) sigmaDelta(); return deltaFin;}
105 // Some kinematics quantities.
106 int id1() const {return sigmaProcessPtr->id(1);}
107 int id2() const {return sigmaProcessPtr->id(2);}
108 double x1() const {return phaseSpacePtr->x1();}
109 double x2() const {return phaseSpacePtr->x2();}
110 double Q2Fac() const {return sigmaProcessPtr->Q2Fac();}
111 double mHat() const {return sqrtpos(phaseSpacePtr->sHat());}
112 double pTHat() const {return phaseSpacePtr->pTHat();}
114 // Tell whether container is for Les Houches events.
115 bool isLHAContainer() const {return isLHA;}
116 int lhaStrategy() const {return lhaStrat;}
118 // Info on Les Houches events.
119 int codeLHASize() const {return codeLHA.size();}
120 int subCodeLHA(int i) const {return codeLHA[i];}
121 long nTriedLHA(int i) const {return nTryLHA[i];}
122 long nSelectedLHA(int i) const {return nSelLHA[i];}
123 long nAcceptedLHA(int i) const {return nAccLHA[i];}
125 // When two hard processes set or get info whether process is matched.
126 void isSame( bool isSameIn) { isSameSave = isSameIn;}
127 bool isSame() const {return isSameSave;}
131 // Constants: could only be changed in the code itself.
132 static const int N12SAMPLE, N3SAMPLE;
134 // Pointer to the subprocess matrix element. Mark if external.
135 SigmaProcess* sigmaProcessPtr;
138 // Pointer to the phase space generator.
139 PhaseSpace* phaseSpacePtr;
141 // Pointer to various information on the generation.
144 // Pointer to the particle data table.
145 ParticleData* particleDataPtr;
147 // Pointer to the random number generator.
150 // Pointer to ResonanceDecays object for sequential resonance decays.
151 ResonanceDecays* resDecaysPtr;
153 // Pointer to userHooks object for user interaction with program.
154 UserHooks* userHooksPtr;
156 // Pointer to LHAup for generating external events.
160 bool isLHA, isMinBias, isResolved, isDiffA, isDiffB, isDiffC, isQCD3body,
161 allowNegSig, hasOctetOnium, isSameSave, increaseMaximum,
163 int lhaStrat, lhaStratAbs;
165 // Statistics on generation process. (Long integers just in case.)
167 long nTry, nSel, nAcc, nTryStat;
168 double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
169 sigmaFin, deltaFin, weightNow, wtAccSum;
171 // Statistics for Les Houches event classification.
173 vector<long> nTryLHA, nSelLHA, nAccLHA;
174 long nTryLast, nSelLast;
176 // Estimate integrated cross section and its uncertainty.
181 //==========================================================================
183 // The SetupContainers class turns the list of user-requested processes
184 // into a vector of ProcessContainer objects, each with a process.
186 class SetupContainers {
193 // Initialization assuming all necessary data already read.
194 bool init(vector<ProcessContainer*>& containerPtrs, Settings& settings,
195 ParticleData* particleDataPtr, Couplings* couplings);
197 // Initialization of a second hard process.
198 bool init2(vector<ProcessContainer*>& container2Ptrs, Settings& settings);
202 //==========================================================================
204 } // end namespace Pythia8
206 #endif // Pythia8_ProcessContainer_H