PYTHIA  8.313
VinciaQED.h
1 // VinciaQED.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 Peter Skands, 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 QED antenna-shower class and auxiliary
7 // classes. Main author is Rob Verheyen.
8 
9 #ifndef Pythia8_VinciaQED_H
10 #define Pythia8_VinciaQED_H
11 
12 // PYTHIA 8 headers.
13 #include "Pythia8/BeamParticle.h"
14 #include "Pythia8/Event.h"
15 #include "Pythia8/StandardModel.h"
16 #include "Pythia8/PartonSystems.h"
17 
18 // VINCIA headers.
19 #include "Pythia8/VinciaCommon.h"
20 #include "Pythia8/VinciaWeights.h"
21 
22 namespace Pythia8 {
23 
24 //==========================================================================
25 
26 // Class for QED emissions.
27 
29 
30 public:
31 
32  // Friends for internal private members.
33  friend class QEDemitSystem;
34 
35  // Constuctor.
36  QEDemitElemental() : rndmPtr(nullptr), partonSystemsPtr(nullptr), q2Sav(0.),
37  zetaSav(0.), phiSav(0.), sxjSav(0.), syjSav(0.), alpha(0.), c(0.),
38  hasTrial(false), x(0), y(0), idx(0), idy(0), mx2(0.), my2(0.),
39  ex(0.), ey(0.), m2Ant(0.), sAnt(0.), QQ(0.), isII(false), isIF(false),
40  isFF(false), isRF(false), isIA(true), isDip(false), shh(0.),
41  isInitPtr(false), isInit(false), verbose(1) {;}
42 
43  // Initialize the pointers.
44  void initPtr(Rndm* rndmPtrIn, PartonSystems* partonSystemsPtrIn);
45  // Initialize.
46  void init(Event &event, int xIn, int yIn, double shhIn,
47  double verboseIn);
48  // Initialize.
49  void init(Event &event, int xIn, vector<int> iRecoilIn, double shhIn,
50  double verboseIn);
51  // Generate a trial point.
52  double generateTrial(Event &event, double q2Start, double q2Low,
53  double alphaIn, double cIn);
54 
55 private:
56 
57  // Random pointer.
58  Rndm* rndmPtr{};
59 
60  // Parton system pointer.
61  PartonSystems* partonSystemsPtr{};
62 
63  // Trial variables.
64  double q2Sav, zetaSav, phiSav, sxjSav, syjSav, alpha, c;
65  bool hasTrial;
66 
67  // Particle indices.
68  int x, y;
69  // Recoiler indices.
70  vector<int> iRecoil;
71  // IDs.
72  int idx, idy;
73  // Number of spin states.
74  int spinTypex, spinTypey;
75  // Particle masses.
76  double mx2, my2;
77  // Particle energies.
78  double ex, ey;
79  // Antenna invariant mass.
80  double m2Ant;
81  // Antenna dot product.
82  double sAnt;
83  // The negative of the product of charges.
84  double QQ;
85 
86  // Type switches.
87  bool isII, isIF, isFF, isRF, isIA, isDip;
88 
89  // Hadronic invariant mass.
90  double shh;
91 
92  // Initialization.
93  bool isInitPtr, isInit;
94  int verbose;
95 
96 };
97 
98 //==========================================================================
99 
100 // Base class for QED systems.
101 
102 class QEDsystem {
103 
104  public:
105 
106  // Constructor.
107  QEDsystem() : infoPtr(nullptr), partonSystemsPtr(nullptr),
108  particleDataPtr(nullptr), rndmPtr(nullptr), settingsPtr(nullptr),
109  loggerPtr(nullptr), vinComPtr(nullptr), isInitPtr(false), iSys(-1),
110  verbose(0), jNew(0), shat(0.) {;}
111 
112  // Destructor.
113  virtual ~QEDsystem() = default;
114 
115  // Initialize pointers.
116  void initPtr(Info* infoPtrIn, ParticleData* particleDataPtrIn,
117  PartonSystems* partonSystemsPtrIn, Rndm* rndmPtrIn,
118  Settings* settingsPtrIn, VinciaCommon* vinComPtrIn);
119 
120  // Initialise settings for current run.
121  virtual void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
122  int verboseIn) = 0;
123  virtual void setVerbose(int verboseIn) { verbose = verboseIn; }
124  // Prepare a parton system for evolution.
125  virtual void prepare(const int iSysIn, Event &event, const double q2CutIn,
126  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
127  AlphaEM alIn) = 0;
128  // Build parton system.
129  virtual void buildSystem(Event &event) = 0;
130  // Generate a trial scale.
131  virtual double q2Next(Event &event, double q2Start) = 0 ;
132  // Generate kinematics and check for veto of branching.
133  virtual bool acceptTrial(Event &event) = 0;
134  // Update the envent.
135  virtual void updateEvent(Event &event) = 0;
136  // Update the parton systems.
137  virtual void updatePartonSystems();
138  // Print information about the system.
139  virtual void print() = 0;
140 
141  // Methods to tell which type of brancher this is.
142  virtual bool isSplitting() {return false;};
143  virtual bool isInitial() {return false;};
144 
145  protected:
146 
147  // Pointers.
148  Info* infoPtr{};
149  PartonSystems* partonSystemsPtr{};
150  ParticleData* particleDataPtr{};
151  Rndm* rndmPtr{};
152  Settings* settingsPtr{};
153  Logger* loggerPtr{};
154  VinciaCommon* vinComPtr{};
155  bool isInitPtr;
156 
157  // Event system.
158  int iSys;
159  vector<Vec4> pNew;
160 
161  // Verbose setting.
162  int verbose;
163 
164  // Information for partonSystems.
165  int jNew;
166  map<int,int> iReplace;
167  double shat;
168 
169 };
170 
171 //==========================================================================
172 
173 // Class for a QED emission system.
174 
175 class QEDemitSystem : public QEDsystem {
176 
177 public:
178 
179  QEDemitSystem() : shh(-1.), cMat(0.), trialIsVec(false), beamAPtr(nullptr),
180  beamBPtr(nullptr), qedMode(-1), qedModeMPI(-1),
181  scaleRegion(0), emitBelowHad(false), q2Cut(-1.), isInit(false),
182  TINYPDF(-1.), kMapTypeFinal(0) {;}
183 
184  // Initialise settings for current run.
185  void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn, int verboseIn)
186  override;
187  // Prepare a parton system for photon emission evolution
188  void prepare(const int iSysIn, Event &event, const double q2CutIn,
189  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
190  AlphaEM alIn) override;
191  // Set up antenna pairing for incoherent mode.
192  void buildSystem(Event &event) override;
193  // Generate a trial scale.
194  double q2Next(Event &event, double q2Start) override;
195  // Generate kinematics and check veto.
196  bool acceptTrial(Event &event) override;
197  // Update the envent
198  void updateEvent(Event &event) override;
199  // Branching type.
200  bool isInitial() override {return eleTrial->isII || eleTrial->isIF;}
201  // Print the QED emit internal system.
202  void print() override;
203 
204  // Trial antenna function.
205  double aTrial(QEDemitElemental* ele, double sxj, double syj, double sxy);
206  // Physical antenna function.
207  double aPhys (QEDemitElemental* ele, double sxj, double syj, double sxy);
208  // Ratio between PDFs.
209  double pdfRatio(bool isA, double eOld, double eNew, int id, double Qt2);
210 
211  private:
212 
213  // Event system.
214  double shh;
215 
216  // Internal storage.
217  vector<vector<QEDemitElemental> > eleMat;
218  vector<int> iCoh;
219  double cMat;
220  vector<QEDemitElemental> eleVec;
221 
222  // AlphaEM.
223  AlphaEM al;
224 
225  // Evolution window.
226  vector<double> evolutionWindows;
227 
228  // Trial pointer.
229  QEDemitElemental* eleTrial{};
230  bool trialIsVec;
231 
232  // Pointers.
233  BeamParticle* beamAPtr{};
234  BeamParticle* beamBPtr{};
235 
236  // Settings.
237  int qedMode, qedModeMPI;
238  vector<bool> useSpinsQEDNow, useSpinsQED, useSpinsQEDHadDec;
239  int scaleRegion, emitBelowHad, isHadronDecay;
240  double q2Cut;
241 
242  // Initialization.
243  bool isInit;
244 
245  // PDF check.
246  double TINYPDF;
247 
248  // Recoil strategy.
249  int kMapTypeFinal;
250 
251  // Global recoil momenta.
252  Vec4 pRecSum;
253  vector<Vec4> pRec;
254  vector<int> iRec;
255 
256 };
257 
258 //==========================================================================
259 
260 // Class for trial QED splittings.
261 
263 
264 public:
265 
266  // Friends for internal private members.
267  friend class QEDsplitSystem;
268 
269  // Constructor.
270  QEDsplitElemental(Event &event, int iPhotIn, int iSpecIn):
271  iPhot(iPhotIn), iSpec(iSpecIn), ariWeight(0) {
272  m2Ant = max(VinciaConstants::PICO,
273  m2(event[iPhotIn], event[iSpecIn]));
274  sAnt = max(VinciaConstants::PICO,
275  2.*event[iPhotIn].p()*event[iSpecIn].p());
276  m2Spec = max(0., event[iSpecIn].m2());
277  }
278 
279  // Kallen function.
280  double getKallen() {return m2Ant/(m2Ant - m2Spec);}
281 
282 private:
283 
284  // Internal members.
285  int iPhot, iSpec;
286  double m2Spec, m2Ant, sAnt;
287  double ariWeight;
288 };
289 
290 //==========================================================================
291 
292 // Class for a QED splitting system.
293 
294 class QEDsplitSystem : public QEDsystem {
295 
296 public:
297 
298  // Constructor.
300  totIdWeight(-1.), hasTrial(false),
301  q2Trial(-1.), zTrial(-1.), phiTrial(-1.), idTrial(0), nQuark(-1),
302  nLepton(-1), q2Max(-1.), q2Cut(-1.), scaleRegion(0),
303  beamAPtr(nullptr), beamBPtr(nullptr), isInit(false), kMapTypeFinal(0) {;}
304 
305  // Initialize.
306  void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn, int verboseIn)
307  override;
308  // Prepare list of final-state photons - with recoilers - for splittings.
309  void prepare(const int iSysIn, Event &event, const double q2CutIn,
310  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
311  AlphaEM alIn) override;
312  // Build the splitting system.
313  void buildSystem(Event &event) override;
314  // Generate a scale for the system.
315  double q2Next(Event &event, double q2Start) override;
316  // Generate kinematics and check veto.
317  bool acceptTrial(Event &event) override;
318  // Update Event.
319  void updateEvent(Event &event) override;
320  // Branching type: isSplitting() = true.
321  bool isSplitting() override { return true;}
322  // Print the system.
323  void print() override;
324 
325 private:
326 
327  // AlphaEM.
328  AlphaEM al;
329 
330  // Evolution window.
331  vector<double> evolutionWindows;
332 
333  // Weights for splitting IDs.
334  vector<int> ids;
335  vector<double> idWeights;
336  double totIdWeight;
337 
338  // Antennae.
339  vector<QEDsplitElemental> eleVec;
340 
341  // Trial variables.
342  bool hasTrial;
343  double q2Trial, zTrial, phiTrial, idTrial;
344  QEDsplitElemental* eleTrial{};
345 
346  // Settings.
347  int nQuark, nLepton;
348  double q2Max, q2Cut;
349  int scaleRegion;
350 
351  // Pointers.
352  BeamParticle* beamAPtr;
353  BeamParticle* beamBPtr;
354 
355  // Initialization.
356  bool isInit;
357 
358  // Recoil strategy.
359  int kMapTypeFinal;
360 
361 };
362 
363 //==========================================================================
364 
365 // Class for a QED conversion system.
366 
367 class QEDconvSystem : public QEDsystem {
368 
369 public:
370 
371  // Constructor.
372  QEDconvSystem() : totIdWeight(-1.), maxIdWeight(-1.), shh(-1.), s(-1.),
373  iA(-1), iB(-1), isAPhot(false), isBPhot(false), hasTrial(false),
374  iPhotTrial(-1), iSpecTrial(-1), q2Trial(-1.), zTrial(-1.), phiTrial(-1.),
375  idTrial(-1), nQuark(-1), q2Cut(-1.),
376  scaleRegion(0), beamAPtr(nullptr), beamBPtr(nullptr),isInit(false),
377  TINYPDF(-1.) {;}
378 
379  // Initialize.
380  void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn, int verboseIn)
381  override;
382  // Prepare for backwards-evolution of photons.
383  void prepare(const int iSysIn, Event &event, const double q2CutIn,
384  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
385  AlphaEM alIn) override;
386  // Build the system.
387  void buildSystem(Event &event) override;
388  // Generate a trial scale.
389  double q2Next(Event &event, double q2Start) override;
390  // Generate kinematics and check veto.
391  bool acceptTrial(Event &event) override;
392  // Update.
393  void updateEvent(Event &event) override;
394  // Branching type: isInitial() = true.
395  bool isInitial() override {return true;};
396  // Print.
397  void print() override;
398 
399 private:
400 
401  // Trial pdf ratios for conversion.
402  map<int,double> Rhat;
403 
404  // AlphaEM.
405  AlphaEM al;
406 
407  // Evolution window.
408  vector<double> evolutionWindows;
409 
410  // Weights for conversion IDs.
411  vector<int> ids;
412  vector<double> idWeights;
413  double totIdWeight, maxIdWeight;
414  double shh;
415 
416  // Antenna parameters.
417  double s;
418  int iA, iB;
419  bool isAPhot, isBPhot;
420 
421  // Trial variables.
422  bool hasTrial;
423  int iPhotTrial, iSpecTrial;
424  double q2Trial, zTrial, phiTrial, idTrial;
425 
426  // Settings.
427  int nQuark;
428  double q2Cut;
429 
430  // Scale Region determines what region and scale we are showering at.
431  // Scale Region = 0 is above the hadronisation scale.
432  // Scale Region = 1 is below the hadronisation scale.
433  // Scale Region = 2 is below the hadronisation scale and contains remanants.
434  int scaleRegion;
435 
436  // Pointers.
437  BeamParticle* beamAPtr;
438  BeamParticle* beamBPtr;
439 
440  // Initialization.
441  bool isInit;
442 
443  // PDF check.
444  double TINYPDF;
445 
446  // Global recoil momenta
447  vector<Vec4> pRec;
448  vector<int> iRec;
449 
450 };
451 
452 //==========================================================================
453 
454 // Base class for Vincia's QED and EW shower modules.
455 
457 
458 public:
459 
460  // Default constructor.
461  VinciaModule(): verbose(-1), isInitPtr(false), isInitSav(false) {;};
462 
463  // Default destructor.
464  virtual ~VinciaModule() = default;
465 
466  // Initialise pointers (called at construction time).
467  virtual void initPtr(Info* infoPtrIn, VinciaCommon* vinComPtrIn) = 0;
468 
469  //Some early initialisation (needed for EW shower).
470  virtual void load() {;}
471 
472  // Initialise settings for current run (called as part of Pythia::init()).
473  virtual void init(BeamParticle* beamAPtrIn = 0, BeamParticle* beamBPtrIn = 0)
474  = 0;
475  bool isInit() {return isInitSav;}
476 
477  // Select helicities for a system of particles.
478  virtual bool polarise(vector<Particle>&) {return false;}
479 
480  // Prepare to shower a system.
481  virtual bool prepare(int iSysIn, Event &event, int scaleRegionIn) = 0;
482 
483  // Update shower system each time something has changed in event.
484  virtual void update(Event &event, int iSys) = 0;
485 
486  // Set verbosity level.
487  virtual void setVerbose(int verboseIn) {verbose = verboseIn;}
488 
489  // Clear everything (optionally for specific system).
490  virtual void clear(int iSys = -1) = 0;
491 
492  // Generate a trial scale.
493  virtual double q2Next(Event &event, double q2Start, double q2End) = 0;
494 
495  // Which system does the winner belong to?
496  virtual int sysWin() = 0;
497 
498  // Get information about the latest branching.
499  virtual bool lastIsSplitting() = 0;
500  virtual bool lastIsInitial() = 0;
501  virtual bool lastIsResonanceDecay() {return false;}
502 
503  // Generate kinematics and check veto.
504  virtual bool acceptTrial(Event &event) = 0;
505 
506  // Update event after branching accepted.
507  virtual void updateEvent(Event &event) = 0;
508 
509  // Update partonSystems after branching accepted.
510  virtual void updatePartonSystems(Event &event) = 0;
511 
512  // End scales.
513  virtual double q2minColoured() = 0;
514  virtual double q2min() = 0;
515 
516  // Get number of branchers / systems.
517  virtual unsigned int nBranchers() = 0;
518  virtual unsigned int nResDec() = 0;
519 
520  // Members.
521  BeamParticle* beamAPtr{};
522  BeamParticle* beamBPtr{};
523  Info* infoPtr{};
524  ParticleData* particleDataPtr{};
525  Logger* loggerPtr{};
526  PartonSystems* partonSystemsPtr{};
527  Rndm* rndmPtr{};
528  Settings* settingsPtr{};
529  VinciaCommon* vinComPtr{};
530 
531  protected:
532 
533  int verbose;
534  bool isInitPtr, isInitSav;
535 
536 };
537 
538 //==========================================================================
539 
540 // Class for performing QED showers.
541 
542 class VinciaQED : public VinciaModule {
543 
544 public:
545 
546  // Friends for internal private members.
547  friend class VinciaFSR;
548 
549  // Constructor.
550  VinciaQED() {;}
551 
552  // Initialise pointers (called at construction time).
553  void initPtr(Info* infoPtrIn, VinciaCommon* vinComPtrIn) override;
554  // Initialise settings for current run (called as part of Pythia::init()).
555  void init(BeamParticle* beamAPtrIn = 0, BeamParticle* beamBPtrIn = 0)
556  override;
557  // Prepare to shower a system.
558  bool prepare(int iSysIn, Event& event, int scaleRegionIn) override;
559  // Update QED shower system(s) each time something has changed in event.
560  void update(Event& event, int iSys) override;
561  // Set or change verbosity level, and propagate to QED systems.
562  virtual void setVerbose(int verboseIn) override {
563  verbose = verboseIn;
564  emptyQEDemitSystem.setVerbose(verboseIn);
565  emptyQEDsplitSystem.setVerbose(verboseIn);
566  emptyQEDconvSystem.setVerbose(verboseIn);
567  }
568  // Clear everything, or specific system.
569  void clear(int iSys = -1) override {
570  if (iSys < 0) {emitSystems.clear(); splitSystems.clear();
571  convSystems.clear();}
572  else {emitSystems.erase(iSys); splitSystems.erase(iSys);
573  convSystems.erase(iSys);}
574  qedTrialSysPtr = nullptr;}
575 
576  // Generate a trial scale.
577  double q2Next(Event& event, double q2Start, double) override;
578  // Return the system window.
579  int sysWin() override {return iSysTrial;}
580  // Information about last branching.
581  bool lastIsSplitting() override {
582  if (qedTrialSysPtr != nullptr) return qedTrialSysPtr->isSplitting();
583  else return false;}
584  bool lastIsInitial() override {
585  if (qedTrialSysPtr != nullptr) return qedTrialSysPtr->isInitial();
586  else return false;}
587  // Generate kinematics and check veto.
588  bool acceptTrial(Event &event) override;
589  // Update event after branching accepted.
590  void updateEvent(Event &event) override;
591  // Update partonSystems after branching accepted.
592  void updatePartonSystems(Event &event) override;
593  // Return scales.
594  double q2minColoured() override {return q2minColouredSav;}
595  double q2min() override {return q2minSav;}
596 
597  // Getter for number of systems.
598  unsigned int nBranchers() override {
599  int sizeNow = emitSystems.size();
600  sizeNow = max(sizeNow, (int)splitSystems.size());
601  sizeNow = max(sizeNow, (int)convSystems.size());
602  return sizeNow;}
603 
604  // This module does not implement resonance decays.
605  unsigned int nResDec() override { return 0; }
606 
607 private:
608 
609  // Get Q2 for QED system.
610  template <class T>
611  void q2NextSystem(map<int, T>& QEDsystemList, Event& event, double q2Start);
612 
613  // Systems.
614  QEDemitSystem emptyQEDemitSystem;
615  QEDsplitSystem emptyQEDsplitSystem;
616  QEDconvSystem emptyQEDconvSystem;
617  map< int, QEDemitSystem> emitSystems;
618  map< int, QEDsplitSystem> splitSystems;
619  map< int, QEDconvSystem> convSystems;
620 
621  // Settings.
622  bool doQED, doEmission;
623  int nGammaToLepton, nGammaToQuark;
624  bool doConvertGamma, doConvertQuark;
625 
626  // Scales.
627  double q2minSav, q2minColouredSav;
628 
629  // Trial information
630  int iSysTrial;
631  double q2Trial;
632  QEDsystem* qedTrialSysPtr{};
633 
634  // AlphaEM.
635  AlphaEM al;
636 
637  // Evolution windows
638  vector<double> evolutionWindows;
639 
640 };
641 
642 //==========================================================================
643 
644 } // end namespace Pythia8
645 
646 #endif // Pythia8_VinciaQED_H
bool isInitial() override
Branching type.
Definition: VinciaQED.h:200
virtual bool isSplitting()
Methods to tell which type of brancher this is.
Definition: VinciaQED.h:142
int jNew
Information for partonSystems.
Definition: VinciaQED.h:165
Class for a QED emission system.
Definition: VinciaQED.h:175
bool isSplitting() override
Branching type: isSplitting() = true.
Definition: VinciaQED.h:321
Definition: Info.h:45
The Event class holds all info on the generated event.
Definition: Event.h:408
int iSys
Event system.
Definition: VinciaQED.h:158
bool isInitial() override
Branching type: isInitial() = true.
Definition: VinciaQED.h:395
Definition: BeamParticle.h:133
QEDsplitSystem()
Constructor.
Definition: VinciaQED.h:299
unsigned int nResDec() override
This module does not implement resonance decays.
Definition: VinciaQED.h:605
void init(Event &event, int xIn, int yIn, double shhIn, double verboseIn)
Initialize.
Definition: VinciaQED.cc:35
unsigned int nBranchers() override
Getter for number of systems.
Definition: VinciaQED.h:598
Base class for Vincia&#39;s QED and EW shower modules.
Definition: VinciaQED.h:456
virtual bool polarise(vector< Particle > &)
Select helicities for a system of particles.
Definition: VinciaQED.h:478
friend class QEDemitSystem
Friends for internal private members.
Definition: VinciaQED.h:33
double generateTrial(Event &event, double q2Start, double q2Low, double alphaIn, double cIn)
Generate a trial point.
Definition: VinciaQED.cc:133
virtual void setVerbose(int verboseIn) override
Set or change verbosity level, and propagate to QED systems.
Definition: VinciaQED.h:562
Class for a QED conversion system.
Definition: VinciaQED.h:367
Definition: StandardModel.h:106
Class for performing QED showers.
Definition: VinciaQED.h:542
Definition: Logger.h:23
QEDconvSystem()
Constructor.
Definition: VinciaQED.h:372
virtual void load()
Some early initialisation (needed for EW shower).
Definition: VinciaQED.h:470
VinciaModule()
Default constructor.
Definition: VinciaQED.h:461
VinciaQED()
Constructor.
Definition: VinciaQED.h:550
Definition: VinciaFSR.h:582
Definition: Basics.h:388
Definition: VinciaCommon.h:494
QEDsplitElemental(Event &event, int iPhotIn, int iSpecIn)
Constructor.
Definition: VinciaQED.h:270
Base class for QED systems.
Definition: VinciaQED.h:102
double q2minColoured() override
Return scales.
Definition: VinciaQED.h:594
QEDemitElemental()
Constuctor.
Definition: VinciaQED.h:36
bool lastIsSplitting() override
Information about last branching.
Definition: VinciaQED.h:581
Class for trial QED splittings.
Definition: VinciaQED.h:262
Class for QED emissions.
Definition: VinciaQED.h:28
double m2(const Vec4 &v1)
The squared invariant mass of one or more four-vectors.
Definition: Basics.cc:598
void clear(int iSys=-1) override
Clear everything, or specific system.
Definition: VinciaQED.h:569
The PartonSystems class describes the whole set of subcollisions.
Definition: PartonSystems.h:42
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
void initPtr(Rndm *rndmPtrIn, PartonSystems *partonSystemsPtrIn)
Initialize the pointers.
Definition: VinciaQED.cc:24
double getKallen()
Kallen function.
Definition: VinciaQED.h:280
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
int verbose
Verbose setting.
Definition: VinciaQED.h:162
virtual void setVerbose(int verboseIn)
Set verbosity level.
Definition: VinciaQED.h:487
Definition: Basics.h:32
int sysWin() override
Return the system window.
Definition: VinciaQED.h:579
Definition: Settings.h:196
QEDsystem()
Constructor.
Definition: VinciaQED.h:107
Class for a QED splitting system.
Definition: VinciaQED.h:294