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