9 #ifndef Pythia8_VinciaQED_H 10 #define Pythia8_VinciaQED_H 13 #include "Pythia8/BeamParticle.h" 14 #include "Pythia8/Event.h" 15 #include "Pythia8/StandardModel.h" 16 #include "Pythia8/PartonSystems.h" 19 #include "Pythia8/VinciaCommon.h" 20 #include "Pythia8/VinciaWeights.h" 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) {;}
46 void init(
Event &event,
int xIn,
int yIn,
double shhIn,
49 void init(
Event &event,
int xIn, vector<int> iRecoilIn,
double shhIn,
53 double alphaIn,
double cIn);
64 double q2Sav, zetaSav, phiSav, sxjSav, syjSav, alpha, c;
85 bool isII, isIF, isFF, isRF, isIA, isDip;
91 bool isInitPtr, isInit;
105 QEDsystem() : infoPtr(nullptr), partonSystemsPtr(nullptr),
106 particleDataPtr(nullptr), rndmPtr(nullptr), settingsPtr(nullptr),
107 loggerPtr(nullptr), vinComPtr(nullptr), isInitPtr(false), iSys(-1),
108 verbose(0), jNew(0), shat(0.) {;}
121 virtual void setVerbose(
int verboseIn) { verbose = verboseIn; }
123 virtual void prepare(
int iSysIn,
Event &event,
double q2CutIn,
124 bool isBelowHadIn, vector<double> evolutionWindowsIn,
AlphaEM alIn) = 0;
126 virtual void buildSystem(
Event &event) = 0;
128 virtual double q2Next(
Event &event,
double q2Start) = 0 ;
130 virtual bool acceptTrial(
Event &event) = 0;
132 virtual void updateEvent(
Event &event) = 0;
134 virtual void updatePartonSystems();
136 virtual void print() = 0;
140 virtual bool isInitial() {
return false;};
163 map<int,int> iReplace;
176 QEDemitSystem() : shh(-1.), cMat(0.), trialIsVec(
false), beamAPtr(
nullptr),
177 beamBPtr(
nullptr), qedMode(-1), qedModeMPI(-1), useFullWkernel(
false),
178 isBelowHad(
false), emitBelowHad(
false), q2Cut(-1.), isInit(
false),
179 TINYPDF(-1.), kMapTypeFinal(0) {;}
185 void prepare(
int iSysIn,
Event &event,
double q2CutIn,
bool isBelowHadIn,
186 vector<double> evolutionWindowsIn,
AlphaEM alIn)
override;
188 void buildSystem(
Event &event)
override;
190 double q2Next(
Event &event,
double q2Start)
override;
192 bool acceptTrial(
Event &event)
override;
194 void updateEvent(
Event &event)
override;
196 bool isInitial()
override {
return eleTrial->isII || eleTrial->isIF;}
198 void print()
override;
205 double pdfRatio(
bool isA,
double eOld,
double eNew,
int id,
double Qt2);
213 vector<vector<QEDemitElemental> > eleMat;
216 vector<QEDemitElemental> eleVec;
222 vector<double> evolutionWindows;
233 int qedMode, qedModeMPI;
234 bool useFullWkernel, isBelowHad, emitBelowHad;
266 iPhot(iPhotIn), iSpec(iSpecIn), ariWeight(0) {
267 m2Ant = max(VinciaConstants::PICO,
268 m2(event[iPhotIn], event[iSpecIn]));
269 sAnt = max(VinciaConstants::PICO,
270 2.*event[iPhotIn].p()*event[iSpecIn].p());
271 m2Spec = max(0., event[iSpecIn].
m2());
281 double m2Spec, m2Ant, sAnt;
295 totIdWeight(-1.), hasTrial(false),
296 q2Trial(-1.), zTrial(-1.), phiTrial(-1.), idTrial(0), nQuark(-1),
297 nLepton(-1), q2Max(-1.), q2Cut(-1.), isBelowHad(false),
298 beamAPtr(nullptr), beamBPtr(nullptr), isInit(false), kMapTypeFinal(0) {;}
304 void prepare(
int iSysIn,
Event &event,
double q2CutIn,
bool isBelowHadIn,
305 vector<double> evolutionWindowsIn,
AlphaEM alIn)
override;
307 void buildSystem(
Event &event)
override;
309 double q2Next(
Event &event,
double q2Start)
override;
311 bool acceptTrial(
Event &event)
override;
313 void updateEvent(
Event &event)
override;
317 void print()
override;
325 vector<double> evolutionWindows;
329 vector<double> idWeights;
333 vector<QEDsplitElemental> eleVec;
337 double q2Trial, zTrial, phiTrial, idTrial;
367 iA(-1), iB(-1), isAPhot(false), isBPhot(false), hasTrial(false),
368 iPhotTrial(-1), iSpecTrial(-1), q2Trial(-1.), zTrial(-1.), phiTrial(-1.),
369 idTrial(-1), nQuark(-1), q2Cut(-1.),
370 isBelowHad(false), beamAPtr(nullptr), beamBPtr(nullptr),isInit(false),
377 void prepare(
int iSysIn,
Event &event,
double q2CutIn,
bool isBelowHadIn,
378 vector<double> evolutionWindowsIn,
AlphaEM alIn)
override;
380 void buildSystem(
Event &event)
override;
382 double q2Next(
Event &event,
double q2Start)
override;
384 bool acceptTrial(
Event &event)
override;
386 void updateEvent(
Event &event)
override;
390 void print()
override;
395 map<int,double> Rhat;
401 vector<double> evolutionWindows;
405 vector<double> idWeights;
406 double totIdWeight, maxIdWeight;
412 bool isAPhot, isBPhot;
416 int iPhotTrial, iSpecTrial;
417 double q2Trial, zTrial, phiTrial, idTrial;
463 bool isInit() {
return isInitSav;}
466 virtual bool polarise(vector<Particle>&) {
return false;}
469 virtual bool prepare(
int iSysIn,
Event &event,
bool isBelowHadIn) = 0;
472 virtual void update(
Event &event,
int iSys) = 0;
475 virtual void setVerbose(
int verboseIn) {verbose = verboseIn;}
478 virtual void clear(
int iSys = -1) = 0;
481 virtual double q2Next(
Event &event,
double q2Start,
double q2End) = 0;
484 virtual int sysWin() = 0;
487 virtual bool lastIsSplitting() = 0;
488 virtual bool lastIsInitial() = 0;
489 virtual bool lastIsResonanceDecay() {
return false;}
492 virtual bool acceptTrial(
Event &event) = 0;
495 virtual void updateEvent(
Event &event) = 0;
498 virtual void updatePartonSystems(
Event &event) = 0;
501 virtual double q2minColoured() = 0;
502 virtual double q2min() = 0;
505 virtual unsigned int nBranchers() = 0;
506 virtual unsigned int nResDec() = 0;
522 bool isInitPtr, isInitSav;
546 bool prepare(
int iSysIn,
Event& event,
bool isBelowHadIn)
override;
548 void update(
Event& event,
int iSys)
override;
552 emptyQEDemitSystem.setVerbose(verboseIn);
553 emptyQEDsplitSystem.setVerbose(verboseIn);
554 emptyQEDconvSystem.setVerbose(verboseIn);
557 void clear(
int iSys = -1)
override {
558 if (iSys < 0) {emitSystems.clear(); splitSystems.clear();
559 convSystems.clear();}
560 else {emitSystems.erase(iSys); splitSystems.erase(iSys);
561 convSystems.erase(iSys);}
562 qedTrialSysPtr =
nullptr;}
565 double q2Next(
Event& event,
double q2Start,
double)
override;
567 int sysWin()
override {
return iSysTrial;}
570 if (qedTrialSysPtr !=
nullptr)
return qedTrialSysPtr->isSplitting();
572 bool lastIsInitial()
override {
573 if (qedTrialSysPtr !=
nullptr)
return qedTrialSysPtr->isInitial();
576 bool acceptTrial(
Event &event)
override;
578 void updateEvent(
Event &event)
override;
580 void updatePartonSystems(
Event &event)
override;
583 double q2min()
override {
return q2minSav;}
587 int sizeNow = emitSystems.size();
588 sizeNow = max(sizeNow, (
int)splitSystems.size());
589 sizeNow = max(sizeNow, (
int)convSystems.size());
599 void q2NextSystem(map<int, T>& QEDsystemList,
Event& event,
double q2Start);
605 map< int, QEDemitSystem> emitSystems;
606 map< int, QEDsplitSystem> splitSystems;
607 map< int, QEDconvSystem> convSystems;
610 bool doQED, doEmission;
611 int nGammaToLepton, nGammaToQuark;
612 bool doConvertGamma, doConvertQuark;
615 double q2minSav, q2minColouredSav;
626 vector<double> evolutionWindows;
bool isInitial() override
Branching type.
Definition: VinciaQED.h:196
virtual bool isSplitting()
Methods to tell which type of brancher this is.
Definition: VinciaQED.h:139
int jNew
Information for partonSystems.
Definition: VinciaQED.h:162
Class for a QED emission system.
Definition: VinciaQED.h:172
bool isSplitting() override
Branching type: isSplitting() = true.
Definition: VinciaQED.h:315
The Event class holds all info on the generated event.
Definition: Event.h:453
int iSys
Event system.
Definition: VinciaQED.h:155
bool isInitial() override
Branching type: isInitial() = true.
Definition: VinciaQED.h:388
Definition: BeamParticle.h:133
QEDsplitSystem()
Constructor.
Definition: VinciaQED.h:294
unsigned int nResDec() override
This module does not implement resonance decays.
Definition: VinciaQED.h:593
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:586
Base class for Vincia's QED and EW shower modules.
Definition: VinciaQED.h:444
virtual bool polarise(vector< Particle > &)
Select helicities for a system of particles.
Definition: VinciaQED.h:466
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:132
virtual void setVerbose(int verboseIn) override
Set or change verbosity level, and propagate to QED systems.
Definition: VinciaQED.h:550
Class for a QED conversion system.
Definition: VinciaQED.h:361
Definition: StandardModel.h:106
Class for performing QED showers.
Definition: VinciaQED.h:530
QEDconvSystem()
Constructor.
Definition: VinciaQED.h:366
virtual void load()
Some early initialisation (needed for EW shower).
Definition: VinciaQED.h:458
VinciaModule()
Default constructor.
Definition: VinciaQED.h:449
VinciaQED()
Constructor.
Definition: VinciaQED.h:538
Definition: VinciaFSR.h:582
Definition: VinciaCommon.h:494
QEDsplitElemental(Event &event, int iPhotIn, int iSpecIn)
Constructor.
Definition: VinciaQED.h:265
Base class for QED systems.
Definition: VinciaQED.h:100
double q2minColoured() override
Return scales.
Definition: VinciaQED.h:582
QEDemitElemental()
Constuctor.
Definition: VinciaQED.h:36
bool lastIsSplitting() override
Information about last branching.
Definition: VinciaQED.h:569
Class for trial QED splittings.
Definition: VinciaQED.h:257
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:605
void clear(int iSys=-1) override
Clear everything, or specific system.
Definition: VinciaQED.h:557
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:275
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
int verbose
Verbose setting.
Definition: VinciaQED.h:159
virtual void setVerbose(int verboseIn)
Set verbosity level.
Definition: VinciaQED.h:475
int sysWin() override
Return the system window.
Definition: VinciaQED.h:567
Definition: Settings.h:195
QEDsystem()
Constructor.
Definition: VinciaQED.h:105
Class for a QED splitting system.
Definition: VinciaQED.h:289