ATLAS Offline Software
PowhegBB4Ltms.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // PowhegHooksBB4L.h
6 // Copyright (C) 2017 Silvia Ferrario Ravasio, Tomas Jezo, Paolo Nason, Markus Seidel
7 // based on PowhegHooks.h by Richard Corke
8 // Adapted to be Athena and Pythia8_i compliant by Simone Amoroso (amoroso@cern.ch)
9 
10 #include "Pythia8/Pythia.h"
11 #include "UserHooksUtils.h"
13 #if PYTHIA_VERSION_INTEGER >= 8310
14 // Disable the pythia 8.310 plugin mechanism.
15 // We use our own mechanism to create these classes, and we get link errors
16 // if these are used in more than one compilation unit.
17 # include "Pythia8/Plugins.h"
18 # undef PYTHIA8_PLUGIN_CLASS
19 # undef PYTHIA8_PLUGIN_VERSIONS
20 # define PYTHIA8_PLUGIN_CLASS(BASE, CLASS, PYTHIA, SETTINGS, LOGGER)
21 # define PYTHIA8_PLUGIN_VERSIONS(...)
22 #endif
23 #include "Pythia8Plugins/PowhegHooks.h"
24 #include "UserSetting.h"
25 #include <iostream>
26 #include <cassert>
27 #include <stdexcept>
28 
29 struct{
30  int radtype;
32 
33 
34 
35 
36 namespace Pythia8 { class PowhegBB4Ltms; }
38 
39 
40 
41 namespace Pythia8 {
42  using namespace std;
43 
44  //==========================================================================
45  // Use userhooks to veto PYTHIA emissions above the POWHEG scale.
46  class PowhegBB4Ltms : public PowhegHooks {
47  public:
48 
49  // Constructor
51  : m_debug("Powheg:bb4l:DEBUG", false),
52  m_vetoFSREmission("Powheg:bb4l:FSREmission:veto", true),
53  m_dryRunFSR("Powheg:bb4l:dryRunFSR", false),
54  m_onlyDistance1("Powheg:bb4l:onlyDistance1", false),
55  m_vetoAtPL("Powheg:bb4l:vetoAtPL", false),
56  m_vetoQED("Powheg:bb4l:vetoQED", false),
57  m_vetoPartonLevel("Powheg:bb4l:PartonLevel:veto", 0),
58  m_wouldVetoFsr(false),
59  m_topresscale(-1.), m_vetoDecScale(-1.), m_atopresscale(-1.),
60  m_nFSRvetoBB4l(0),
61  m_pTmin("Powheg:bb4l:pTminVeto", 0.5),
62  m_vetoTopCharge(-1),
63  m_vetoProduction("Powheg:veto", 0),
64  m_pTpythiaVeto("Powheg:bb4l:pTpythiaVeto", 0),
65  m_vetoDipoleFrame("Powheg:bb4l:FSREmission:vetoDipoleFrame", 0),
66  m_scaleResonanceVeto("Powheg:bb4l:ScaleResonance:veto", 0.)
67  {
68  std::cout << "**********************************************************" << std::endl;
69  std::cout << "* *" << std::endl;
70  std::cout << "* Applying Powheg BB4L UserHook! *" << std::endl;
71  std::cout << "* Must run on dedicated Powheg LHE input file *" << std::endl;
72  std::cout << "* (on your own responsibility) *" << std::endl;
73  std::cout << "* *" << std::endl;
74  std::cout << "**********************************************************" << std::endl;
75 
76  }
77 
78 
79  ~PowhegBB4Ltms() override { std::cout << "Number of FSR vetoed in BB4l = " << m_nFSRvetoBB4l << std::endl; }
80 
81 
82  //Initialization
83  virtual bool initAfterBeams() override {
84  // settings of the parent class
85  PowhegHooks::initAfterBeams();
86  // settings of this class
87  return true;
88  }
89 
90 
91  //--- PROCESS LEVEL HOOK ---------------------------------------------------
92  // called at the LHE level
93  inline virtual bool canVetoProcessLevel() override { return true; }
94  inline virtual bool doVetoProcessLevel(Event &e) override {
96  // count myself final particles ....
97  //int count = 0;
98  //int nFinal = settingsPtr->mode("POWHEG:nFinal");
99  //double pT1 = 0., pTsum = 0.;
100  //for (int i = e.size() - 1; i > 0; i--) {
101  // if (e[i].isFinal()) {
102  // count++;
103  // pT1 = e[i].pT();
104  // pTsum += e[i].pT();
105  // } else break;
106  //}
108  //std::cout << "INFO: counted final particles is "<< count<< "; nFinal is "<<nFinal << std::endl;
109  // }
111 
112 
113  // extract the radtype from the event comment
114  stringstream ss;
115  ss << infoPtr->getEventComments();
116  string temp;
117  ss >> temp >> radtype_.radtype;
118  assert (temp == "#rwgt");
119 
120  // find last top and the last anti-top in the record
121  int i_top = -1, i_atop = -1;
122  for (int i = 0; i < e.size(); i++) {
123  if (e[i].id() == 6) i_top = i;
124  if (e[i].id() == -6) i_atop = i;
125  }
126  if (i_top != -1)
127  m_topresscale = findresscale(i_top, e);
128  else
129  m_topresscale = 1e30;
130  if (i_top != -1)
131  m_atopresscale = findresscale(i_atop, e);
132  else
133  m_atopresscale = 1e30;
134  // initialize stuff
135  doVetoFSRInit();
136  // do not veto, ever
137  return false;
138  }
139 
140 
141  //--- PARTON LEVEL HOOK ----------------------------------------------------
142  // called after shower
143  virtual bool retryPartonLevel() override { return m_vetoPartonLevel(settingsPtr) || m_vetoAtPL(settingsPtr); }
144  inline virtual bool canVetoPartonLevel() override { return m_vetoPartonLevel(settingsPtr) || m_vetoAtPL(settingsPtr); }
145  inline virtual bool doVetoPartonLevel(const Event &e) override {
146  if(radtype_.radtype==2)
147  return false;
148  if(m_debug(settingsPtr)){
149  if (m_dryRunFSR(settingsPtr) && m_wouldVetoFsr) {
150  double scale = getdechardness(m_vetoTopCharge, e);
151  std::cout << "FSRdecScale = " << m_vetoDecScale << ", PLdecScale = " << scale << ", ratio = " << m_vetoDecScale/scale << std::endl;
152  }
153  }
154  if (m_vetoPartonLevel(settingsPtr)) {
155  double topdecscale = getdechardness(1, e);
156  double atopdecscale = getdechardness(-1, e);
157  if ((topdecscale > m_topresscale) || (atopdecscale > m_atopresscale)) {
158  return true;
159  }
160  else
161  return false;
162  }
163  if (m_vetoAtPL(settingsPtr)) {
164  if (m_dryRunFSR(settingsPtr) && m_wouldVetoFsr) return true;
165  else return false;
166  }
167  return false;
168  }
169 
170 
171 
172  //--- Internal helper functions --------------------------------------------
173  // Calculates the scale of the hardest emission from within the resonance system
174  // translated by Markus Seidel modified by Tomas Jezo
175  inline double findresscale( const int iRes, const Event& event) {
176  double scale = 0.;
177 
178  int nDau = event[iRes].daughterList().size();
179 
180  if (nDau == 0) {
181  // No resonance found, set scale to high value
182  // Pythia will shower any MC generated resonance unrestricted
183  scale = 1e30;
184  }
185  else if (nDau < 3) {
186  // No radiating resonance found
187  scale = m_pTmin(settingsPtr);
188  }
189  else if (std::abs(event[iRes].id()) == 6) {
190  // Find top daughters
191  int idw = -1, idb = -1, idg = -1;
192 
193  for (int i = 0; i < nDau; i++) {
194  int iDau = event[iRes].daughterList()[i];
195  if (std::abs(event[iDau].id()) == 24) idw = iDau;
196  if (std::abs(event[iDau].id()) == 5) idb = iDau;
197  if (std::abs(event[iDau].id()) == 21) idg = iDau;
198  }
199 
200  // Get daughter 4-vectors in resonance frame
201  Vec4 pw(event[idw].p());
202  pw.bstback(event[iRes].p());
203 
204  Vec4 pb(event[idb].p());
205  pb.bstback(event[iRes].p());
206 
207  Vec4 pg(event[idg].p());
208  pg.bstback(event[iRes].p());
209 
210  // Calculate scale
211  scale = std::sqrt(2*pg*pb*pg.e()/pb.e());
212  }
213  else {
214  scale = 1e30;
215  }
216 
217  return scale;
218  }
219 
220  // The following routine will match daughters of particle `e[iparticle]`
221  // to an expected pattern specified via the list of expected particle PDG ID's `ids`,
222  // id wildcard is specified as 0 if match is obtained, the positions and the momenta
223  // of these particles are returned in vectors `positions` and `momenta` respectively
224  // if exitOnExtraLegs==true, it will exit if the decay has more particles than expected, but not less
225  inline bool match_decay(int iparticle, const Event &e, const vector<int> &ids,
226  vector<int> &positions, vector<Vec4> &momenta, bool exitOnExtraLegs = true){
227  // compare sizes
228  if (e[iparticle].daughterList().size() != ids.size()) {
229  if (exitOnExtraLegs && e[iparticle].daughterList().size() > ids.size()) {
230  std::cout << "extra leg" << std::endl;
231  exit(-1);
232  }
233  return false;
234  }
235 
236  // compare content
237  for (unsigned int i = 0; i < e[iparticle].daughterList().size(); i++) {
238  int di = e[iparticle].daughterList()[i];
239  if (ids[i] != 0 && e[di].id() != ids[i])
240  return false;
241  }
242  // reset the positions and momenta vectors (because they may be reused)
243  positions.clear();
244  momenta.clear();
245  // construct the array of momenta
246  for (unsigned int i = 0; i < e[iparticle].daughterList().size(); i++) {
247  int di = e[iparticle].daughterList()[i];
248  positions.push_back(di);
249  momenta.push_back(e[di].p());
250  }
251  return true;
252  }
253 
254  inline double qSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2){
255  p1.bstback(pt);
256  p2.bstback(pt);
257  return std::sqrt( 2*p1*p2*p2.e()/p1.e() );
258  }
259 
260  inline double gSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2){
261  p1.bstback(pt);
262  p2.bstback(pt);
263  return std::sqrt( 2*p1*p2*p1.e()*p2.e()/(std::pow(p1.e()+p2.e(),2)) );
264  }
265  // Routines to calculate the pT (according to pTdefMode) in a FS splitting:
266  // i (radiator before) -> j (emitted after) k (radiator after)
267  // For the Pythia pT definition, a recoiler (after) must be specified.
268  // (INSPIRED BY pythia8F77_31.cc double pTpythia)
269  inline double pTpythia(const Event &e, int RadAfterBranch, int EmtAfterBranch,
270  int RecAfterBranch)
271  {
272 
273  // Convenient shorthands for later
274  Vec4 radVec = e[RadAfterBranch].p();
275  Vec4 emtVec = e[EmtAfterBranch].p();
276  Vec4 recVec = e[RecAfterBranch].p();
277  int radID = e[RadAfterBranch].id();
278 
279  // Calculate virtuality of splitting
280  Vec4 Q(radVec + emtVec);
281  double Qsq = Q.m2Calc();
282 
283 
284  // Mass term of radiator
285  double m2Rad = (std::abs(radID) >= 4 && std::abs(radID) < 7) ?
286  pow2(particleDataPtr->m0(radID)) : 0.;
287 
288  // z values for FSR
289  double z, pTnow;
290  // Construct 2 -> 3 variables
291  Vec4 sum = radVec + recVec + emtVec;
292  double m2Dip = sum.m2Calc();
293  double x1 = 2. * (sum * radVec) / m2Dip;
294  double x3 = 2. * (sum * emtVec) / m2Dip;
295  z = x1 / (x1 + x3);
296  pTnow = z * (1. - z);
297 
298  // Virtuality
299  pTnow *= (Qsq - m2Rad);
300 
301  if (pTnow < 0.) {
302 #if PYTHIA_VERSION_INTEGER >= 8310
303  loggerPtr->ERROR_MSG("Warning: pTpythia was negative");
304 #else
305  infoPtr->errorMsg("Warning: pTpythia was negative");
306 #endif
307  return -1.;
308  }
309  else
310  return(std::sqrt(pTnow));
311  }
312 
313 
314  inline double getdechardness(int topcharge, const Event &e){
315  int tid = 6*topcharge, wid = 24*topcharge, bid = 5*topcharge, gid = 21, wildcard = 0;
316  // find last top in the record
317  int i_top = -1;
318  Vec4 p_top, p_b, p_g, p_g1, p_g2;
319  for (int i = 0; i < e.size(); i++)
320  if (e[i].id() == tid) {
321  i_top = i;
322  p_top = e[i].p();
323  }
324  if (i_top == -1) return -1.0;
325 
326  // summary of cases
327  // 1.) t > W b
328  // a.) b > 3 ... error
329  // b.) b > b g ... h = std::sqrt(2*p_g*p_b*p_g.e()/p_b.e())
330  // c.) b > other ... h = -1
331  // return h
332  // 2.) t > W b g
333  // a.) b > 3 ... error
334  // b.) b > b g ... h1 = std::sqrt(2*p_g*p_b*p_g.e()/p_b.e())
335  // c.) b > other ... h1 = -1
336  // i.) g > 3 ... error
337  // ii.) g > 2 ... h2 = std::sqrt(2*p_g1*p_g2*p_g1.e()*p_g2.e()/(pow(p_g1.e(),2)+pow(p_g2.e(),2))) );
338  // iii.) g > other ... h2 = -1
339  // return max(h1,h2)
340  // 3.) else ... error
341 
342  vector<Vec4> momenta;
343  vector<int> positions;
344 
345  // 1.) t > b W
346  if ( match_decay(i_top, e, vector<int> {wid, bid}, positions, momenta, false) ) {
347  double h;
348  int i_b = positions[1];
349  // a.+b.) b > 3 or b > b g
350  if ( match_decay(i_b, e, vector<int> {bid, gid}, positions, momenta) )
351  h = qSplittingScale(e[i_top].p(), momenta[0], momenta[1]);
352  // c.) b > other
353  else
354  h = -1;
355  return h;
356  }
357  // 2.) t > b W g
358  else if ( match_decay(i_top, e, vector<int> {wid, bid, gid}, positions, momenta, false) ) {
359  double h1, h2;
360  int i_b = positions[1], i_g = positions[2];
361  // a.+b.) b > 3 or b > b g
362  if ( match_decay(i_b, e, vector<int> {bid, gid}, positions, momenta) )
363  h1 = qSplittingScale(e[i_top].p(), momenta[0], momenta[1]);
364  // c.) b > other
365  else
366  h1 = -1;
367  // i.+ii.) g > 3 or g > 2
368  if ( match_decay(i_g, e, vector<int> {wildcard, wildcard}, positions, momenta) )
369  h2 = gSplittingScale(e[i_top].p(), momenta[0], momenta[1]);
370  // c.) b > other
371  else
372  h2 = -1;
373  return max(h1, h2);
374  }
375  // 3.) else
376  else {
377  std::cout << "getdechardness" << std::endl;
378  std::cout << "top at position " << i_top << std::endl;
379  std::cout << "with " << e[i_top].daughterList().size() << " daughters " << std::endl;
380  for (unsigned int i = 0; i < e[i_top].daughterList().size(); i++) {
381  int di = e[i_top].daughterList()[i];
382  std::cout << "with daughter " << di << ": " << e[di].id() << std::endl;
383  }
384  exit(-1);
385  }
386  }
387 
388 
389 
390  //--- SCALE RESONANCE HOOK -------------------------------------------------
391  // called before each resonance decay shower
392  inline virtual bool canSetResonanceScale() override { return m_scaleResonanceVeto(settingsPtr); }
393  // if the resonance is the (anti)top set the scale to:
394  // ---> (anti)top virtuality if radtype=2
395  // ---> (a)topresscale otherwise
396  // if is not the top, set it to a big number
397  inline virtual double scaleResonance(int iRes, const Event &e) override {
398  if (e[iRes].id() == 6){
399  if(radtype_.radtype == 2)
400  return std::sqrt(e[iRes].m2Calc());
401  else
402  return m_topresscale;
403  }
404  else if (e[iRes].id() == -6){
405  if(radtype_.radtype == 2)
406  return std::sqrt(e[iRes].m2Calc());
407  else
408  return m_atopresscale;
409  }
410  else
411  return std::pow(10.0,30.);
412  }
413 
414 
415  //--- FSR EMISSION LEVEL HOOK ----------------------------------------------
416 
417  // FSR veto: this should be true if we want PowhegHooksBB4l veto in decay
418  // OR PowhegHooks veto in production. (The virtual method
419  // PowhegHooks::canVetoFSREmission has been replaced by
420  // PowhegHooksBB4L::canVetoFSREmission, so FSR veto in production
421  // must be handled here. ISR and MPI veto are instead still
422  // handled by PowhegHooks.)
423  virtual inline bool canVetoFSREmission() override { return m_vetoFSREmission(settingsPtr) || m_vetoProduction(settingsPtr); }
424  virtual inline bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance) override {
426  //VETO INSIDE THE RESONANCE //
428  if (inResonance && m_vetoFSREmission(settingsPtr)) {
429 
430  // get the participants of the splitting: the recoiler, the radiator and the emitted
431  int iRecAft = e.size() - 1;
432  int iEmt = e.size() - 2;
433  int iRadAft = e.size() - 3;
434  int iRadBef = e[iEmt].mother1();
435 
436  // find the top resonance the radiator originates from
437  int iTop = e[iRadBef].mother1();
438  int distance = 1;
439  while (std::abs(e[iTop].id()) != 6 && iTop > 0) {
440  iTop = e[iTop].mother1();
441  distance ++;
442  }
443  if (iTop == 0) {
444 #if PYTHIA_VERSION_INTEGER >= 8310
445  loggerPtr->ERROR_MSG("Warning in PowhegHooksBB4L::doVetoFSREmission: emission in resonance not from top quark, not vetoing");
446 #else
447  infoPtr->errorMsg("Warning in PowhegHooksBB4L::doVetoFSREmission: emission in resonance not from top quark, not vetoing");
448 #endif
449  return doVetoFSR(false,0,0);
450  }
451  int iTopCharge = (e[iTop].id()>0)?1:-1;
452 
453  // calculate the scale of the emission
454  double scale;
455  //using pythia pT definition ...
456  if(m_pTpythiaVeto(settingsPtr))
457  scale = pTpythia(e, iRadAft, iEmt, iRecAft);
458  else{ //.. or using POWHEG pT definition
459  Vec4 pr(e[iRadAft].p()), pe(e[iEmt].p()), pt(e[iTop].p()), prec(e[iRecAft].p()), psystem;
460  // The computation of the POWHEG pT can be done in the top rest frame or in the dipole one.
461  // pdipole = pemt +prec +prad (after the emission)
462  // For the first emission off the top resonance pdipole = pw +pb (before the emission) = ptop
463  if(m_vetoDipoleFrame(settingsPtr))
464  psystem = pr+pe+prec;
465  else
466  psystem = pt;
467 
468  // gluon splitting into two partons
469  if (e[iRadBef].id() == 21)
470  scale = gSplittingScale(psystem, pr, pe);
471  // quark emitting a gluon (or a photon)
472  else if (std::abs(e[iRadBef].id()) == 5 && ((e[iEmt].id() == 21) && ! m_vetoQED(settingsPtr)) )
473  scale = qSplittingScale(psystem, pr, pe);
474  // other stuff (which we should not veto)
475  else {
476  scale = 0;
477  }
478  }
479 
480  if (iTopCharge > 0) {
481  if (m_onlyDistance1(settingsPtr)) {
482  if ( m_debug(settingsPtr) && (distance == 1) && scale > m_topresscale && ! m_wouldVetoFsr)
483  std::cout << e[iTop].id() << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << std::endl;
484  return doVetoFSR((distance == 1) && scale > m_topresscale,scale,iTopCharge);
485  }
486  else {
487  if ( m_debug(settingsPtr) && scale > m_topresscale && ! m_wouldVetoFsr)
488  std::cout << e[iTop].id() << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << std::endl;
489  return doVetoFSR(scale > m_topresscale,scale,iTopCharge);
490  }
491  }
492  else if (iTopCharge < 0){
493  if (m_onlyDistance1(settingsPtr)){
494  if ( m_debug(settingsPtr) && (distance == 1) && scale > m_atopresscale && ! m_wouldVetoFsr)
495  std::cout << e[iTop].id() << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << std::endl;
496  return doVetoFSR((distance == 1) && scale > m_atopresscale,scale,iTopCharge);
497  }
498  else {
499  if ( m_debug(settingsPtr) && scale > m_topresscale && ! m_wouldVetoFsr)
500  std::cout << e[iTop].id() << ": " << e[iRadBef].id() << " > " << e[iRadAft].id() << " + " << e[iEmt].id() << "; " << scale << std::endl;
501  return doVetoFSR(scale > m_atopresscale,scale,iTopCharge);
502  }
503  }
504  else {
505  std::cout << "Bug in PohwgeHooksBB4l" << std::endl;
506  }
507  }
508 
510  // VETO THE PRODUCTION PROCESS //
512  else if(!inResonance && m_vetoProduction(settingsPtr)){
513  return PowhegHooks::doVetoFSREmission(sizeOld, e, iSys, inResonance);
514  }
515  return 0;
516  }
517 
518 
519 
520  inline bool doVetoFSR(bool condition, double scale, int iTopCharge) {
521  if(radtype_.radtype==2)
522  return false;
523  if (condition) {
524  if (!m_wouldVetoFsr) {
525  m_wouldVetoFsr = true;
526  m_vetoDecScale = scale;
527  m_vetoTopCharge = iTopCharge;
528  }
529  if (m_dryRunFSR(settingsPtr)) return false;
530  else {
531  m_nFSRvetoBB4l++;
532  return true;
533  }
534  }
535  else return false;
536  }
537 
538  inline void doVetoFSRInit() {
539  m_wouldVetoFsr = false;
540  m_vetoDecScale = -1;
541  m_vetoTopCharge = 0;
542  }
543 
544 
545  private:
546 
547 
549  Pythia8_UserHooks::UserSetting<bool> m_vetoFSREmission, m_dryRunFSR, m_onlyDistance1, m_vetoAtPL, m_vetoQED;
552  double m_topresscale, m_vetoDecScale, m_atopresscale;
553  unsigned long int m_nFSRvetoBB4l;
554 
557  Pythia8_UserHooks::UserSetting<int> m_vetoProduction, m_pTpythiaVeto, m_vetoDipoleFrame;
559 
560  };
561 
562  //==========================================================================
563 
564 } // end namespace Pythia8
temp
Definition: JetEventDict.h:21
plotBeamSpotCompare.x1
x1
Definition: plotBeamSpotCompare.py:216
Pythia8::PowhegBB4Ltms::m_nFSRvetoBB4l
unsigned long int m_nFSRvetoBB4l
Definition: PowhegBB4Ltms.cxx:553
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
max
#define max(a, b)
Definition: cfImp.cxx:41
PowhegControl_ttHplus_NLO.ss
ss
Definition: PowhegControl_ttHplus_NLO.py:83
Pythia8::PowhegBB4Ltms::m_vetoDecScale
double m_vetoDecScale
Definition: PowhegBB4Ltms.cxx:552
keylayer_zslicemap.pb
pb
Definition: keylayer_zslicemap.py:188
Pythia8::PowhegBB4Ltms::qSplittingScale
double qSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2)
Definition: PowhegBB4Ltms.cxx:254
Pythia8::PowhegBB4Ltms::doVetoFSREmission
virtual bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance) override
Definition: PowhegBB4Ltms.cxx:424
Pythia8::PowhegBB4Ltms::m_scaleResonanceVeto
Pythia8_UserHooks::UserSetting< double > m_scaleResonanceVeto
Definition: PowhegBB4Ltms.cxx:558
Pythia8::PowhegBB4Ltms::m_wouldVetoFsr
bool m_wouldVetoFsr
Definition: PowhegBB4Ltms.cxx:551
Pythia8::PowhegBB4Ltms::canVetoPartonLevel
virtual bool canVetoPartonLevel() override
Definition: PowhegBB4Ltms.cxx:144
Pythia8::PowhegBB4Ltms::doVetoFSRInit
void doVetoFSRInit()
Definition: PowhegBB4Ltms.cxx:538
Event
Definition: trigbs_orderedMerge.cxx:42
Pythia8::PowhegBB4Ltms::scaleResonance
virtual double scaleResonance(int iRes, const Event &e) override
Definition: PowhegBB4Ltms.cxx:397
Pythia8::PowhegBB4Ltms::pTpythia
double pTpythia(const Event &e, int RadAfterBranch, int EmtAfterBranch, int RecAfterBranch)
Definition: PowhegBB4Ltms.cxx:269
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
Pythia8::PowhegBB4Ltms::m_vetoPartonLevel
Pythia8_UserHooks::UserSetting< int > m_vetoPartonLevel
Definition: PowhegBB4Ltms.cxx:550
test_pyathena.pt
pt
Definition: test_pyathena.py:11
LArG4GenerateShowerLib.condition
condition
Definition: LArG4GenerateShowerLib.py:19
Pythia8::PowhegBB4Ltms::canVetoFSREmission
virtual bool canVetoFSREmission() override
Definition: PowhegBB4Ltms.cxx:423
Pythia8::PowhegBB4Ltms::findresscale
double findresscale(const int iRes, const Event &event)
Definition: PowhegBB4Ltms.cxx:175
Pythia8_UserHooks::UserHooksFactory::Creator
Definition: UserHooksFactory.h:54
read_hist_ntuple.h1
h1
Definition: read_hist_ntuple.py:21
yodamerge_tmp.scale
scale
Definition: yodamerge_tmp.py:138
LArG4AODNtuplePlotter.pe
pe
Definition: LArG4AODNtuplePlotter.py:116
UserHooksFactory.h
Pythia8
Author: James Monk (jmonk@cern.ch)
Definition: IPythia8Custom.h:13
HitType::wildcard
@ wildcard
Pythia8_UserHooks::UserSetting< bool >
Pythia8::PowhegBB4Ltms::m_vetoProduction
Pythia8_UserHooks::UserSetting< int > m_vetoProduction
Definition: PowhegBB4Ltms.cxx:557
python.setupRTTAlg.size
int size
Definition: setupRTTAlg.py:39
Pythia8::PowhegBB4Ltms
Definition: PowhegBB4Ltms.cxx:46
Pythia8::PowhegBB4Ltms::doVetoProcessLevel
virtual bool doVetoProcessLevel(Event &e) override
Definition: PowhegBB4Ltms.cxx:94
convertTimingResiduals.sum
sum
Definition: convertTimingResiduals.py:55
Pythia8::PowhegBB4Ltms::PowhegBB4Ltms
PowhegBB4Ltms()
Definition: PowhegBB4Ltms.cxx:50
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
Pythia8::PowhegBB4Ltms::m_debug
Pythia8_UserHooks::UserSetting< bool > m_debug
Definition: PowhegBB4Ltms.cxx:548
lumiFormat.i
int i
Definition: lumiFormat.py:92
z
#define z
h
extractSporadic.h
list h
Definition: extractSporadic.py:97
Pythia8::PowhegBB4Ltms::gSplittingScale
double gSplittingScale(Vec4 pt, Vec4 p1, Vec4 p2)
Definition: PowhegBB4Ltms.cxx:260
Pythia8::PowhegBB4Ltms::m_vetoQED
Pythia8_UserHooks::UserSetting< bool > m_vetoQED
Definition: PowhegBB4Ltms.cxx:549
Pythia8::PowhegBB4Ltms::m_pTmin
Pythia8_UserHooks::UserSetting< double > m_pTmin
Definition: PowhegBB4Ltms.cxx:555
TileDCSDataPlotter.pr
pr
Definition: TileDCSDataPlotter.py:922
UserSetting.h
Pythia8::PowhegBB4Ltms::retryPartonLevel
virtual bool retryPartonLevel() override
Definition: PowhegBB4Ltms.cxx:143
gid
Acts::GeometryIdentifier gid
Definition: ActsWriteTrackingGeometryTransforms.cxx:27
Pythia8::PowhegBB4Ltms::doVetoPartonLevel
virtual bool doVetoPartonLevel(const Event &e) override
Definition: PowhegBB4Ltms.cxx:145
calibdata.exit
exit
Definition: calibdata.py:236
Pythia8::PowhegBB4Ltms::canVetoProcessLevel
virtual bool canVetoProcessLevel() override
Definition: PowhegBB4Ltms.cxx:93
Pythia8::PowhegBB4Ltms::match_decay
bool match_decay(int iparticle, const Event &e, const vector< int > &ids, vector< int > &positions, vector< Vec4 > &momenta, bool exitOnExtraLegs=true)
Definition: PowhegBB4Ltms.cxx:225
Pythia8::PowhegBB4Ltms::getdechardness
double getdechardness(int topcharge, const Event &e)
Definition: PowhegBB4Ltms.cxx:314
ParticleGun_EoverP_Config.pg
pg
Definition: ParticleGun_EoverP_Config.py:61
python.subdetectors.mmg.ids
ids
Definition: mmg.py:8
radtype
int radtype
Definition: PowhegBB4Ltms.cxx:30
Pythia8::PowhegBB4Ltms::canSetResonanceScale
virtual bool canSetResonanceScale() override
Definition: PowhegBB4Ltms.cxx:392
Pythia8::PowhegBB4Ltms::~PowhegBB4Ltms
~PowhegBB4Ltms() override
Definition: PowhegBB4Ltms.cxx:79
UserHooksUtils.h
Pythia8::PowhegBB4Ltms::initAfterBeams
virtual bool initAfterBeams() override
Definition: PowhegBB4Ltms.cxx:83
powhegBB4LtmsCreator
Pythia8_UserHooks::UserHooksFactory::Creator< Pythia8::PowhegBB4Ltms > powhegBB4LtmsCreator("PowhegBB4Ltms")
Pythia8::PowhegBB4Ltms::m_vetoTopCharge
int m_vetoTopCharge
Definition: PowhegBB4Ltms.cxx:556
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
Pythia8::PowhegBB4Ltms::doVetoFSR
bool doVetoFSR(bool condition, double scale, int iTopCharge)
Definition: PowhegBB4Ltms.cxx:520
jobOptions.prec
prec
Definition: jobOptions.Superchic_UPC_yyMuMu.py:20
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
radtype_
struct @60 radtype_
BchCleanup.idb
idb
Definition: BchCleanup.py:152