CMS 3D CMS Logo

MultiUserHook.h
Go to the documentation of this file.
1 //to allow combining multiple user hooks
2 
3 #include "Pythia8/UserHooks.h"
4 #include "Pythia8/StringFragmentation.h"
5 
6 class MultiUserHook : public Pythia8::UserHooks {
7 
8 public:
9 
10  // Constructor and destructor.
12 
13  std::vector<Pythia8::UserHooks*> hooks() const { return hooks_; }
14  unsigned int nHooks() const { return hooks_.size(); }
15  void addHook(Pythia8::UserHooks *hook) { hooks_.push_back(hook); }
16 
17 // Initialisation after beams have been set by Pythia::init().
18  bool initAfterBeams() override {
19  bool test = true;
20  for (Pythia8::UserHooks *hook : hooks_) {
21  hook->initPtr(infoPtr, settingsPtr,
22  particleDataPtr, rndmPtr,
23  beamAPtr, beamBPtr,
24  beamPomAPtr, beamPomBPtr,
25  coupSMPtr, partonSystemsPtr,
26  sigmaTotPtr);
27  test &= hook->initAfterBeams();
28  }
29 
30  //Certain types of hooks can only be handled by one UserHook at a time.
31  //Check for this and return an error if needed
32  int nCanVetoPT = 0;
33  int nCanVetoStep = 0;
34  int nCanVetoMPIStep = 0;
35  int nCanSetResonanceScale = 0;
36  int nCanEnhance = 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;
43  }
44 
45  if (nCanVetoPT>1) {
46  infoPtr->errorMsg("Error in MultiUserHook::initAfterBeams "
47  "multiple UserHooks with canVetoPT() not allowed");
48  test = false;
49  }
50  if (nCanVetoStep>1) {
51  infoPtr->errorMsg("Error in MultiUserHook::initAfterBeams "
52  "multiple UserHooks with canVetoStep() not allowed");
53  test = false;
54  }
55  if (nCanVetoMPIStep>1) {
56  infoPtr->errorMsg("Error in MultiUserHook::initAfterBeams "
57  "multiple UserHooks with canVetoMPIStep() not allowed");
58  test = false;
59  }
60  if (nCanSetResonanceScale>1) {
61  infoPtr->errorMsg("Error in MultiUserHook::initAfterBeams "
62  "multiple UserHooks with canSetResonanceScale() not allowed");
63  test = false;
64  }
65  if (nCanEnhance>1) {
66  infoPtr->errorMsg("Error in MultiUserHook::initAfterBeams "
67  "multiple UserHooks with canEnhanceEmission() or canEnhanceTrial() not allowed");
68  test = false;
69  }
70 
71  return test;
72  }
73 
74  // Possibility to modify cross section of process.
75  bool canModifySigma() override {
76  bool test = false;
77  for (Pythia8::UserHooks *hook : hooks_) {
78  test |= hook->canModifySigma();
79  }
80  return test;
81  }
82 
83  // Multiplicative factor modifying the cross section of a hard process.
84  double multiplySigmaBy(const Pythia8::SigmaProcess* sigmaProcessPtr,
85  const Pythia8::PhaseSpace* phaseSpacePtr, bool inEvent) override {
86  double factor = 1.;
87  for (Pythia8::UserHooks *hook : hooks_) {
88  if (hook->canModifySigma()) factor *= hook->multiplySigmaBy(sigmaProcessPtr,phaseSpacePtr,inEvent);
89  }
90  return factor;
91  }
92 
93  // Possibility to bias selection of events, compensated by a weight.
94  bool canBiasSelection() override {
95  bool test = false;
96  for (Pythia8::UserHooks *hook : hooks_) {
97  test |= hook->canBiasSelection();
98  }
99  return test;
100  }
101 
102  // Multiplicative factor in the phase space selection of a hard process.
103  double biasSelectionBy(const Pythia8::SigmaProcess* sigmaProcessPtr,
104  const Pythia8::PhaseSpace* phaseSpacePtr, bool inEvent) override {
105  double factor = 1.;
106  for (Pythia8::UserHooks *hook : hooks_) {
107  if (hook->canBiasSelection()) factor *= hook->biasSelectionBy(sigmaProcessPtr,phaseSpacePtr,inEvent);
108  }
109  return factor;
110  };
111 
112  // Event weight to compensate for selection weight above.
113  double biasedSelectionWeight() override {
114  double factor = 1.;
115  for (Pythia8::UserHooks *hook : hooks_) {
116  if (hook->canBiasSelection()) factor *= hook->biasedSelectionWeight();
117  }
118  return factor;
119  };
120 
121  // Possibility to veto event after process-level selection.
122  bool canVetoProcessLevel() override {
123  bool test = false;
124  for (Pythia8::UserHooks *hook : hooks_) {
125  test |= hook->canVetoProcessLevel();
126  }
127  return test;
128  }
129 
130  // Decide whether to veto current process or not, based on process record.
131  // Usage: doVetoProcessLevel( process).
132  bool doVetoProcessLevel(Pythia8::Event& event) override {
133  bool test = false;
134  for (Pythia8::UserHooks *hook : hooks_) {
135  if (hook->canVetoProcessLevel()) test |= hook->doVetoProcessLevel(event);
136  }
137  return test;
138  }
139 
140  // Possibility to veto resonance decay chain.
141  bool canVetoResonanceDecays() override {
142  bool test = false;
143  for (Pythia8::UserHooks *hook : hooks_) {
144  test |= hook->canVetoResonanceDecays();
145  }
146  return test;
147  }
148 
149  // Decide whether to veto current resonance decay chain or not, based on
150  // process record. Usage: doVetoProcessLevel( process).
151  bool doVetoResonanceDecays(Pythia8::Event& event) override {
152  bool test = false;
153  for (Pythia8::UserHooks *hook : hooks_) {
154  if (hook->canVetoResonanceDecays()) test |= hook->doVetoResonanceDecays(event);
155  }
156  return test;
157  }
158 
159  // Possibility to veto MPI + ISR + FSR evolution and kill event,
160  // making decision at a fixed pT scale. Useful for MLM-style matching.
161  bool canVetoPT() override {
162  bool test = false;
163  for (Pythia8::UserHooks *hook : hooks_) {
164  test |= hook->canVetoPT();
165  }
166  return test;
167  }
168 
169  // Transverse-momentum scale for veto test. (only one hook allowed)
170  double scaleVetoPT() override {
171  for (Pythia8::UserHooks *hook : hooks_) {
172  if (hook->canVetoPT()) return hook->scaleVetoPT();
173  }
174  return 0.;
175  };
176 
177  // Decide whether to veto current event or not, based on event record.
178  // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
179  // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
180  // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
181  bool doVetoPT( int iPos, const Pythia8::Event& event) override {
182  for (Pythia8::UserHooks *hook : hooks_) {
183  if (hook->canVetoPT()) return hook->doVetoPT(iPos,event);
184  }
185  return false;
186  }
187 
188  // Possibility to veto MPI + ISR + FSR evolution and kill event,
189  // making decision after fixed number of ISR or FSR steps.
190  bool canVetoStep() override {
191  bool test = false;
192  for (Pythia8::UserHooks *hook : hooks_) {
193  test |= hook->canVetoStep();
194  }
195  return test;
196  }
197 
198  // Up to how many ISR + FSR steps of hardest interaction should be checked.
199  int numberVetoStep() override {
200  for (Pythia8::UserHooks *hook : hooks_) {
201  if (hook->canVetoStep()) return hook->numberVetoStep();
202  }
203  return 1;
204  };
205 
206  // Decide whether to veto current event or not, based on event record.
207  // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
208  // nISR and nFSR number of emissions so far for hard interaction only.
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);
212  }
213  return false;
214  }
215 
216  // Possibility to veto MPI + ISR + FSR evolution and kill event,
217  // making decision after fixed number of MPI steps.
218  bool canVetoMPIStep() override {
219  bool test = false;
220  for (Pythia8::UserHooks *hook : hooks_) {
221  test |= hook->canVetoMPIStep();
222  }
223  return test;
224  }
225 
226  // Up to how many MPI steps should be checked.
227  int numberVetoMPIStep() override {
228  for (Pythia8::UserHooks *hook : hooks_) {
229  if (hook->canVetoMPIStep()) return hook->numberVetoMPIStep();
230  }
231  return 1;
232  };
233 
234  // Decide whether to veto current event or not, based on event record.
235  // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
236  bool doVetoMPIStep( int nMPI, const Pythia8::Event& event) override {
237  for (Pythia8::UserHooks *hook : hooks_) {
238  if (hook->canVetoMPIStep()) return hook->doVetoMPIStep(nMPI,event);
239  }
240  return false;
241  }
242 
243  // Possibility to veto event after ISR + FSR + MPI in parton level,
244  // but before beam remnants and resonance decays.
245  bool canVetoPartonLevelEarly() override {
246  bool test = false;
247  for (Pythia8::UserHooks *hook : hooks_) {
248  test |= hook->canVetoPartonLevelEarly();
249  }
250  return test;
251  }
252 
253  // Decide whether to veto current partons or not, based on event record.
254  // Usage: doVetoPartonLevelEarly( event).
255  bool doVetoPartonLevelEarly( const Pythia8::Event& event) override {
256  bool test = false;
257  for (Pythia8::UserHooks *hook : hooks_) {
258  if (hook->canVetoPartonLevelEarly()) test |= hook->doVetoPartonLevelEarly(event);
259  }
260  return test;
261  }
262 
263  // Retry same ProcessLevel with a new PartonLevel after a veto in
264  // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
265  // if you overload this method to return true.
266  bool retryPartonLevel() override {
267  bool test = false;
268  for (Pythia8::UserHooks *hook : hooks_) {
269  test |= hook->retryPartonLevel();
270  }
271  return test;
272  }
273 
274  // Possibility to veto event after parton-level selection.
275  bool canVetoPartonLevel() override {
276  bool test = false;
277  for (Pythia8::UserHooks *hook : hooks_) {
278  test |= hook->canVetoPartonLevel();
279  }
280  return test;
281  }
282 
283  // Decide whether to veto current partons or not, based on event record.
284  // Usage: doVetoPartonLevel( event).
285  bool doVetoPartonLevel( const Pythia8::Event& event) override {
286  bool test = false;
287  for (Pythia8::UserHooks *hook : hooks_) {
288  if (hook->canVetoPartonLevel()) test |= hook->doVetoPartonLevel(event);
289  }
290  return test;
291  }
292 
293  // Possibility to set initial scale in TimeShower for resonance decay.
294  bool canSetResonanceScale() override {
295  bool test = false;
296  for (Pythia8::UserHooks *hook : hooks_) {
297  test |= hook->canSetResonanceScale();
298  }
299  return test;
300  }
301 
302  // Initial scale for TimeShower evolution.
303  // Usage: scaleResonance( iRes, event), where iRes is location
304  // of decaying resonance in the event record.
305  // Only one UserHook setting the resonance scale is allowed
306  double scaleResonance( int iRes, const Pythia8::Event& event) override {
307  for (Pythia8::UserHooks *hook : hooks_) {
308  if (hook->canSetResonanceScale()) return hook->scaleResonance(iRes,event);
309  }
310  return 0.;
311  };
312 
313  // Possibility to veto an emission in the ISR machinery.
314  bool canVetoISREmission() override {
315  bool test = false;
316  for (Pythia8::UserHooks *hook : hooks_) {
317  test |= hook->canVetoISREmission();
318  }
319  return test;
320  }
321 
322  // Decide whether to veto current emission or not, based on event record.
323  // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
324  // of event record before current emission-to-be-scrutinized was added,
325  // and iSys is the system of the radiation (according to PartonSystems).
326  bool doVetoISREmission( int sizeOld, const Pythia8::Event& event, int iSys) override {
327  bool test = false;
328  for (Pythia8::UserHooks *hook : hooks_) {
329  if (hook->canVetoISREmission()) test |= hook->doVetoISREmission(sizeOld,event,iSys);
330  }
331  return test;
332  }
333 
334  // Possibility to veto an emission in the FSR machinery.
335  bool canVetoFSREmission() override {
336  bool test = false;
337  for (Pythia8::UserHooks *hook : hooks_) {
338  test |= hook->canVetoFSREmission();
339  }
340  return test;
341  }
342 
343  // Decide whether to veto current emission or not, based on event record.
344  // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
345  // sizeOld is size of event record before current emission-to-be-scrutinized
346  // was added, iSys is the system of the radiation (according to
347  // PartonSystems), and inResonance is true if the emission takes place in a
348  // resonance decay.
349  bool doVetoFSREmission( int sizeOld, const Pythia8::Event& event, int iSys, bool inResonance = false ) override {
350  bool test = false;
351  for (Pythia8::UserHooks *hook : hooks_) {
352  if (hook->canVetoFSREmission()) test |= hook->doVetoFSREmission(sizeOld,event,iSys,inResonance);
353  }
354  return test;
355  }
356 
357  // Possibility to veto an MPI.
358  bool canVetoMPIEmission() override {
359  bool test = false;
360  for (Pythia8::UserHooks *hook : hooks_) {
361  test |= hook->canVetoMPIEmission();
362  }
363  return test;
364  }
365 
366  // Decide whether to veto an MPI based on event record.
367  // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
368  // is size of event record before the current MPI.
369  bool doVetoMPIEmission( int sizeOld, const Pythia8::Event &event) override {
370  bool test = false;
371  for (Pythia8::UserHooks *hook : hooks_) {
372  if (hook->canVetoMPIEmission()) test |= hook->doVetoMPIEmission(sizeOld,event);
373  }
374  return test;
375  }
376 
377  // Possibility to reconnect colours from resonance decay systems.
379  bool test = false;
380  for (Pythia8::UserHooks *hook : hooks_) {
381  test |= hook->canReconnectResonanceSystems();
382  }
383  return test;
384  }
385 
386  // Do reconnect colours from resonance decay systems.
387  // Usage: doVetoFSREmission( oldSizeEvt, event)
388  // where oldSizeEvent is the event size before resonance decays.
389  // Should normally return true, while false means serious failure.
390  // Value of PartonLevel:earlyResDec determines where method is called.
391  bool doReconnectResonanceSystems( int oldSizeEvt, Pythia8::Event &event) override {
392  bool test = true;
393  for (Pythia8::UserHooks *hook : hooks_) {
394  if (hook->canReconnectResonanceSystems()) test &= hook->doReconnectResonanceSystems(oldSizeEvt,event);
395  }
396  return test;
397  }
398 
399  // Enhance emission rates (sec. 4 in EPJC (2013) 73).
400  bool canEnhanceEmission() override {
401  bool test = false;
402  for (Pythia8::UserHooks *hook : hooks_) {
403  test |= hook->canEnhanceEmission();
404  }
405  return test;
406  }
407 
408  // Bookkeeping of weights for enhanced actual or trial emissions
409  // (sec. 3 in EPJC (2013) 73).
410  bool canEnhanceTrial() override {
411  bool test = false;
412  for (Pythia8::UserHooks *hook : hooks_) {
413  test |= hook->canEnhanceTrial();
414  }
415  return test;
416  }
417 
418  double enhanceFactor( std::string str) override {
419  for (Pythia8::UserHooks *hook : hooks_) {
420  if (hook->canEnhanceEmission() || hook->canEnhanceTrial()) return hook->enhanceFactor(str);
421  }
422  return 1.;
423  };
424 
425  double vetoProbability( std::string str ) override {
426  for (Pythia8::UserHooks *hook : hooks_) {
427  if (hook->canEnhanceEmission() || hook->canEnhanceTrial()) return hook->vetoProbability(str);
428  }
429  return 0.;
430  };
431 
432  // Can change fragmentation parameters.
433  bool canChangeFragPar() override {
434  bool test = false;
435  for (Pythia8::UserHooks *hook : hooks_) {
436  test |= hook->canChangeFragPar();
437  }
438  return test;
439  }
440 
441  // Do change fragmentation parameters.
442  // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton.
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 {
445  bool test = true;
446  for (Pythia8::UserHooks *hook : hooks_) {
447  if (hook->canChangeFragPar()) test &= hook->doChangeFragPar(flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton, SE);
448  }
449  return test;
450  }
451 
452  // Do a veto on a hadron just before it is added to the final state.
453  bool doVetoFragmentation( Pythia8::Particle part, const Pythia8::StringEnd* SE) override {
454  bool test = false;
455  for (Pythia8::UserHooks *hook : hooks_) {
456  if (hook->canChangeFragPar()) test |= hook->doVetoFragmentation(part, SE);
457  }
458  return test;
459  }
460 
461 
462 
463 //--------------------------------------------------------------------------
464 
465 private:
466  std::vector<Pythia8::UserHooks*> hooks_;
467 
468 };
double vetoProbability(std::string str) override
bool doVetoPartonLevel(const Pythia8::Event &event) override
std::vector< Pythia8::UserHooks * > hooks() const
Definition: MultiUserHook.h:13
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)
Definition: MultiUserHook.h:15
bool initAfterBeams() override
Definition: MultiUserHook.h:18
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
Definition: MultiUserHook.h:75
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
part
Definition: HCALResponse.h:20
bool doVetoResonanceDecays(Pythia8::Event &event) override
int numberVetoStep() override
bool canVetoMPIEmission() override
unsigned int nHooks() const
Definition: MultiUserHook.h:14
bool retryPartonLevel() override
bool canBiasSelection() override
Definition: MultiUserHook.h:94
#define str(s)
bool canEnhanceEmission() override
double multiplySigmaBy(const Pythia8::SigmaProcess *sigmaProcessPtr, const Pythia8::PhaseSpace *phaseSpacePtr, bool inEvent) override
Definition: MultiUserHook.h:84
bool doVetoPartonLevelEarly(const Pythia8::Event &event) override
Definition: event.py:1