3 #include "Pythia8/UserHooks.h" 4 #include "Pythia8/StringFragmentation.h" 13 std::vector<Pythia8::UserHooks*>
hooks()
const {
return hooks_; }
20 for (Pythia8::UserHooks *hook :
hooks_) {
21 hook->initPtr(infoPtr, settingsPtr,
22 particleDataPtr, rndmPtr,
24 beamPomAPtr, beamPomBPtr,
25 coupSMPtr, partonSystemsPtr,
27 test &= hook->initAfterBeams();
34 int nCanVetoMPIStep = 0;
35 int nCanSetResonanceScale = 0;
37 for (Pythia8::UserHooks *hook : hooks_) {
38 if (hook->canVetoPT()) ++nCanVetoPT;
39 if (hook->canVetoStep()) ++nCanVetoStep;
40 if (hook->canVetoMPIStep()) ++nCanVetoMPIStep;
41 if (hook->canSetResonanceScale()) ++nCanSetResonanceScale;
42 if (hook->canEnhanceEmission() || hook->canEnhanceTrial()) ++nCanEnhance;
46 infoPtr->errorMsg(
"Error in MultiUserHook::initAfterBeams " 47 "multiple UserHooks with canVetoPT() not allowed");
51 infoPtr->errorMsg(
"Error in MultiUserHook::initAfterBeams " 52 "multiple UserHooks with canVetoStep() not allowed");
55 if (nCanVetoMPIStep>1) {
56 infoPtr->errorMsg(
"Error in MultiUserHook::initAfterBeams " 57 "multiple UserHooks with canVetoMPIStep() not allowed");
60 if (nCanSetResonanceScale>1) {
61 infoPtr->errorMsg(
"Error in MultiUserHook::initAfterBeams " 62 "multiple UserHooks with canSetResonanceScale() not allowed");
66 infoPtr->errorMsg(
"Error in MultiUserHook::initAfterBeams " 67 "multiple UserHooks with canEnhanceEmission() or canEnhanceTrial() not allowed");
77 for (Pythia8::UserHooks *hook :
hooks_) {
78 test |= hook->canModifySigma();
85 const Pythia8::PhaseSpace* phaseSpacePtr,
bool inEvent)
override {
87 for (Pythia8::UserHooks *hook :
hooks_) {
88 if (hook->canModifySigma()) factor *= hook->multiplySigmaBy(sigmaProcessPtr,phaseSpacePtr,inEvent);
96 for (Pythia8::UserHooks *hook :
hooks_) {
97 test |= hook->canBiasSelection();
104 const Pythia8::PhaseSpace* phaseSpacePtr,
bool inEvent)
override {
106 for (Pythia8::UserHooks *hook :
hooks_) {
107 if (hook->canBiasSelection()) factor *= hook->biasSelectionBy(sigmaProcessPtr,phaseSpacePtr,inEvent);
115 for (Pythia8::UserHooks *hook :
hooks_) {
116 if (hook->canBiasSelection()) factor *= hook->biasedSelectionWeight();
124 for (Pythia8::UserHooks *hook :
hooks_) {
125 test |= hook->canVetoProcessLevel();
134 for (Pythia8::UserHooks *hook :
hooks_) {
135 if (hook->canVetoProcessLevel()) test |= hook->doVetoProcessLevel(event);
143 for (Pythia8::UserHooks *hook :
hooks_) {
144 test |= hook->canVetoResonanceDecays();
153 for (Pythia8::UserHooks *hook :
hooks_) {
154 if (hook->canVetoResonanceDecays()) test |= hook->doVetoResonanceDecays(event);
163 for (Pythia8::UserHooks *hook :
hooks_) {
164 test |= hook->canVetoPT();
171 for (Pythia8::UserHooks *hook :
hooks_) {
172 if (hook->canVetoPT())
return hook->scaleVetoPT();
182 for (Pythia8::UserHooks *hook :
hooks_) {
183 if (hook->canVetoPT())
return hook->doVetoPT(iPos,event);
192 for (Pythia8::UserHooks *hook :
hooks_) {
193 test |= hook->canVetoStep();
200 for (Pythia8::UserHooks *hook :
hooks_) {
201 if (hook->canVetoStep())
return hook->numberVetoStep();
209 bool doVetoStep(
int iPos,
int nISR,
int nFSR,
const Pythia8::Event&
event)
override {
210 for (Pythia8::UserHooks *hook :
hooks_) {
211 if (hook->canVetoStep())
return hook->doVetoStep(iPos,nISR,nFSR,event);
220 for (Pythia8::UserHooks *hook :
hooks_) {
221 test |= hook->canVetoMPIStep();
228 for (Pythia8::UserHooks *hook :
hooks_) {
229 if (hook->canVetoMPIStep())
return hook->numberVetoMPIStep();
237 for (Pythia8::UserHooks *hook :
hooks_) {
238 if (hook->canVetoMPIStep())
return hook->doVetoMPIStep(nMPI,event);
247 for (Pythia8::UserHooks *hook :
hooks_) {
248 test |= hook->canVetoPartonLevelEarly();
257 for (Pythia8::UserHooks *hook :
hooks_) {
258 if (hook->canVetoPartonLevelEarly()) test |= hook->doVetoPartonLevelEarly(event);
268 for (Pythia8::UserHooks *hook :
hooks_) {
269 test |= hook->retryPartonLevel();
277 for (Pythia8::UserHooks *hook :
hooks_) {
278 test |= hook->canVetoPartonLevel();
287 for (Pythia8::UserHooks *hook :
hooks_) {
288 if (hook->canVetoPartonLevel()) test |= hook->doVetoPartonLevel(event);
296 for (Pythia8::UserHooks *hook :
hooks_) {
297 test |= hook->canSetResonanceScale();
307 for (Pythia8::UserHooks *hook :
hooks_) {
308 if (hook->canSetResonanceScale())
return hook->scaleResonance(iRes,event);
316 for (Pythia8::UserHooks *hook :
hooks_) {
317 test |= hook->canVetoISREmission();
328 for (Pythia8::UserHooks *hook :
hooks_) {
329 if (hook->canVetoISREmission()) test |= hook->doVetoISREmission(sizeOld,event,iSys);
337 for (Pythia8::UserHooks *hook :
hooks_) {
338 test |= hook->canVetoFSREmission();
351 for (Pythia8::UserHooks *hook :
hooks_) {
352 if (hook->canVetoFSREmission()) test |= hook->doVetoFSREmission(sizeOld,event,iSys,inResonance);
360 for (Pythia8::UserHooks *hook :
hooks_) {
361 test |= hook->canVetoMPIEmission();
371 for (Pythia8::UserHooks *hook :
hooks_) {
372 if (hook->canVetoMPIEmission()) test |= hook->doVetoMPIEmission(sizeOld,event);
380 for (Pythia8::UserHooks *hook :
hooks_) {
381 test |= hook->canReconnectResonanceSystems();
393 for (Pythia8::UserHooks *hook :
hooks_) {
394 if (hook->canReconnectResonanceSystems()) test &= hook->doReconnectResonanceSystems(oldSizeEvt,event);
402 for (Pythia8::UserHooks *hook :
hooks_) {
403 test |= hook->canEnhanceEmission();
412 for (Pythia8::UserHooks *hook :
hooks_) {
413 test |= hook->canEnhanceTrial();
419 for (Pythia8::UserHooks *hook :
hooks_) {
420 if (hook->canEnhanceEmission() || hook->canEnhanceTrial())
return hook->enhanceFactor(str);
426 for (Pythia8::UserHooks *hook :
hooks_) {
427 if (hook->canEnhanceEmission() || hook->canEnhanceTrial())
return hook->vetoProbability(str);
435 for (Pythia8::UserHooks *hook :
hooks_) {
436 test |= hook->canChangeFragPar();
443 bool doChangeFragPar( Pythia8::StringFlav* flavPtr, Pythia8::StringZ* zPtr, Pythia8::StringPT* pTPtr,
int idEnd,
444 double m2Had, std::vector<int> iParton,
const Pythia8::StringEnd* SE)
override {
446 for (Pythia8::UserHooks *hook :
hooks_) {
447 if (hook->canChangeFragPar()) test &= hook->doChangeFragPar(flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton, SE);
455 for (Pythia8::UserHooks *hook :
hooks_) {
456 if (hook->canChangeFragPar()) test |= hook->doVetoFragmentation(part, SE);
double vetoProbability(std::string str) override
bool doVetoPartonLevel(const Pythia8::Event &event) override
std::vector< Pythia8::UserHooks * > hooks() const
bool doVetoFSREmission(int sizeOld, const Pythia8::Event &event, int iSys, bool inResonance=false) override
std::vector< Pythia8::UserHooks * > hooks_
double scaleResonance(int iRes, const Pythia8::Event &event) override
double biasedSelectionWeight() override
bool canVetoMPIStep() override
void addHook(Pythia8::UserHooks *hook)
bool initAfterBeams() override
bool doVetoStep(int iPos, int nISR, int nFSR, const Pythia8::Event &event) override
bool canVetoISREmission() override
bool canVetoPT() override
bool canEnhanceTrial() override
bool canSetResonanceScale() override
bool canVetoFSREmission() override
bool doVetoPT(int iPos, const Pythia8::Event &event) override
bool canVetoPartonLevel() override
bool doVetoFragmentation(Pythia8::Particle part, const Pythia8::StringEnd *SE) override
bool canVetoPartonLevelEarly() override
bool canVetoProcessLevel() override
bool canModifySigma() override
bool canVetoResonanceDecays() override
bool doVetoProcessLevel(Pythia8::Event &event) override
bool doVetoMPIStep(int nMPI, const Pythia8::Event &event) override
bool doVetoISREmission(int sizeOld, const Pythia8::Event &event, int iSys) override
bool canChangeFragPar() override
bool doChangeFragPar(Pythia8::StringFlav *flavPtr, Pythia8::StringZ *zPtr, Pythia8::StringPT *pTPtr, int idEnd, double m2Had, std::vector< int > iParton, const Pythia8::StringEnd *SE) override
double enhanceFactor(std::string str) override
double scaleVetoPT() override
int numberVetoMPIStep() override
double biasSelectionBy(const Pythia8::SigmaProcess *sigmaProcessPtr, const Pythia8::PhaseSpace *phaseSpacePtr, bool inEvent) override
bool canVetoStep() override
bool doReconnectResonanceSystems(int oldSizeEvt, Pythia8::Event &event) override
bool canReconnectResonanceSystems() override
bool doVetoMPIEmission(int sizeOld, const Pythia8::Event &event) override
bool doVetoResonanceDecays(Pythia8::Event &event) override
int numberVetoStep() override
bool canVetoMPIEmission() override
unsigned int nHooks() const
bool retryPartonLevel() override
bool canBiasSelection() override
bool canEnhanceEmission() override
double multiplySigmaBy(const Pythia8::SigmaProcess *sigmaProcessPtr, const Pythia8::PhaseSpace *phaseSpacePtr, bool inEvent) override
bool doVetoPartonLevelEarly(const Pythia8::Event &event) override