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