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) 2020 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/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"
31 
32 namespace Pythia8 {
33 
34 //==========================================================================
35 
36 // The ProcessContainer class combines pointers to matrix element and
37 // phase space generator with general generation info.
38 
39 class ProcessContainer : public PhysicsBase {
40 
41 public:
42 
43  // Constructor.
44  ProcessContainer(SigmaProcess* sigmaProcessPtrIn = 0,
45  bool externalPtrIn = false, PhaseSpace* phaseSpacePtrIn = 0) :
46  sigmaProcessPtr(sigmaProcessPtrIn), externalPtr(externalPtrIn),
47  phaseSpacePtr(phaseSpacePtrIn), resDecaysPtr(), gammaKinPtr(),
48  matchInOut(), idRenameBeams(), setLifetime(), setQuarkMass(),
49  setLeptonMass(), idNewM(), mRecalculate(), mNewM(), isLHA(), isNonDiff(),
50  isResolved(), isDiffA(), isDiffB(), isDiffC(), isQCD3body(),
51  allowNegSig(), isSameSave(), increaseMaximum(), canVetoResDecay(),
52  lhaStrat(), lhaStratAbs(), processCode(), useStrictLHEFscales(),
53  isAsymLHA(), betazLHA(), newSigmaMx(), nTry(), nSel(), nAcc(),
54  nTryStat(), sigmaMx(), sigmaSgn(), sigmaSum(), sigma2Sum(), sigmaNeg(),
55  sigmaAvg(), sigmaFin(), deltaFin(), weightNow(), wtAccSum(),
56  beamAhasResGamma(), beamBhasResGamma(), beamHasResGamma(),
57  beamHasGamma(), beamAgammaMode(), beamBgammaMode(), gammaModeEvent(),
58  approximatedGammaFlux(), nTryRequested(), nSelRequested(),
59  nAccRequested(), sigmaTemp(), sigma2Temp(), normVar3() {}
60 
61  // Destructor. Do not destroy external sigmaProcessPtr.
62  ~ProcessContainer() {delete phaseSpacePtr;
63  if (!externalPtr) delete sigmaProcessPtr;}
64 
65  // Initialize phase space and counters.
66  bool init(bool isFirst, ResonanceDecays* resDecaysPtrIn,
67  SLHAinterface* slhaInterfacePtr, GammaKinematics* gammaKinPtrIn);
68 
69  // Store or replace Les Houches pointer.
70  void setLHAPtr( LHAupPtr lhaUpPtrIn, ParticleData* particleDataPtrIn = 0,
71  Settings* settingsPtrIn = 0, Rndm* rndmPtrIn = 0)
72  {lhaUpPtr = lhaUpPtrIn; setLifetime = 0;
73  if (settingsPtrIn && rndmPtrIn) {
74  rndmPtr = rndmPtrIn;
75  setLifetime = settingsPtrIn->mode("LesHouches:setLifetime");
76  }
77  if (particleDataPtrIn != 0) particleDataPtr = particleDataPtrIn;
78  if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
79  if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
80 
81  // Update the CM energy of the event.
82  void newECM(double eCM) {phaseSpacePtr->newECM(eCM);}
83 
84  // Generate a trial event; accepted or not.
85  bool trialProcess();
86 
87  // Pick flavours and colour flow of process.
88  bool constructState();
89 
90  // Give the hard subprocess (with option for a second hard subprocess).
91  bool constructProcess( Event& process, bool isHardest = true);
92 
93  // Give the Les Houches decay chain for external resonance.
94  bool constructDecays( Event& process);
95 
96  // Do resonance decays.
97  bool decayResonances( Event& process);
98 
99  // Accumulate statistics after user veto.
100  void accumulate();
101 
102  // Reset statistics on events generated so far.
103  void reset();
104 
105  // Set whether (photon) beam is resolved or unresolved.
106  // Method propagates the choice of photon process type to beam pointers.
107  void setBeamModes(bool setVMD = false, bool isSampled = true);
108 
109  // Process name and code, and the number of final-state particles.
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);}
117 
118  // Member functions for info on generation process.
119  bool newSigmaMax() const {return newSigmaMx;}
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;}
131 
132  // Some kinematics quantities.
133  int id1() const {return sigmaProcessPtr->id(1);}
134  int id2() const {return sigmaProcessPtr->id(2);}
135  double x1() const {return phaseSpacePtr->x1();}
136  double x2() const {return phaseSpacePtr->x2();}
137  double Q2Fac() const {return sigmaProcessPtr->Q2Fac();}
138  double mHat() const {return sqrtpos(phaseSpacePtr->sHat());}
139  double pTHat() const {return phaseSpacePtr->pTHat();}
140 
141  // Tell whether container is for Les Houches events.
142  bool isLHAContainer() const {return isLHA;}
143  int lhaStrategy() const {return lhaStrat;}
144 
145  // Info on Les Houches events.
146  int codeLHASize() const {return codeLHA.size();}
147  int subCodeLHA(int i) const {return codeLHA[i];}
148  long nTriedLHA(int i) const {return nTryLHA[i];}
149  long nSelectedLHA(int i) const {return nSelLHA[i];}
150  long nAcceptedLHA(int i) const {return nAccLHA[i];}
151 
152  // When two hard processes set or get info whether process is matched.
153  void isSame( bool isSameIn) { isSameSave = isSameIn;}
154  bool isSame() const {return isSameSave;}
155 
156 private:
157 
158  // Constants: could only be changed in the code itself.
159  static const int N12SAMPLE, N3SAMPLE;
160 
161  // Pointer to the subprocess matrix element. Mark if external.
162  SigmaProcess* sigmaProcessPtr;
163  bool externalPtr;
164 
165  // Pointer to the phase space generator.
166  PhaseSpace* phaseSpacePtr;
167 
168  // Pointer to ResonanceDecays object for sequential resonance decays.
169  ResonanceDecays* resDecaysPtr;
170 
171  // Pointer to LHAup for generating external events.
172  LHAupPtr lhaUpPtr;
173 
174  // Pointer to the phase space generator of photons from leptons.
175  GammaKinematics* gammaKinPtr;
176 
177  // Possibility to modify Les Houches input.
178  bool matchInOut;
179  int idRenameBeams, setLifetime, setQuarkMass, setLeptonMass, idNewM[9];
180  double mRecalculate, mNewM[9];
181 
182  // Info on process.
183  bool isLHA, isNonDiff, isResolved, isDiffA, isDiffB, isDiffC, isQCD3body,
184  allowNegSig, isSameSave, increaseMaximum, canVetoResDecay;
185  int lhaStrat, lhaStratAbs, processCode;
186  bool useStrictLHEFscales;
187 
188  // Boost Les Houches events to CM frame (when originally asymmetric).
189  bool isAsymLHA;
190  double betazLHA;
191 
192  // Statistics on generation process. (Long integers just in case.)
193  bool newSigmaMx;
194  long nTry, nSel, nAcc, nTryStat;
195  double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
196  sigmaFin, deltaFin, weightNow, wtAccSum;
197 
198  // Flags to store whether beam has a (un)resolved photon.
199  bool beamAhasResGamma, beamBhasResGamma, beamHasResGamma, beamHasGamma;
200  int beamAgammaMode, beamBgammaMode, gammaModeEvent;
201 
202  // Use approximated photon flux for process sampling.
203  bool approximatedGammaFlux;
204 
205  // Statistics for Les Houches event classification.
206  vector<int> codeLHA;
207  vector<long> nTryLHA, nSelLHA, nAccLHA;
208 
209  // More fine-grained event counting.
210  long nTryRequested, nSelRequested, nAccRequested;
211 
212  // Temporary summand for handling (weighted) events when vetoes are applied.
213  double sigmaTemp, sigma2Temp, normVar3;
214 
215  // Estimate integrated cross section and its uncertainty.
216  void sigmaDelta();
217 
218 };
219 
220 //==========================================================================
221 
222 // The SetupContainers class turns the list of user-requested processes
223 // into a vector of ProcessContainer objects, each with a process.
224 
225 class SetupContainers {
226 
227 public:
228 
229  // Constructor.
230  SetupContainers() : nVecA(), nVecB() {}
231 
232  // Initialization assuming all necessary data already read.
233  bool init(vector<ProcessContainer*>& containerPtrs, Info* infoPtr);
234 
235  // Initialization of a second hard process.
236  bool init2(vector<ProcessContainer*>& container2Ptrs, Info* infoPtr);
237 
238 private:
239 
240  // Methods to check that outgoing SUSY particles are allowed ones.
241  void setupIdVecs( Settings& settings);
242  bool allowIdVals( int idCheck1, int idCheck2);
243 
244  // Arrays of allowed outgoing SUSY particles and their lengths.
245  vector<int> idVecA, idVecB;
246  int nVecA, nVecB;
247 
248  // Helper class to setup onia production.
249  SigmaOniaSetup charmonium, bottomonium;
250 
251 };
252 
253 //==========================================================================
254 
255 } // end namespace Pythia8
256 
257 #endif // Pythia8_ProcessContainer_H
Definition: AgUStep.h:26