CMS 3D CMS Logo

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