StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ProcessContainer.h
1 // ProcessContainer.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
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.
9 
10 #ifndef Pythia8_ProcessContainer_H
11 #define Pythia8_ProcessContainer_H
12 
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/PythiaStdlib.h"
21 #include "Pythia8/ResonanceDecays.h"
22 #include "Pythia8/Settings.h"
23 #include "Pythia8/SigmaProcess.h"
24 #include "Pythia8/SigmaTotal.h"
25 #include "Pythia8/SigmaOnia.h"
26 #include "Pythia8/StandardModel.h"
27 #include "Pythia8/SusyCouplings.h"
28 #include "Pythia8/SLHAinterface.h"
29 #include "Pythia8/UserHooks.h"
30 
31 namespace Pythia8 {
32 
33 //==========================================================================
34 
35 // The ProcessContainer class combines pointers to matrix element and
36 // phase space generator with general generation info.
37 
38 class ProcessContainer {
39 
40 public:
41 
42  // Constructor.
43  ProcessContainer(SigmaProcess* sigmaProcessPtrIn = 0,
44  bool externalPtrIn = false, PhaseSpace* phaseSpacePtrIn = 0) :
45  sigmaProcessPtr(sigmaProcessPtrIn),
46  externalPtr(externalPtrIn), phaseSpacePtr(phaseSpacePtrIn) {}
47 
48  // Destructor. Do not destroy external sigmaProcessPtr.
49  ~ProcessContainer() {delete phaseSpacePtr;
50  if (!externalPtr) delete sigmaProcessPtr;}
51 
52  // Initialize phase space and counters.
53  bool init(bool isFirst, Info* infoPtrIn, Settings& settings,
54  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn, BeamParticle* beamAPtr,
55  BeamParticle* beamBPtr, Couplings* couplings, SigmaTotal* sigmaTotPtrIn,
56  ResonanceDecays* resDecaysPtrIn, SLHAinterface* slhaInterfacePtr,
57  UserHooks* userHooksPtr, GammaKinematics* gammaKinPtrIn);
58 
59  // Store or replace Les Houches pointer.
60  void setLHAPtr( LHAup* lhaUpPtrIn, ParticleData* particleDataPtrIn = 0,
61  Settings* settingsPtrIn = 0, Rndm* rndmPtrIn = 0)
62  {lhaUpPtr = lhaUpPtrIn; setLifetime = 0;
63  if (settingsPtrIn && rndmPtrIn) {
64  rndmPtr = rndmPtrIn;
65  setLifetime = settingsPtrIn->mode("LesHouches:setLifetime");
66  }
67  if (particleDataPtrIn != 0) particleDataPtr = particleDataPtrIn;
68  if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
69  if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
70 
71  // Update the CM energy of the event.
72  void newECM(double eCM) {phaseSpacePtr->newECM(eCM);}
73 
74  // Generate a trial event; accepted or not.
75  bool trialProcess();
76 
77  // Pick flavours and colour flow of process.
78  bool constructState();
79 
80  // Give the hard subprocess (with option for a second hard subprocess).
81  bool constructProcess( Event& process, bool isHardest = true);
82 
83  // Give the Les Houches decay chain for external resonance.
84  bool constructDecays( Event& process);
85 
86  // Do resonance decays.
87  bool decayResonances( Event& process);
88 
89  // Accumulate statistics after user veto.
90  void accumulate();
91 
92  // Reset statistics on events generated so far.
93  void reset();
94 
95  // Set whether (photon) beam is resolved or unresolved.
96  // Method propagates the choice of photon process type to beam pointers.
97  void setBeamModes();
98 
99  // If resolved photons, then choose VMD states. Method propagates
100  // states to beam pointers and info pointer.
101  pair<int, int> chooseVMDstates(int idA, int idB);
102 
103  // Process name and code, and the number of final-state particles.
104  string name() const {return sigmaProcessPtr->name();}
105  int code() const {return sigmaProcessPtr->code();}
106  int nFinal() const {return sigmaProcessPtr->nFinal();}
107  bool isSUSY() const {return sigmaProcessPtr->isSUSY();}
108  bool isNonDiffractive() const {return isNonDiff;}
109  bool isSoftQCD() const {return (code() > 100 && code() < 107);}
110 
111  // Member functions for info on generation process.
112  bool newSigmaMax() const {return newSigmaMx;}
113  double sigmaMax() const {return sigmaMx;}
114  long nTried() const {return nTry;}
115  long nSelected() const {return nSel;}
116  long nAccepted() const {return nAcc;}
117  double weightSum() const {return wtAccSum;}
118  double sigmaSelMC( bool doAccumulate = true)
119  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaAvg;}
120  double sigmaMC( bool doAccumulate = true)
121  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaFin;}
122  double deltaMC( bool doAccumulate = true)
123  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return deltaFin;}
124 
125  // Some kinematics quantities.
126  int id1() const {return sigmaProcessPtr->id(1);}
127  int id2() const {return sigmaProcessPtr->id(2);}
128  double x1() const {return phaseSpacePtr->x1();}
129  double x2() const {return phaseSpacePtr->x2();}
130  double Q2Fac() const {return sigmaProcessPtr->Q2Fac();}
131  double mHat() const {return sqrtpos(phaseSpacePtr->sHat());}
132  double pTHat() const {return phaseSpacePtr->pTHat();}
133 
134  // Tell whether container is for Les Houches events.
135  bool isLHAContainer() const {return isLHA;}
136  int lhaStrategy() const {return lhaStrat;}
137 
138  // Info on Les Houches events.
139  int codeLHASize() const {return codeLHA.size();}
140  int subCodeLHA(int i) const {return codeLHA[i];}
141  long nTriedLHA(int i) const {return nTryLHA[i];}
142  long nSelectedLHA(int i) const {return nSelLHA[i];}
143  long nAcceptedLHA(int i) const {return nAccLHA[i];}
144 
145  // When two hard processes set or get info whether process is matched.
146  void isSame( bool isSameIn) { isSameSave = isSameIn;}
147  bool isSame() const {return isSameSave;}
148 
149 private:
150 
151  // Constants: could only be changed in the code itself.
152  static const int N12SAMPLE, N3SAMPLE;
153 
154  // Pointer to the subprocess matrix element. Mark if external.
155  SigmaProcess* sigmaProcessPtr;
156  bool externalPtr;
157 
158  // Pointer to the phase space generator.
159  PhaseSpace* phaseSpacePtr;
160 
161  // Pointer to various information on the generation.
162  Info* infoPtr;
163 
164  // Pointer to the particle data table.
165  ParticleData* particleDataPtr;
166 
167  // Pointer to the random number generator.
168  Rndm* rndmPtr;
169 
170  // Pointer to ResonanceDecays object for sequential resonance decays.
171  ResonanceDecays* resDecaysPtr;
172 
173  // Pointer to SigmaTotal.
174  SigmaTotal* sigmaTotPtr;
175 
176  // Pointer to userHooks object for user interaction with program.
177  UserHooks* userHooksPtr;
178 
179  // Pointer to LHAup for generating external events.
180  LHAup* lhaUpPtr;
181 
182  // Pointers to the two incoming beams.
183  BeamParticle* beamAPtr;
184  BeamParticle* beamBPtr;
185 
186  // Pointer to the phase space generator of photons from leptons.
187  GammaKinematics* gammaKinPtr;
188 
189  // Possibility to modify Les Houches input.
190  bool matchInOut;
191  int idRenameBeams, setLifetime, setQuarkMass, setLeptonMass, idNewM[9];
192  double mRecalculate, mNewM[9];
193 
194  // Info on process.
195  bool isLHA, isNonDiff, isResolved, isDiffA, isDiffB, isDiffC,
196  isDiff, isSingleDiff, isDoubleDiff, isCentralDiff, isQCD3body,
197  allowNegSig, isSameSave, increaseMaximum, canVetoResDecay;
198  int lhaStrat, lhaStratAbs;
199  bool useStrictLHEFscales;
200 
201  // Statistics on generation process. (Long integers just in case.)
202  bool newSigmaMx;
203  long nTry, nSel, nAcc, nTryStat;
204  double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
205  sigmaFin, deltaFin, weightNow, wtAccSum;
206 
207  // Flags to store whether beam has a (un)resolved photon.
208  bool beamAhasResGamma, beamBhasResGamma, beamHasResGamma, beamHasGamma;
209  int beamAgammaMode, beamBgammaMode, gammaModeEvent;
210 
211  // Use external photon flux.
212  bool externalFlux;
213 
214  // Statistics for Les Houches event classification.
215  vector<int> codeLHA;
216  vector<long> nTryLHA, nSelLHA, nAccLHA;
217 
218  // More fine-grained event counting.
219  long nTryRequested, nSelRequested, nAccRequested;
220 
221  // Temporary summand for handling (weighted) events when vetoes are applied.
222  double sigmaTemp, sigma2Temp;
223 
224  // Estimate integrated cross section and its uncertainty.
225  void sigmaDelta();
226 
227 };
228 
229 //==========================================================================
230 
231 // The SetupContainers class turns the list of user-requested processes
232 // into a vector of ProcessContainer objects, each with a process.
233 
234 class SetupContainers {
235 
236 public:
237 
238  // Constructor.
239  SetupContainers() {}
240 
241  // Initialization assuming all necessary data already read.
242  bool init(vector<ProcessContainer*>& containerPtrs, Info* infoPtr,
243  Settings& settings, ParticleData* particleDataPtr, Couplings* couplings);
244 
245  // Initialization of a second hard process.
246  bool init2(vector<ProcessContainer*>& container2Ptrs, Settings& settings);
247 
248 private:
249 
250  // Methods to check that outgoing SUSY particles are allowed ones.
251  void setupIdVecs( Settings& settings);
252  bool allowIdVals( int idCheck1, int idCheck2);
253 
254  // Arrays of allowed outgoing SUSY particles and their lengths.
255  vector<int> idVecA, idVecB;
256  int nVecA, nVecB;
257 
258  // Helper class to setup onia production.
259  SigmaOniaSetup charmonium, bottomonium;
260 
261 };
262 
263 //==========================================================================
264 
265 } // end namespace Pythia8
266 
267 #endif // Pythia8_ProcessContainer_H
Definition: AgUStep.h:26