ATLAS Offline Software
Loading...
Searching...
No Matches
PowhegHooksBB4Ldlsl.cxx
Go to the documentation of this file.
1// PowhegHooksBB4L.h
2// Rewritten by T. Jezo in 2021. With various contributions from S. Ferrario
3// Ravasio, B. Nachman, P. Nason and M. Seidel. Inspired by
4// ttb_NLO_dec/main-PYTHIA8.f by P. Nason and E. Re and by PowhegHooks.h by R.
5// Corke.
6//
7// Adapted to be Athena and Pythia8_i compliant by Katharina Voß (katharina.voss@cern.ch)
8// following the example from Simone Amoroso in PowhegBB4Ltms.cxx
9
10// # Introduction
11//
12// This hook is intended for use together with POWHEG-BOX-RES/b_bbar_4l NLO LHE
13// events. This also includes events in which one of the W bosons was
14// re-decayed hadronically. (Note that LHE format version larger than 3 may not
15// require this hook).
16//
17// The hook inherits from PowhegHooks and as such it indirectly implements the
18// following:
19// - doVetoMPIStep
20// - doVetoISREmission
21// - doVetoMPIEmission
22// and it overloads:
23// - doVetoFSREmission, which works as follows (if POWHEG:veto==1):
24// - if inResonance==true it vetoes all the emission that is harder than
25// the scale of its parent (anti-)top quark or W^+(-)
26// - if inResonance==false, it calls PowhegHooks::doVetoISREmission
27// and it also implements:
28// - doVetoProcessLevel, which is never used for vetoing (i.e. it always
29// returns false). Instead it is used for the calculation of reconance scales
30// using LHE kinematics.
31//
32// This version of the hooks is only suitable for use with fully compatible
33// POWHEG-BOX Les Houches readers (such the one in main-PYTHIA82-lhef but
34// not the one in main31.cc.)
35//
36//
37// # Basic use
38//
39// In order to use this hook you must replace all the declarations and
40// constructor calls of PowhegHooks to PowhegHooksBB4L:
41//
42// PowhegHooks *powhegHooks; -> PowhegHooksBB4L *powhegHooks;
43// *powhegHooks = new PowhegHooks(); -> *powhegHooks = new PowhegHooksBB4L();
44//
45// In order to switch it on set POWHEG:veto = 1 and
46// POWHEG:bb4l:FSREmission:veto = 1. This will lead to a veto in ISR, FSR and
47// MPI steps of pythia as expected using PowhegHooks in all the cases other than
48// the case of FSR emission in resonance decay. Within resonance decays
49// PowhegHooksBB4L takes over the control.
50//
51// Furthermore, this hook can also be used standalone without PowhegHooks, i.e.
52// only the FSR emission from resonances will be vetoed (the default use in
53// 1801.03944 and 1906.09166). In order to do that set
54// POWHEG:bb4l:FSREmission:veto = 1 and POWHEG:veto = 0. Note that the this is
55// not recommended for comparing against data but because it is easier to
56// interpret it is often done in theoretical studies.
57//
58// Note that this version of the hook relies on the event "radtype" (1 for
59// btilde, 2 for remnant) to be set by an external program, such as
60// main-PYTHIA82-lhef in the radtype_ common block.
61// There also exists a version of this hook in which the event "radtype" is
62// read in from the .lhe file using pythia built in functionality. You need
63// that version if you want to use this hook with main31.cc.
64//
65//
66// # Expert use
67//
68// This hook also implements an alternative veto procedure which allows to
69// assign a "SCALUP" type of scale to a resonance using the scaleResonance
70// method. This is a much simpler veto but it is also clearly inferior as
71// compared to the one implemented using the doVetoFSREmission method because
72// the definition of the scale of the emission does not match the
73// POWHEG-BOX-RES definition. Nevertheless, it can be activated using
74// POWHEG:bb4l:ScaleResonance:veto = 1. Additionally one MUST switch off the
75// other veto by calling on POWHEG:bb4l:FSREmission:veto = 0.
76//
77// The following settings are at the disposal of the user to control the
78// behaviour of the hook
79// - On/off switches for the veto:
80// - POWHEG:bb4l:FSREmission:veto
81// on/off switch for the default veto based on doFSREmission
82// - POWHEG:bb4l:ScaleResonance:veto
83// on/off switch for the alternative veto based on scaleResonance (only
84// for expert users)
85// - Important settings:
86// - POWHEG:bb4l:ptMinVeto: MUST be set to the same value as the
87// corresponding flag in POWHEG-BOX-RES
88// - Alternatives for scale calculations
89// - default: emission scale is calculated using POWHEG definitions and in
90// the resonance rest frame
91// - POWHEG:bb4l:FSREmission:vetoDipoleFrame: emission scale is calculated
92// using POWHEG definitions in the dipole frame
93// - POWHEG:bb4l:FSREmission:pTpythiaVeto: emission scale is calculated
94// using Pythia definitions
95// - Other flags:
96// - POWHEG:bb4l:FSREmission:vetoQED: decides whether or not QED emission
97// off quarks should also be vetoed (not implemented in the case of
98// the ScaleResonance:veto)
99// - POWHEG:bb4l:DEBUG: enables debug printouts on standard output
100
101#include "Pythia8/Pythia.h"
102#include "UserHooksUtils.h"
104#if PYTHIA_VERSION_INTEGER >= 8310
105// Disable the pythia 8.310 plugin mechanism.
106// We use our own mechanism to create these classes, and we get link errors
107// if these are used in more than one compilation unit.
108# include "Pythia8/Plugins.h"
109# undef PYTHIA8_PLUGIN_CLASS
110# undef PYTHIA8_PLUGIN_VERSIONS
111# define PYTHIA8_PLUGIN_CLASS(BASE, CLASS, PYTHIA, SETTINGS, LOGGER)
112# define PYTHIA8_PLUGIN_VERSIONS(...)
113#endif
114
116
117#include "UserSetting.h"
118#include <iostream>
119#include <cassert>
120#include <stdexcept>
121
122struct{
125
126
127namespace Pythia8 { class PowhegBB4Ldlsl; }
129
130
131namespace Pythia8 {
132 using namespace std;
133
135
136 public:
137
139 : m_debug("Powheg:bb4l:DEBUG", false),
140 m_vetoFSREmission("Powheg:bb4l:FSREmission:veto", true),
141 m_vetoQED("Powheg:bb4l:vetoQED", false),
142 m_topresscale(-1.), m_atopresscale(-1.),
143 m_wpresscale(-1.),m_wmresscale(-1.),
145 m_pTmin("Powheg:bb4l:pTminVeto", 0.5),
146 m_vetoProduction("Powheg:veto", 1),
147 m_pTpythiaVeto("Powheg:bb4l:pTpythiaVeto", 0),
148 m_vetoDipoleFrame("Powheg:bb4l:FSREmission:vetoDipoleFrame", 0),
149 m_scaleResonanceVeto("Powheg:bb4l:ScaleResonance:veto", 0.),
150 m_pThardMode("POWHEG:pThard",0)
151 {
152 std::cout << "**********************************************************" << std::endl;
153 std::cout << "* *" << std::endl;
154 std::cout << "* Applying Powheg BB4L UserHook (PowhegBB4Ldlsl)! *" << std::endl;
155 std::cout << "* Must run on dedicated Powheg LHE input file *" << std::endl;
156 std::cout << "* (on your own responsibility) *" << std::endl;
157 std::cout << "* *" << std::endl;
158 std::cout << "**********************************************************" << std::endl;
159
160 }
161
162 ~PowhegBB4Ldlsl() override { }
163
164 //--- Initialization -------------------------------------------------------
165 virtual bool initAfterBeams() override {
166 // initialize settings of the parent class
168 // initialize settings of this class already initialized in constructor
169 return true;
170 }
171
172 //--- PROCESS LEVEL HOOK ---------------------------------------------------
173 // This hook gets triggered for each event before the shower starts, i.e. at
174 // the LHE level. We use it to calculate the scales of resonances.
175 inline virtual bool canVetoProcessLevel() override { return true; }
176 inline virtual bool doVetoProcessLevel(Event &e) override {
177
178 // extract the radtype from the event comment
179 stringstream ss;
180 ss << infoPtr->getEventComments();
181 string temp;
182 ss >> temp >> radtype_.radtype;
183 assert (temp == "#rwgt");
184
185 // we only calculate resonance scales for btilde events (radtype_.radtype==1)
186 // remnant events are not vetoed
187 if (radtype_.radtype==2) return false;
188 // find last top and the last anti-top in the record
189 int i_top = -1, i_atop = -1, i_wp = -1, i_wm = -1;
190 for (int i = 0; i < e.size(); i++) {
191 if (e[i].id() == 6) i_top = i;
192 if (e[i].id() == -6) i_atop = i;
193 if (e[i].id() == 24) i_wp = i;
194 if (e[i].id() == -24) i_wm = i;
195 }
196 // if found calculate the resonance scale
197 m_topresscale = findresscale(i_top, e);
198 // similarly for anti-top
199 m_atopresscale = findresscale(i_atop, e);
200 // and for W^+ and W^-
201 m_wpresscale = findresscale(i_wp, e);
202 m_wmresscale = findresscale(i_wm, e);
203
204 // this is from old UserHook -> not needed anymore?
205 // // initialize stuff
206 // doVetoFSRInit();
207
208
209 // do not veto, ever
210 return false;
211 }
212
213 //--- FSR EMISSION LEVEL HOOK ----------------------------------------------
214 // This hook gets triggered everytime the parton shower attempts to attach
215 // a FSR emission.
216 inline virtual bool canVetoFSREmission() override { return m_vetoFSREmission(settingsPtr) || m_vetoProduction(settingsPtr); }
217 inline virtual bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance) override {
218
219 // FSR VETO INSIDE THE RESONANCE (if it is switched on)
220 if (inResonance && m_vetoFSREmission(settingsPtr)) {
221
222 // get the participants of the splitting: the recoiler, the radiator and the emitted
223 int iRecAft = e.size() - 1;
224 int iEmt = e.size() - 2;
225 int iRadAft = e.size() - 3;
226 int iRadBef = e[iEmt].mother1();
227
228 // find the resonance the radiator originates from
229 int iRes = e[iRadBef].mother1();
230 while ( iRes > 0 && (std::abs(e[iRes].id()) !=6 && std::abs(e[iRes].id()) != 24) ) {
231 iRes = e[iRes].mother1();
232 }
233 if (iRes == 0) {
234 #if PYTHIA_VERSION_INTEGER >= 8310
235 loggerPtr->ERROR_MSG("Warning in PowhegHooksBB4L::doVetoFSREmission: emission in resonance not from the top quark or from the W boson, not vetoing");
236 #else
237 infoPtr->errorMsg("Warning in PowhegHooksBB4L::doVetoFSREmission: emission in resonance not from the top quark or from the W boson, not vetoing");
238 #endif
239 return doVetoFSR(false);
240 }
241 int iResId = e[iRes].id();
242
243 // calculate the scale of the emission
244 double scale;
245 //using pythia pT definition ...
246 if(m_pTpythiaVeto(settingsPtr))
247 scale = pTpythia(e, iRadAft, iEmt, iRecAft);
248 //.. or using POWHEG pT definition
249 else{
250 Vec4 pr(e[iRadAft].p()), pe(e[iEmt].p()), pres(e[iRes].p()), prec(e[iRecAft].p()), psystem;
251 // The computation of the POWHEG pT can be done in the top rest frame or in the diple one.
252 // pdipole = pemt +prec +prad (after the emission)
253 // For the first emission off the top resonance pdipole = pw +pb (before the emission) = ptop
254 if(m_vetoDipoleFrame(settingsPtr))
255 psystem = pr+pe+prec;
256 else
257 psystem = pres;
258
259 // gluon splitting into two partons
260 if (e[iRadBef].id() == 21)
261 scale = gSplittingScale(psystem, pr, pe);
262 // quark emitting a gluon (or a photon)
263 else if (std::abs(e[iRadBef].id()) <= 5 && ((e[iEmt].id() == 21) && ! m_vetoQED(settingsPtr)) )
264 scale = qSplittingScale(psystem, pr, pe);
265 // other stuff (which we should not veto)
266 else {
267 scale = 0;
268 }
269 }
270
271 // compare the current splitting scale to the correct resonance scale
272 if (iResId == 6) {
273 if ( m_debug(settingsPtr) && scale > m_topresscale)
274 cout << iResId << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << endl;
275 return doVetoFSR(scale > m_topresscale);
276 }
277 else if (iResId == -6){
278 if ( m_debug(settingsPtr) && scale > m_atopresscale)
279 cout << iResId << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << endl;
280 return doVetoFSR(scale > m_atopresscale);
281 }
282 else if (iResId == 24) {
283 if ( m_debug(settingsPtr) && scale > m_wpresscale)
284 cout << iResId << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << endl;
285 return doVetoFSR(scale > m_wpresscale);
286 }
287 else if (iResId == -24){
288 if ( m_debug(settingsPtr) && scale > m_wmresscale)
289 cout << iResId << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << endl;
290 return doVetoFSR(scale > m_wmresscale);
291 }
292 else {
293 #if PYTHIA_VERSION_INTEGER >= 8310
294 loggerPtr->ERROR_MSG("Error in PowhegHooksBB4L::doVetoFSREmission: unimplemented case");
295 #else
296 infoPtr->errorMsg("Error in PowhegHooksBB4L::doVetoFSREmission: unimplemented case");
297 #endif
298 exit(-1);
299 }
300 }
301 // FSR VETO THE PRODUCTION PROCESS, i.e. OUTSIDE RESONANCE (if it is switched on)
302 else if(!inResonance && m_vetoProduction(settingsPtr)){
303 return PowhegHooksSetterMethod::doVetoFSREmission(sizeOld, e, iSys, inResonance);
304 }
305 // OTHERWISE DON'T VETO
306 else {
307 return false;
308 }
309 }
310
311 inline bool doVetoFSR(bool condition) {
312 if (radtype_.radtype==2) return false;
313 if (condition) {
315 return true;
316 }
317 return false;
318 }
319
320 //--- SCALE RESONANCE HOOK -------------------------------------------------
321 // called before each resonance decay shower
322 inline virtual bool canSetResonanceScale() override { return m_scaleResonanceVeto(settingsPtr); }
323 // if the resonance is the (anti)top or W+/W- set the scale to:
324 // - if radtype=2 (remnant): resonance virtuality
325 // - if radtype=1 (btilde):
326 // - (a)topresscale/wp(m)resscale for tops and Ws
327 // - a large number otherwise
328 // if is not the top, set it to a big number
329 inline virtual double scaleResonance(int iRes, const Event &e) override {
330 if(radtype_.radtype == 2)
331 return sqrt(e[iRes].m2Calc());
332 else {
333 if (e[iRes].id() == 6)
334 return m_topresscale;
335 else if (e[iRes].id() == -6)
336 return m_atopresscale;
337 else if (e[iRes].id() == 24)
338 return m_wpresscale;
339 else if (e[iRes].id() == -24)
340 return m_wmresscale;
341 else
342 return 1e30;
343 }
344 }
345
346 //--- Internal helper functions --------------------------------------------
347 // Calculates the scale of the hardest emission from within the resonance system
348 // translated by Markus Seidel modified by Tomas Jezo
349 inline double findresscale( const int iRes, const Event& event) {
350
351 // return large scale if the resonance position is ill defined
352 if (iRes < 0) return 1e30;
353
354 // get number of resonance decay products
355 int nDau = event[iRes].daughterList().size();
356
357 // iRes is not decayed, return high scale equivalent to
358 // unrestricted shower
359 if (nDau == 0) {
360 return 1e30;
361 }
362 // iRes did not radiate, this means that POWHEG pt scale has
363 // evolved all the way down to pTmin
364 else if (nDau < 3) {
365 return m_pTmin(settingsPtr);
366 }
367 // iRes is a (anti-)top quark
368 else if (std::abs(event[iRes].id()) == 6) {
369 // find top daughters
370 int idw = -1, idb = -1, idg = -1;
371 for (int i = 0; i < nDau; i++) {
372 int iDau = event[iRes].daughterList()[i];
373 if (std::abs(event[iDau].id()) == 24) idw = iDau;
374 if (std::abs(event[iDau].id()) == 5) idb = iDau;
375 if (std::abs(event[iDau].id()) == 21) idg = iDau;
376 }
377
378 // Get daughter 4-vectors in resonance frame
379 Vec4 pw(event[idw].p());
380 pw.bstback(event[iRes].p());
381 Vec4 pb(event[idb].p());
382 pb.bstback(event[iRes].p());
383 Vec4 pg(event[idg].p());
384 pg.bstback(event[iRes].p());
385
386 // Calculate scale and return it
387 return sqrt(2*pg*pb*pg.e()/pb.e());
388 }
389 // iRes is a W+(-) boson
390 else if (std::abs(event[iRes].id()) == 24) {
391 // Find W daughters
392 int idq = -1, ida = -1, idg = -1;
393 for (int i = 0; i < nDau; i++) {
394 int iDau = event[iRes].daughterList()[i];
395 if (event[iDau].id() == 21) idg = iDau;
396 else if (event[iDau].id() > 0) idq = iDau;
397 else if (event[iDau].id() < 0) ida = iDau;
398 }
399 if (idq<0 or ida <0){
400 throw std::out_of_range("idq or ida out of range in PowhegHooksBB4Ldlsl.cxx");
401 }
402 // Get daughter 4-vectors in resonance frame
403 Vec4 pq(event[idq].p());
404 pq.bstback(event[iRes].p());
405 Vec4 pa(event[ida].p());
406 pa.bstback(event[iRes].p());
407 Vec4 pg(event[idg].p());
408 pg.bstback(event[iRes].p());
409
410 // Calculate scale
411 Vec4 pw = pq + pa + pg;
412 double q2 = pw*pw;
413 double csi = 2*pg.e()/sqrt(q2);
414 double yq = 1 - pg*pq/(pg.e()*pq.e());
415 double ya = 1 - pg*pa/(pg.e()*pa.e());
416 // and return it
417 return sqrt(min(1-yq,1-ya)*pow2(csi)*q2/2);
418 }
419 // in any other case just return a high scale equivalent to
420 // unrestricted shower
421 return 1e30;
422 }
423
424 // The following routine will match daughters of particle `e[iparticle]` to an expected pattern specified via the list of expected particle PDG ID's `ids`,
425 // id wildcard is specified as 0 if match is obtained, the positions and the momenta of these particles are returned in vectors `positions` and `momenta`
426 // respectively
427 // if exitOnExtraLegs==true, it will exit if the decay has more particles than expected, but not less
428 inline bool match_decay(int iparticle, const Event &e, const vector<int> &ids, vector<int> &positions, vector<Vec4> &momenta, bool exitOnExtraLegs = true){
429 // compare sizes
430 if (e[iparticle].daughterList().size() != ids.size()) {
431 if (exitOnExtraLegs && e[iparticle].daughterList().size() > ids.size()) {
432 cout << "extra leg" << endl;
433 exit(-1);
434 }
435 return false;
436 }
437 // compare content
438 for (long unsigned int i = 0; i < e[iparticle].daughterList().size(); i++) {
439 int di = e[iparticle].daughterList()[i];
440 if (ids[i] != 0 && e[di].id() != ids[i])
441 return false;
442 }
443 // reset the positions and momenta vectors (because they may be reused)
444 positions.clear();
445 momenta.clear();
446 // construct the array of momenta
447 for (long unsigned int i = 0; i < e[iparticle].daughterList().size(); i++) {
448 int di = e[iparticle].daughterList()[i];
449 positions.push_back(di);
450 momenta.push_back(e[di].p());
451 }
452 return true;
453 }
454
455 inline double qSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2){
456 p1.bstback(pt);
457 p2.bstback(pt);
458 return sqrt( 2*p1*p2*p2.e()/p1.e() );
459 }
460
461 inline double gSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2){
462 p1.bstback(pt);
463 p2.bstback(pt);
464 return sqrt( 2*p1*p2*p1.e()*p2.e()/(pow(p1.e()+p2.e(),2)) );
465 }
466
467 // Routines to calculate the pT (according to pTdefMode) in a FS splitting:
468 // i (radiator before) -> j (emitted after) k (radiator after)
469 // For the Pythia pT definition, a recoiler (after) must be specified.
470 // (INSPIRED BY pythia8F77_31.cc double pTpythia)
471 inline double pTpythia(const Event &e, int RadAfterBranch, int EmtAfterBranch,
472 int RecAfterBranch)
473 {
474
475 // Convenient shorthands for later
476 Vec4 radVec = e[RadAfterBranch].p();
477 Vec4 emtVec = e[EmtAfterBranch].p();
478 Vec4 recVec = e[RecAfterBranch].p();
479 int radID = e[RadAfterBranch].id();
480
481 // Calculate virtuality of splitting
482 Vec4 Q(radVec + emtVec);
483 double Qsq = Q.m2Calc();
484
485 // Mass term of radiator
486 double m2Rad = (std::abs(radID) >= 4 && std::abs(radID) < 7) ?
487 pow2(particleDataPtr->m0(radID)) : 0.;
488
489 // z values for FSR
490 double z, pTnow;
491 // Construct 2 -> 3 variables
492 Vec4 sum = radVec + recVec + emtVec;
493 double m2Dip = sum.m2Calc();
494
495 double x1 = 2. * (sum * radVec) / m2Dip;
496 double x3 = 2. * (sum * emtVec) / m2Dip;
497 z = x1 / (x1 + x3);
498 pTnow = z * (1. - z);
499
500
501 // Virtuality
502 pTnow *= (Qsq - m2Rad);
503
504 if (pTnow < 0.) {
505 cout << "Warning: pTpythia was negative" << endl;
506 return -1.;
507 }
508 else
509 return(sqrt(pTnow));
510 }
511
512 // Functions to return statistics about the veto
514
515
516 //--------------------------------------------------------------------------
517
518 // Extraction of pThard based on the incoming event.
519 // Assume that all the final-state particles are in a continuous block
520 // at the end of the event and the final entry is the POWHEG emission.
521 // If there is no POWHEG emission, then pThard is set to SCALUP.
522
523 inline bool canVetoMPIStep() override { return true; }
524 inline int numberVetoMPIStep() override { return 1; }
525 inline bool doVetoMPIStep(int nMPI, const Event &e) override {
526 // let the original PowhegHook intialise all necessary variables
527 // and set shower starting scale for pThard=0 mode
528 // consequence: always need to set nFinal = -1 in job options!
530 // only reset pThard value if using pThardMode =1 or = 2
531 if (m_pThardMode(settingsPtr) == 1 || m_pThardMode(settingsPtr) == 2){
532 // Extra check on nMPI
533 if (nMPI > 1) return false;
534 int count = 0;
535 //check if there is a top and an antitop in the event
536 bool top_exists = false;
537 bool antitop_exists = false;
538 for (int i = e.size() - 1; i > 0; i--) {
539 if (e[i].isFinal()) {
540 count++;
541 //check if there is a top and an antitop in the event
542 if (e[i].id() == 6) top_exists = true;
543 if (e[i].id() == -6) antitop_exists = true;
544 } else break;
545 }
546
547 //adjust nFinal depending on the number of tops
548 int nCurrentFinal = -1;
549 // if top and antitop exist, then pp -> t tbar (+ pwhg emissions) -> WbWb (+X)
550 if (top_exists && antitop_exists) nCurrentFinal = 2;
551 // if only top or antitop exist, then pp -> t W b (+ pwhg emissions) -> WbWb (+X)
552 else if (top_exists || antitop_exists) nCurrentFinal = 3;
553 else{
554 std::cout << "Error: neither top nor anti-top found - something went wrong" << std::endl;
555 exit(1);
556 }
557 // other cases not considered currently
558 if (count != nCurrentFinal && count != nCurrentFinal + 1) {
559 cout << "Error: wrong number of final state particles in event: count " << count << " nCurrentFinal: "<< nCurrentFinal << endl;
560 for (int i = e.size() - 1; i > 0; i--) {
561 if (e[i].isFinal()) { std:: cout << "FS particle " << i << " PID: " << e[i].id() << std::endl;
562 }
563 }
564 exit(1);
565 }
566 // Flag if POWHEG radiation present and index
567 bool isEmt = (count == nCurrentFinal) ? false : true;
568 int iEmt = (isEmt) ? e.size() - 1 : -1;
569
570 double pThard = 0;
571 if (m_pThardMode(settingsPtr) == 1) {
572 pThard = PowhegHooksSetterMethod::pTcalc(e, -1, iEmt, -1, -1, -1);
573 // If pThardMode is 2, then the pT of all final-state partons is checked
574 // against all other incoming and outgoing partons, with the minimal value
575 // taken.
576 } else if (m_pThardMode(settingsPtr) == 2) {
577 pThard = PowhegHooksSetterMethod::pTcalc(e, -1, -1, -1, -1, -1);
578 } else {
579 std::cout << "Error: m_pThardMode neither 1 or 2 - something went wrong" << std::endl;
580 exit(1);
581 }
582
583 // set pThard value in base class UserHook
585 }
586
587 // Do not veto the event
588 return false;
589 }
590
591 //--------------------------------------------------------------------------
592
593 private:
597 unsigned long int m_nInResonanceFSRveto;
602
603 };
604
605 //==========================================================================
606
607} // end namespace Pythia8
608
609
static Double_t ss
struct @265362120060103244115120172162007233333034374360 radtype_
Pythia8_UserHooks::UserHooksFactory::Creator< Pythia8::PowhegBB4Ldlsl > powhegBB4LdlslCreator("PowhegBB4Ldlsl")
static const std::map< unsigned int, unsigned int > pow2
#define z
constexpr int pow(int base, int exp) noexcept
#define min(a, b)
Definition cfImp.cxx:40
virtual bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance) override
virtual bool canVetoProcessLevel() override
Pythia8_UserHooks::UserSetting< bool > m_vetoFSREmission
Pythia8_UserHooks::UserSetting< int > m_vetoProduction
double findresscale(const int iRes, const Event &event)
bool doVetoMPIStep(int nMPI, const Event &e) override
virtual double scaleResonance(int iRes, const Event &e) override
Pythia8_UserHooks::UserSetting< double > m_pTmin
bool doVetoFSR(bool condition)
virtual bool initAfterBeams() override
Pythia8_UserHooks::UserSetting< bool > m_vetoQED
bool match_decay(int iparticle, const Event &e, const vector< int > &ids, vector< int > &positions, vector< Vec4 > &momenta, bool exitOnExtraLegs=true)
virtual bool canVetoFSREmission() override
Pythia8_UserHooks::UserSetting< int > m_vetoDipoleFrame
virtual bool doVetoProcessLevel(Event &e) override
Pythia8_UserHooks::UserSetting< bool > m_debug
unsigned long int m_nInResonanceFSRveto
double pTpythia(const Event &e, int RadAfterBranch, int EmtAfterBranch, int RecAfterBranch)
virtual bool canSetResonanceScale() override
Pythia8_UserHooks::UserSetting< double > m_scaleResonanceVeto
Pythia8_UserHooks::UserSetting< int > m_pTpythiaVeto
double qSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2)
Pythia8_UserHooks::UserSetting< int > m_pThardMode
double gSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2)
double pTcalc(const Event &e, int i, int j, int k, int r, int xSRin)
bool doVetoFSREmission(int, const Event &e, int iSys, bool)
bool doVetoMPIStep(int nMPI, const Event &e)
int count(std::string s, const std::string &regx)
count how many occurances of a regx are in a string
Definition hcg.cxx:146
Author: James Monk (jmonk@cern.ch)
STL namespace.