CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EmissionVetoHook1.cc
Go to the documentation of this file.
2 
5  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
6  << "EmissionVeto: " << message << endl;
7 }
8 
9 //--------------------------------------------------------------------------
10 
11 // Routines to calculate the pT (according to pTdefMode) in a splitting:
12 // ISR: i (radiator after) -> j (emitted after) k (radiator before)
13 // FSR: i (radiator before) -> j (emitted after) k (radiator after)
14 // For the Pythia pT definition, a recoiler (after) must be specified.
15 
16 // Compute the Pythia pT separation. Based on pTLund function in History.cc
17 double EmissionVetoHook1::pTpythia(const Pythia8::Event &e, int RadAfterBranch,
18  int EmtAfterBranch, int RecAfterBranch, bool FSR) {
19 
20  // Convenient shorthands for later
21  Pythia8::Vec4 radVec = e[RadAfterBranch].p();
22  Pythia8::Vec4 emtVec = e[EmtAfterBranch].p();
23  Pythia8::Vec4 recVec = e[RecAfterBranch].p();
24  int radID = e[RadAfterBranch].id();
25 
26  // Calculate virtuality of splitting
27  double sign = (FSR) ? 1. : -1.;
28  Pythia8::Vec4 Q(radVec + sign * emtVec);
29  double Qsq = sign * Q.m2Calc();
30 
31  // Mass term of radiator
32  double m2Rad = (abs(radID) >= 4 && abs(radID) < 7) ?
33  Pythia8::pow2(particleDataPtr->m0(radID)) : 0.;
34 
35  // z values for FSR and ISR
36  double z, pTnow;
37  if (FSR) {
38  // Construct 2 -> 3 variables
39  Pythia8::Vec4 sum = radVec + recVec + emtVec;
40  double m2Dip = sum.m2Calc();
41  double x1 = 2. * (sum * radVec) / m2Dip;
42  double x3 = 2. * (sum * emtVec) / m2Dip;
43  z = x1 / (x1 + x3);
44  pTnow = z * (1. - z);
45 
46  } else {
47  // Construct dipoles before/after splitting
48  Pythia8::Vec4 qBR(radVec - emtVec + recVec);
49  Pythia8::Vec4 qAR(radVec + recVec);
50  z = qBR.m2Calc() / qAR.m2Calc();
51  pTnow = (1. - z);
52  }
53 
54  // Virtuality with correct sign
55  pTnow *= (Qsq - sign * m2Rad);
56 
57  // Can get negative pT for massive splittings
58  if (pTnow < 0.) {
59  cout << "Warning: pTpythia was negative" << endl;
60  return -1.;
61  }
62 
63 #ifdef DBGOUTPUT
64  cout << "pTpythia: rad = " << RadAfterBranch << ", emt = "
65  << EmtAfterBranch << ", rec = " << RecAfterBranch
66  << ", pTnow = " << sqrt(pTnow) << endl;
67 #endif
68 
69  // Return pT
70  return sqrt(pTnow);
71 }
72 
73 // Compute the POWHEG pT separation between i and j
74 double EmissionVetoHook1::pTpowheg(const Pythia8::Event &e, int i, int j, bool FSR) {
75 
76  // pT value for FSR and ISR
77  double pTnow = 0.;
78  if (FSR) {
79  // POWHEG d_ij (in CM frame). Note that the incoming beams have not
80  // been updated in the parton systems pointer yet (i.e. prior to any
81  // potential recoil).
82  int iInA = partonSystemsPtr->getInA(0);
83  int iInB = partonSystemsPtr->getInB(0);
84  double betaZ = - ( e[iInA].pz() + e[iInB].pz() ) /
85  ( e[iInA].e() + e[iInB].e() );
86  Pythia8::Vec4 iVecBst(e[i].p()), jVecBst(e[j].p());
87  iVecBst.bst(0., 0., betaZ);
88  jVecBst.bst(0., 0., betaZ);
89  pTnow = sqrt( (iVecBst + jVecBst).m2Calc() *
90  iVecBst.e() * jVecBst.e() /
91  Pythia8::pow2(iVecBst.e() + jVecBst.e()) );
92 
93  } else {
94  // POWHEG pT_ISR is just kinematic pT
95  pTnow = e[j].pT();
96  }
97 
98  // Check result
99  if (pTnow < 0.) {
100  cout << "Warning: pTpowheg was negative" << endl;
101  return -1.;
102  }
103 
104 #ifdef DBGOUTPUT
105  cout << "pTpowheg: i = " << i << ", j = " << j
106  << ", pTnow = " << pTnow << endl;
107 #endif
108 
109  return pTnow;
110 }
111 
112 // Calculate pT for a splitting based on pTdefMode.
113 // If j is -1, all final-state partons are tried.
114 // If i, k, r and xSR are -1, then all incoming and outgoing
115 // partons are tried.
116 // xSR set to 0 means ISR, while xSR set to 1 means FSR
117 double EmissionVetoHook1::pTcalc(const Pythia8::Event &e, int i, int j, int k, int r, int xSRin) {
118 
119  // Loop over ISR and FSR if necessary
120  double pTemt = -1., pTnow;
121  int xSR1 = (xSRin == -1) ? 0 : xSRin;
122  int xSR2 = (xSRin == -1) ? 2 : xSRin + 1;
123  for (int xSR = xSR1; xSR < xSR2; xSR++) {
124  // FSR flag
125  bool FSR = (xSR == 0) ? false : true;
126 
127  // If all necessary arguments have been given, then directly calculate.
128  // POWHEG ISR and FSR, need i and j.
129  if ((pTdefMode == 0 || pTdefMode == 1) && i > 0 && j > 0) {
130  pTemt = pTpowheg(e, i, j, (pTdefMode == 0) ? false : FSR);
131 
132  // Pythia ISR, need i, j and r.
133  } else if (!FSR && pTdefMode == 2 && i > 0 && j > 0 && r > 0) {
134  pTemt = pTpythia(e, i, j, r, FSR);
135 
136  // Pythia FSR, need k, j and r.
137  } else if (FSR && pTdefMode == 2 && j > 0 && k > 0 && r > 0) {
138  pTemt = pTpythia(e, k, j, r, FSR);
139 
140  // Otherwise need to try all possible combintations.
141  } else {
142  // Start by finding incoming legs to the hard system after
143  // branching (radiator after branching, i for ISR).
144  // Use partonSystemsPtr to find incoming just prior to the
145  // branching and track mothers.
146  int iInA = partonSystemsPtr->getInA(0);
147  int iInB = partonSystemsPtr->getInB(0);
148  while (e[iInA].mother1() != 1) { iInA = e[iInA].mother1(); }
149  while (e[iInB].mother1() != 2) { iInB = e[iInB].mother1(); }
150 
151  // If we do not have j, then try all final-state partons
152  int jNow = (j > 0) ? j : 0;
153  int jMax = (j > 0) ? j + 1 : e.size();
154  for (; jNow < jMax; jNow++) {
155 
156  // Final-state and coloured jNow or photon only
157  if (!e[jNow].isFinal()) continue;
158  if (e[jNow].colType() == 0 && e[jNow].id() != 22) continue;
159 
160  // POWHEG
161  if (pTdefMode == 0 || pTdefMode == 1) {
162 
163  // ISR - only done once as just kinematical pT
164  if (!FSR) {
165  pTnow = pTpowheg(e, iInA, jNow, (pTdefMode == 0) ? false : FSR);
166  if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
167 
168  // FSR - try all outgoing partons from system before branching
169  // as i. Note that for the hard system, there is no
170  // "before branching" information.
171  } else {
172 
173  int outSize = partonSystemsPtr->sizeOut(0);
174  for (int iMem = 0; iMem < outSize; iMem++) {
175  int iNow = partonSystemsPtr->getOut(0, iMem);
176 
177  // Coloured only, i != jNow and no carbon copies
178  if (iNow == jNow || e[iNow].colType() == 0) continue;
179  if (jNow == e[iNow].daughter1()
180  && jNow == e[iNow].daughter2()) continue;
181 
182  pTnow = pTpowheg(e, iNow, jNow, (pTdefMode == 0)
183  ? false : FSR);
184  if (pTnow > 0.) pTemt = (pTemt < 0)
185  ? pTnow : min(pTemt, pTnow);
186  } // for (iMem)
187 
188  } // if (!FSR)
189 
190  // Pythia
191  } else if (pTdefMode == 2) {
192 
193  // ISR - other incoming as recoiler
194  if (!FSR) {
195  pTnow = pTpythia(e, iInA, jNow, iInB, FSR);
196  if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
197  pTnow = pTpythia(e, iInB, jNow, iInA, FSR);
198  if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
199 
200  // FSR - try all final-state coloured partons as radiator
201  // after emission (k).
202  } else {
203  for (int kNow = 0; kNow < e.size(); kNow++) {
204  if (kNow == jNow || !e[kNow].isFinal() ||
205  e[kNow].colType() == 0) continue;
206 
207  // For this kNow, need to have a recoiler.
208  // Try two incoming.
209  pTnow = pTpythia(e, kNow, jNow, iInA, FSR);
210  if (pTnow > 0.) pTemt = (pTemt < 0)
211  ? pTnow : min(pTemt, pTnow);
212  pTnow = pTpythia(e, kNow, jNow, iInB, FSR);
213  if (pTnow > 0.) pTemt = (pTemt < 0)
214  ? pTnow : min(pTemt, pTnow);
215 
216  // Try all other outgoing.
217  for (int rNow = 0; rNow < e.size(); rNow++) {
218  if (rNow == kNow || rNow == jNow ||
219  !e[rNow].isFinal() || e[rNow].colType() == 0) continue;
220  pTnow = pTpythia(e, kNow, jNow, rNow, FSR);
221  if (pTnow > 0.) pTemt = (pTemt < 0)
222  ? pTnow : min(pTemt, pTnow);
223  } // for (rNow)
224 
225  } // for (kNow)
226  } // if (!FSR)
227  } // if (pTdefMode)
228  } // for (j)
229  }
230  } // for (xSR)
231 
232 #ifdef DBGOUTPUT
233  cout << "pTcalc: i = " << i << ", j = " << j << ", k = " << k
234  << ", r = " << r << ", xSR = " << xSRin
235  << ", pTemt = " << pTemt << endl;
236 #endif
237 
238  return pTemt;
239 }
240 
241 //--------------------------------------------------------------------------
242 
243 // Extraction of pThard based on the incoming event.
244 // Assume that all the final-state particles are in a continuous block
245 // at the end of the event and the final entry is the POWHEG emission.
246 // If there is no POWHEG emission, then pThard is set to QRen.
247 
248 bool EmissionVetoHook1::doVetoMPIStep(int nMPI, const Pythia8::Event &e) {
249  // Extra check on nMPI
250  if (nMPI > 1) return false;
251 
252  // Find if there is a POWHEG emission. Go backwards through the
253  // event record until there is a non-final particle. Also sum pT and
254  // find pT_1 for possible MPI vetoing
255  int count = 0, inonfinal = 0;
256  double pT1 = 0., pTsum = 0.;
257  for (int i = e.size() - 1; i > 0; i--) {
258  inonfinal = i;
259  if (e[i].isFinal()) {
260  count++;
261  pT1 = e[i].pT();
262  pTsum += e[i].pT();
263  } else break;
264  }
265 
266  nFinal = nFinalExt;
267 
268  if (nFinal < 0) { // nFinal is not specified from external, try to find out
269  int first = -1, myid;
270  int last = -1;
271  for(int ip = 2; ip < e.size(); ip++) {
272  myid = e[ip].id();
273  if(abs(myid) < 6 || abs(myid) == 21) continue;
274  first = ip;
275  break;
276  }
277  if(first < 0) fatalEmissionVeto(string("signal particles not found"));
278  for(int ip = first; ip < e.size(); ip++) {
279  myid = e[ip].id();
280  if(abs(myid) < 6 || abs(myid) == 21) continue;
281  last = ip;
282  }
283  nFinal = last - inonfinal;
284  }
285 
286  // Extra check that we have the correct final state
287  if (count != nFinal && count != nFinal + 1)
288  fatalEmissionVeto(string("Wrong number of final state particles in event"));
289 
290  // Flag if POWHEG radiation present and index
291  bool isEmt = (count == nFinal) ? false : true;
292  int iEmt = (isEmt) ? e.size() - 1 : -1;
293 
294  // If there is no radiation or if pThardMode is 0 then set pThard to QRen.
295  if (!isEmt || pThardMode == 0) {
296  pThard = infoPtr->QRen();
297 
298  // If pThardMode is 1 then the pT of the POWHEG emission is checked against
299  // all other incoming and outgoing partons, with the minimal value taken
300  } else if (pThardMode == 1) {
301  pThard = pTcalc(e, -1, iEmt, -1, -1, -1);
302 
303  // If pThardMode is 2, then the pT of all final-state partons is checked
304  // against all other incoming and outgoing partons, with the minimal value
305  // taken
306  } else if (pThardMode == 2) {
307  pThard = pTcalc(e, -1, -1, -1, -1, -1);
308 
309  }
310 
311  // Find MPI veto pT if necessary
312  if (MPIvetoOn) {
313  pTMPI = (isEmt) ? pTsum / 2. : pT1;
314  }
315 
316  if(Verbosity)
317  cout << "doVetoMPIStep: QFac = " << infoPtr->QFac()
318  << ", QRen = " << infoPtr->QRen()
319  << ", pThard = " << pThard << endl << endl;
320 
321  // Initialise other variables
322  accepted = false;
323  nAcceptSeq = 0;
324 
325  // Do not veto the event
326  return false;
327 }
328 
329 //--------------------------------------------------------------------------
330 
331 // ISR veto
332 
333 bool EmissionVetoHook1::doVetoISREmission(int, const Pythia8::Event &e, int iSys) {
334  // Must be radiation from the hard system
335  if (iSys != 0) return false;
336 
337  // If we already have accepted 'vetoCount' emissions in a row, do nothing
338  if (vetoOn && nAcceptSeq >= vetoCount) return false;
339 
340  // Pythia radiator after, emitted and recoiler after.
341  int iRadAft = -1, iEmt = -1, iRecAft = -1;
342  for (int i = e.size() - 1; i > 0; i--) {
343  if (iRadAft == -1 && e[i].status() == -41) iRadAft = i;
344  else if (iEmt == -1 && e[i].status() == 43) iEmt = i;
345  else if (iRecAft == -1 && e[i].status() == -42) iRecAft = i;
346  if (iRadAft != -1 && iEmt != -1 && iRecAft != -1) break;
347  }
348  if (iRadAft == -1 || iEmt == -1 || iRecAft == -1) {
349  e.list();
350  fatalEmissionVeto(string("Couldn't find Pythia ISR emission"));
351  }
352 
353  // pTemtMode == 0: pT of emitted w.r.t. radiator
354  // pTemtMode == 1: min(pT of emitted w.r.t. all incoming/outgoing)
355  // pTemtMode == 2: min(pT of all outgoing w.r.t. all incoming/outgoing)
356  int xSR = (pTemtMode == 0) ? 0 : -1;
357  int i = (pTemtMode == 0) ? iRadAft : -1;
358  int j = (pTemtMode != 2) ? iEmt : -1;
359  int k = -1;
360  int r = (pTemtMode == 0) ? iRecAft : -1;
361  double pTemt = pTcalc(e, i, j, k, r, xSR);
362 
363 #ifdef DBGOUTPUT
364  cout << "doVetoISREmission: pTemt = " << pTemt << endl << endl;
365 #endif
366 
367  // Veto if pTemt > pThard
368  if (pTemt > pThard) {
369  nAcceptSeq = 0;
370  nISRveto++;
371  return true;
372  }
373 
374  // Else mark that an emission has been accepted and continue
375  nAcceptSeq++;
376  accepted = true;
377  return false;
378 }
379 
380 //--------------------------------------------------------------------------
381 
382 // FSR veto
383 
384 bool EmissionVetoHook1::doVetoFSREmission(int, const Pythia8::Event &e, int iSys, bool) {
385  // Must be radiation from the hard system
386  if (iSys != 0) return false;
387 
388  // If we already have accepted 'vetoCount' emissions in a row, do nothing
389  if (vetoOn && nAcceptSeq >= vetoCount) return false;
390 
391  // Pythia radiator (before and after), emitted and recoiler (after)
392  int iRecAft = e.size() - 1;
393  int iEmt = e.size() - 2;
394  int iRadAft = e.size() - 3;
395  int iRadBef = e[iEmt].mother1();
396  if ( (e[iRecAft].status() != 52 && e[iRecAft].status() != -53) ||
397  e[iEmt].status() != 51 || e[iRadAft].status() != 51) {
398  e.list();
399  fatalEmissionVeto(string("Couldn't find Pythia FSR emission"));
400  }
401 
402  // Behaviour based on pTemtMode:
403  // 0 - pT of emitted w.r.t. radiator before
404  // 1 - min(pT of emitted w.r.t. all incoming/outgoing)
405  // 2 - min(pT of all outgoing w.r.t. all incoming/outgoing)
406  int xSR = (pTemtMode == 0) ? 1 : -1;
407  int i = (pTemtMode == 0) ? iRadBef : -1;
408  int k = (pTemtMode == 0) ? iRadAft : -1;
409  int r = (pTemtMode == 0) ? iRecAft : -1;
410 
411  // When pTemtMode is 0 or 1, iEmt has been selected
412  double pTemt = -1.;
413  if (pTemtMode == 0 || pTemtMode == 1) {
414  // Which parton is emitted, based on emittedMode:
415  // 0 - Pythia definition of emitted
416  // 1 - Pythia definition of radiated after emission
417  // 2 - Random selection of emitted or radiated after emission
418  // 3 - Try both emitted and radiated after emission
419  int j = iRadAft;
420  if (emittedMode == 0 || (emittedMode == 2 && rndmPtr->flat() < 0.5)) j++;
421 
422  for (int jLoop = 0; jLoop < 2; jLoop++) {
423  if (jLoop == 0) pTemt = pTcalc(e, i, j, k, r, xSR);
424  else if (jLoop == 1) pTemt = min(pTemt, pTcalc(e, i, j, k, r, xSR));
425 
426  // For emittedMode == 3, have tried iRadAft, now try iEmt
427  if (emittedMode != 3) break;
428  if (k != -1) swap(j, k); else j = iEmt;
429  }
430 
431  // If pTemtMode is 2, then try all final-state partons as emitted
432  } else if (pTemtMode == 2) {
433  pTemt = pTcalc(e, i, -1, k, r, xSR);
434 
435  }
436 
437 #ifdef DBGOUTPUT
438  cout << "doVetoFSREmission: pTemt = " << pTemt << endl << endl;
439 #endif
440 
441  // Veto if pTemt > pThard
442  if (pTemt > pThard) {
443  nAcceptSeq = 0;
444  nFSRveto++;
445  return true;
446  }
447 
448  // Else mark that an emission has been accepted and continue
449  nAcceptSeq++;
450  accepted = true;
451  return false;
452 }
453 
454 //--------------------------------------------------------------------------
455 
456 // MPI veto
457 
458 bool EmissionVetoHook1::doVetoMPIEmission(int, const Pythia8::Event &e) {
459  if (MPIvetoOn) {
460  if (e[e.size() - 1].pT() > pTMPI) {
461 #ifdef DBGOUTPUT
462  cout << "doVetoMPIEmission: pTnow = " << e[e.size() - 1].pT()
463  << ", pTMPI = " << pTMPI << endl << endl;
464 #endif
465  return true;
466  }
467  }
468  return false;
469 }
470 
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
int i
Definition: DBlmapReader.cc:9
bool doVetoMPIEmission(int, const Pythia8::Event &e)
bool doVetoISREmission(int, const Pythia8::Event &e, int iSys)
ExtVec< T, 4 > Vec4
Definition: ExtVec.h:49
double sign(double x)
void fatalEmissionVeto(string message)
float float float z
unsigned long int nISRveto
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
double pTpowheg(const Pythia8::Event &e, int i, int j, bool FSR)
unsigned long int nFSRveto
double pTcalc(const Pythia8::Event &e, int i, int j, int k, int r, int xSRin)
bool doVetoMPIStep(int nMPI, const Pythia8::Event &e)
double pTpythia(const Pythia8::Event &e, int RadAfterBranch, int EmtAfterBranch, int RecAfterBranch, bool FSR)
bool doVetoFSREmission(int, const Pythia8::Event &e, int iSys, bool)
tuple cout
Definition: gather_cfg.py:121
tuple status
Definition: ntuplemaker.py:245