PYTHIA  8.314
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(BeamParticlePtr beamAPtrIn, BeamParticlePtr 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(BeamParticlePtr beamAPtrIn, BeamParticlePtr beamBPtrIn,
186  int verboseIn) 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  // Beam pointers.
233  BeamParticlePtr beamAPtr{}, beamBPtr{};
234 
235  // Settings.
236  int qedMode, qedModeMPI;
237  vector<bool> useSpinsQEDNow, useSpinsQED, useSpinsQEDHadDec;
238  int scaleRegion, emitBelowHad, isHadronDecay;
239  double q2Cut;
240 
241  // Initialization.
242  bool isInit;
243 
244  // PDF check.
245  double TINYPDF;
246 
247  // Recoil strategy.
248  int kMapTypeFinal;
249 
250  // Global recoil momenta.
251  Vec4 pRecSum;
252  vector<Vec4> pRec;
253  vector<int> iRec;
254 
255 };
256 
257 //==========================================================================
258 
259 // Class for trial QED splittings.
260 
262 
263 public:
264 
265  // Friends for internal private members.
266  friend class QEDsplitSystem;
267 
268  // Constructor.
269  QEDsplitElemental(Event &event, int iPhotIn, int iSpecIn):
270  iPhot(iPhotIn), iSpec(iSpecIn), ariWeight(0) {
271  m2Ant = max(VinciaConstants::PICO,
272  m2(event[iPhotIn], event[iSpecIn]));
273  sAnt = max(VinciaConstants::PICO,
274  2.*event[iPhotIn].p()*event[iSpecIn].p());
275  m2Spec = max(0., event[iSpecIn].m2());
276  }
277 
278  // Kallen function.
279  double getKallen() {return m2Ant/(m2Ant - m2Spec);}
280 
281 private:
282 
283  // Internal members.
284  int iPhot, iSpec;
285  double m2Spec, m2Ant, sAnt;
286  double ariWeight;
287 };
288 
289 //==========================================================================
290 
291 // Class for a QED splitting system.
292 
293 class QEDsplitSystem : public QEDsystem {
294 
295 public:
296 
297  // Constructor.
299  totIdWeight(-1.), hasTrial(false),
300  q2Trial(-1.), zTrial(-1.), phiTrial(-1.), idTrial(0), nQuark(-1),
301  nLepton(-1), q2Max(-1.), q2Cut(-1.), scaleRegion(0),
302  beamAPtr(nullptr), beamBPtr(nullptr), isInit(false), kMapTypeFinal(0) {;}
303 
304  // Initialize.
305  void init(BeamParticlePtr beamAPtrIn, BeamParticlePtr beamBPtrIn,
306  int verboseIn) override;
307  // Prepare list of final-state photons - with recoilers - for splittings.
308  void prepare(const int iSysIn, Event &event, const double q2CutIn,
309  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
310  AlphaEM alIn) override;
311  // Build the splitting system.
312  void buildSystem(Event &event) override;
313  // Generate a scale for the system.
314  double q2Next(Event &event, double q2Start) override;
315  // Generate kinematics and check veto.
316  bool acceptTrial(Event &event) override;
317  // Update Event.
318  void updateEvent(Event &event) override;
319  // Branching type: isSplitting() = true.
320  bool isSplitting() override { return true;}
321  // Print the system.
322  void print() override;
323 
324 private:
325 
326  // AlphaEM.
327  AlphaEM al;
328 
329  // Evolution window.
330  vector<double> evolutionWindows;
331 
332  // Weights for splitting IDs.
333  vector<int> ids;
334  vector<double> idWeights;
335  double totIdWeight;
336 
337  // Antennae.
338  vector<QEDsplitElemental> eleVec;
339 
340  // Trial variables.
341  bool hasTrial;
342  double q2Trial, zTrial, phiTrial, idTrial;
343  QEDsplitElemental* eleTrial{};
344 
345  // Settings.
346  int nQuark, nLepton;
347  double q2Max, q2Cut;
348  int scaleRegion;
349 
350  // Pointers.
351  BeamParticlePtr beamAPtr;
352  BeamParticlePtr beamBPtr;
353 
354  // Initialization.
355  bool isInit;
356 
357  // Recoil strategy.
358  int kMapTypeFinal;
359 
360 };
361 
362 //==========================================================================
363 
364 // Class for a QED conversion system.
365 
366 class QEDconvSystem : public QEDsystem {
367 
368 public:
369 
370  // Constructor.
371  QEDconvSystem() : totIdWeight(-1.), maxIdWeight(-1.), shh(-1.), s(-1.),
372  iA(-1), iB(-1), isAPhot(false), isBPhot(false), hasTrial(false),
373  iPhotTrial(-1), iSpecTrial(-1), q2Trial(-1.), zTrial(-1.), phiTrial(-1.),
374  idTrial(-1), nQuark(-1), q2Cut(-1.),
375  scaleRegion(0), beamAPtr(nullptr), beamBPtr(nullptr),isInit(false),
376  TINYPDF(-1.) {;}
377 
378  // Initialize.
379  void init(BeamParticlePtr beamAPtrIn, BeamParticlePtr beamBPtrIn,
380  int verboseIn) override;
381  // Prepare for backwards-evolution of photons.
382  void prepare(const int iSysIn, Event &event, const double q2CutIn,
383  const int scaleRegionIn, const vector<double> evolutionWindowsIn,
384  AlphaEM alIn) override;
385  // Build the system.
386  void buildSystem(Event &event) override;
387  // Generate a trial scale.
388  double q2Next(Event &event, double q2Start) override;
389  // Generate kinematics and check veto.
390  bool acceptTrial(Event &event) override;
391  // Update.
392  void updateEvent(Event &event) override;
393  // Branching type: isInitial() = true.
394  bool isInitial() override {return true;};
395  // Print.
396  void print() override;
397 
398 private:
399 
400  // Trial pdf ratios for conversion.
401  map<int,double> Rhat;
402 
403  // AlphaEM.
404  AlphaEM al;
405 
406  // Evolution window.
407  vector<double> evolutionWindows;
408 
409  // Weights for conversion IDs.
410  vector<int> ids;
411  vector<double> idWeights;
412  double totIdWeight, maxIdWeight;
413  double shh;
414 
415  // Antenna parameters.
416  double s;
417  int iA, iB;
418  bool isAPhot, isBPhot;
419 
420  // Trial variables.
421  bool hasTrial;
422  int iPhotTrial, iSpecTrial;
423  double q2Trial, zTrial, phiTrial, idTrial;
424 
425  // Settings.
426  int nQuark;
427  double q2Cut;
428 
429  // Scale Region determines what region and scale we are showering at.
430  // Scale Region = 0 is above the hadronisation scale.
431  // Scale Region = 1 is below the hadronisation scale.
432  // Scale Region = 2 is below the hadronisation scale and contains remanants.
433  int scaleRegion;
434 
435  // Pointers.
436  BeamParticlePtr beamAPtr;
437  BeamParticlePtr beamBPtr;
438 
439  // Initialization.
440  bool isInit;
441 
442  // PDF check.
443  double TINYPDF;
444 
445  // Global recoil momenta
446  vector<Vec4> pRec;
447  vector<int> iRec;
448 
449 };
450 
451 //==========================================================================
452 
453 // Base class for Vincia's QED and EW shower modules.
454 
456 
457 public:
458 
459  // Default constructor.
460  VinciaModule(): verbose(-1), isInitPtr(false), isInitSav(false) {;};
461 
462  // Default destructor.
463  virtual ~VinciaModule() = default;
464 
465  // Initialise pointers (called at construction time).
466  virtual void initPtr(Info* infoPtrIn, VinciaCommon* vinComPtrIn) = 0;
467 
468  //Some early initialisation (needed for EW shower).
469  virtual void load() {;}
470 
471  // Initialise settings for current run (called as part of Pythia::init()).
472  virtual void init(BeamParticlePtr beamAPtrIn = nullptr,
473  BeamParticlePtr beamBPtrIn = nullptr) = 0;
474  bool isInit() {return isInitSav;}
475 
476  // Select helicities for a system of particles.
477  virtual bool polarise(vector<Particle>&) {return false;}
478 
479  // Prepare to shower a system.
480  virtual bool prepare(int iSysIn, Event &event, int scaleRegionIn) = 0;
481 
482  // Update shower system each time something has changed in event.
483  virtual void update(Event &event, int iSys) = 0;
484 
485  // Set verbosity level.
486  virtual void setVerbose(int verboseIn) {verbose = verboseIn;}
487 
488  // Clear everything (optionally for specific system).
489  virtual void clear(int iSys = -1) = 0;
490 
491  // Generate a trial scale.
492  virtual double q2Next(Event &event, double q2Start, double q2End) = 0;
493 
494  // Which system does the winner belong to?
495  virtual int sysWin() = 0;
496 
497  // Get information about the latest branching.
498  virtual bool lastIsSplitting() = 0;
499  virtual bool lastIsInitial() = 0;
500  virtual bool lastIsResonanceDecay() {return false;}
501 
502  // Generate kinematics and check veto.
503  virtual bool acceptTrial(Event &event) = 0;
504 
505  // Update event after branching accepted.
506  virtual void updateEvent(Event &event) = 0;
507 
508  // Update partonSystems after branching accepted.
509  virtual void updatePartonSystems(Event &event) = 0;
510 
511  // End scales.
512  virtual double q2minColoured() = 0;
513  virtual double q2min() = 0;
514 
515  // Get number of branchers / systems.
516  virtual unsigned int nBranchers() = 0;
517  virtual unsigned int nResDec() = 0;
518 
519  // Members.
520  BeamParticlePtr beamAPtr{};
521  BeamParticlePtr beamBPtr{};
522  Info* infoPtr{};
523  ParticleData* particleDataPtr{};
524  Logger* loggerPtr{};
525  PartonSystems* partonSystemsPtr{};
526  Rndm* rndmPtr{};
527  Settings* settingsPtr{};
528  VinciaCommon* vinComPtr{};
529 
530  protected:
531 
532  int verbose;
533  bool isInitPtr, isInitSav;
534 
535 };
536 
537 //==========================================================================
538 
539 // Class for performing QED showers.
540 
541 class VinciaQED : public VinciaModule {
542 
543 public:
544 
545  // Friends for internal private members.
546  friend class VinciaFSR;
547 
548  // Constructor.
549  VinciaQED() {;}
550 
551  // Initialise pointers (called at construction time).
552  void initPtr(Info* infoPtrIn, VinciaCommon* vinComPtrIn) override;
553  // Initialise settings for current run (called as part of Pythia::init()).
554  void init(BeamParticlePtr beamAPtrIn = 0, BeamParticlePtr beamBPtrIn = 0)
555  override;
556  // Prepare to shower a system.
557  bool prepare(int iSysIn, Event& event, int scaleRegionIn) override;
558  // Update QED shower system(s) each time something has changed in event.
559  void update(Event& event, int iSys) override;
560  // Set or change verbosity level, and propagate to QED systems.
561  virtual void setVerbose(int verboseIn) override {
562  verbose = verboseIn;
563  emptyQEDemitSystem.setVerbose(verboseIn);
564  emptyQEDsplitSystem.setVerbose(verboseIn);
565  emptyQEDconvSystem.setVerbose(verboseIn);
566  }
567  // Clear everything, or specific system.
568  void clear(int iSys = -1) override {
569  if (iSys < 0) {emitSystems.clear(); splitSystems.clear();
570  convSystems.clear();}
571  else {emitSystems.erase(iSys); splitSystems.erase(iSys);
572  convSystems.erase(iSys);}
573  qedTrialSysPtr = nullptr;}
574 
575  // Generate a trial scale.
576  double q2Next(Event& event, double q2Start, double) override;
577  // Return the system window.
578  int sysWin() override {return iSysTrial;}
579  // Information about last branching.
580  bool lastIsSplitting() override {
581  if (qedTrialSysPtr != nullptr) return qedTrialSysPtr->isSplitting();
582  else return false;}
583  bool lastIsInitial() override {
584  if (qedTrialSysPtr != nullptr) return qedTrialSysPtr->isInitial();
585  else return false;}
586  // Generate kinematics and check veto.
587  bool acceptTrial(Event &event) override;
588  // Update event after branching accepted.
589  void updateEvent(Event &event) override;
590  // Update partonSystems after branching accepted.
591  void updatePartonSystems(Event &event) override;
592  // Return scales.
593  double q2minColoured() override {return q2minColouredSav;}
594  double q2min() override {return q2minSav;}
595 
596  // Getter for number of systems.
597  unsigned int nBranchers() override {
598  int sizeNow = emitSystems.size();
599  sizeNow = max(sizeNow, (int)splitSystems.size());
600  sizeNow = max(sizeNow, (int)convSystems.size());
601  return sizeNow;}
602 
603  // This module does not implement resonance decays.
604  unsigned int nResDec() override { return 0; }
605 
606 private:
607 
608  // Get Q2 for QED system.
609  template <class T>
610  void q2NextSystem(map<int, T>& QEDsystemList, Event& event, double q2Start);
611 
612  // Systems.
613  QEDemitSystem emptyQEDemitSystem;
614  QEDsplitSystem emptyQEDsplitSystem;
615  QEDconvSystem emptyQEDconvSystem;
616  map< int, QEDemitSystem> emitSystems;
617  map< int, QEDsplitSystem> splitSystems;
618  map< int, QEDconvSystem> convSystems;
619 
620  // Settings.
621  bool doQED, doEmission;
622  int nGammaToLepton, nGammaToQuark;
623  bool doConvertGamma, doConvertQuark;
624 
625  // Scales.
626  double q2minSav, q2minColouredSav;
627 
628  // Trial information
629  int iSysTrial;
630  double q2Trial;
631  QEDsystem* qedTrialSysPtr{};
632 
633  // AlphaEM.
634  AlphaEM al;
635 
636  // Evolution windows
637  vector<double> evolutionWindows;
638 
639 };
640 
641 //==========================================================================
642 
643 } // end namespace Pythia8
644 
645 #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:320
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:394
QEDsplitSystem()
Constructor.
Definition: VinciaQED.h:298
unsigned int nResDec() override
This module does not implement resonance decays.
Definition: VinciaQED.h:604
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:597
Base class for Vincia&#39;s QED and EW shower modules.
Definition: VinciaQED.h:455
virtual bool polarise(vector< Particle > &)
Select helicities for a system of particles.
Definition: VinciaQED.h:477
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:561
Class for a QED conversion system.
Definition: VinciaQED.h:366
Definition: StandardModel.h:106
Class for performing QED showers.
Definition: VinciaQED.h:541
Definition: Logger.h:23
QEDconvSystem()
Constructor.
Definition: VinciaQED.h:371
virtual void load()
Some early initialisation (needed for EW shower).
Definition: VinciaQED.h:469
VinciaModule()
Default constructor.
Definition: VinciaQED.h:460
VinciaQED()
Constructor.
Definition: VinciaQED.h:549
Definition: VinciaFSR.h:582
Definition: Basics.h:388
Definition: VinciaCommon.h:494
QEDsplitElemental(Event &event, int iPhotIn, int iSpecIn)
Constructor.
Definition: VinciaQED.h:269
Base class for QED systems.
Definition: VinciaQED.h:102
double q2minColoured() override
Return scales.
Definition: VinciaQED.h:593
QEDemitElemental()
Constuctor.
Definition: VinciaQED.h:36
bool lastIsSplitting() override
Information about last branching.
Definition: VinciaQED.h:580
Class for trial QED splittings.
Definition: VinciaQED.h:261
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:568
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:279
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:486
Definition: Basics.h:32
int sysWin() override
Return the system window.
Definition: VinciaQED.h:578
Definition: Settings.h:196
QEDsystem()
Constructor.
Definition: VinciaQED.h:107
Class for a QED splitting system.
Definition: VinciaQED.h:293