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