gambit is hosted by Hepforge, IPPP Durham
GAMBIT  v1.5.0-2191-ga4742ac
a Global And Modular Bsm Inference Tool
Analysis_ATLAS_8TeV_2LEPEW_20invfb.cpp
Go to the documentation of this file.
1 #include <vector>
2 #include <cmath>
3 #include <memory>
4 #include <iomanip>
5 
9 
10 /* The ATLAS 2 lepton EW analysis (20fb^-1)
11 
12  based on: arXiv: 1403.5294
13 
14  Code by Martin White
15 
16  Known features: Signal leptons in the paper have certain isolation plus ID cuts (these are ignored here by default)
17 
18  a) The trigger efficiencies are important! We need to multiply by these to get the correct yields (after which the cutflows agree rather well).
19 
20 */
21 
22 namespace Gambit {
23  namespace ColliderBit {
24 
25  using namespace std;
26 
27  // double phi_mpi_pi_2lep(double x){
28  // while (x >= pi) x -= 2*M_PI;
29  // while (x < -pi) x += 2*M_PI;
30  // }
31 
32 
33  bool sortByPT_2lep(const HEPUtils::Particle* lep1, const HEPUtils::Particle* lep2) { return (lep1->pT() > lep2->pT()); }
34 
35  class Analysis_ATLAS_8TeV_2LEPEW_20invfb : public Analysis {
36  private:
37 
38  // Numbers passing cuts (doubles because we will use the trigger eff)
39  double _num_MT2_90_SF;
40  double _num_MT2_90_DF;
41  double _num_MT2_120_SF;
42  double _num_MT2_120_DF;
43  double _num_MT2_150_SF;
44  double _num_MT2_150_DF;
45  double _num_WWa_SF;
46  double _num_WWa_DF;
47  double _num_WWb_SF;
48  double _num_WWb_DF;
49  double _num_WWc_SF;
50  double _num_WWc_DF;
51  double _num_Zjets;
52 
53  vector<double> cutFlowVector;
54  vector<double> cutFlowIncrements;
55  vector<string> cutFlowVector_str;
56  const static int NCUTS=90;
57 
58  public:
59 
60  // Required detector sim
61  static constexpr const char* detector = "ATLAS";
62 
63  Analysis_ATLAS_8TeV_2LEPEW_20invfb() {
64 
65  set_analysis_name("ATLAS_8TeV_2LEPEW_20invfb");
66  set_luminosity(20.3);
67 
68  _num_MT2_90_SF=0;
69  _num_MT2_90_DF=0;
70  _num_MT2_120_SF=0;
71  _num_MT2_120_DF=0;
72  _num_MT2_150_SF=0;
73  _num_MT2_150_DF=0;
74  _num_WWa_SF=0;
75  _num_WWa_DF=0;
76  _num_WWb_SF=0;
77  _num_WWb_DF=0;
78  _num_WWc_SF=0;
79  _num_WWc_DF=0;
80  _num_Zjets=0;
81 
82  for(int i=0;i<NCUTS;i++){
83  cutFlowVector.push_back(0);
84  cutFlowVector_str.push_back("");
85  cutFlowIncrements.push_back(0.);
86  }
87 
88  }
89 
90  void EleEleOverlapRemoval(vector<const HEPUtils::Particle*> &vec1, vector<const HEPUtils::Particle*> &vec2, double DeltaRMax) {
91  //Routine to do electron-electron overlap check
92  //Discard lowest energy electron if two are found overlapping
93  vector<const HEPUtils::Particle*> Survivors;
94 
95  for(unsigned int it1 = 0; it1 < vec1.size(); it1++) {
96  bool overlap = false;
97  HEPUtils::P4 lep1mom=vec1.at(it1)->mom();
98  for(unsigned int it2 = 0; it2 < vec2.size(); it2++) {
99  if(it1==it2)continue;
100  HEPUtils::P4 lep2mom=vec2.at(it2)->mom();
101  double dR;
102 
103  dR=lep1mom.deltaR_eta(lep2mom);
104 
105  if(fabs(dR) <= DeltaRMax && lep1mom.E()<lep2mom.E()) overlap=true;
106  }
107  if(overlap) continue;
108  Survivors.push_back(vec1.at(it1));
109  }
110  vec1=Survivors;
111 
112  return;
113  }
114 
115  void LepLepOverlapRemoval(vector<const HEPUtils::Particle*> &vec1, vector<const HEPUtils::Particle*> &vec2, double DeltaRMax) {
116  //Routine to do lepton-lepton overlap check
117  //Discard first lepton if overlap is found
118  vector<const HEPUtils::Particle*> Survivors;
119 
120  for(unsigned int it1 = 0; it1 < vec1.size(); it1++) {
121  bool overlap = false;
122  HEPUtils::P4 lep1mom=vec1.at(it1)->mom();
123  for(unsigned int it2 = 0; it2 < vec2.size(); it2++) {
124  if(it1==it2)continue;
125  HEPUtils::P4 lep2mom=vec2.at(it2)->mom();
126  double dR;
127 
128  dR=lep1mom.deltaR_eta(lep2mom);
129 
130  if(fabs(dR) <= DeltaRMax)overlap=true;
131  }
132  if(overlap) continue;
133  Survivors.push_back(vec1.at(it1));
134  }
135  vec1=Survivors;
136 
137  return;
138  }
139 
140  void JetLeptonOverlapRemoval(vector<const HEPUtils::Jet*> &jetvec, vector<const HEPUtils::Particle*> &lepvec, double DeltaRMax) {
141  //Routine to do jet-lepton check
142  //Discards jets if they are within DeltaRMax of a lepton
143 
144  vector<const HEPUtils::Jet*> Survivors;
145 
146  for(unsigned int itjet = 0; itjet < jetvec.size(); itjet++) {
147  bool overlap = false;
148  HEPUtils::P4 jetmom=jetvec.at(itjet)->mom();
149  for(unsigned int itlep = 0; itlep < lepvec.size(); itlep++) {
150  HEPUtils::P4 lepmom=lepvec.at(itlep)->mom();
151  double dR=jetmom.deltaR_eta(lepmom);
152 
153  if(fabs(dR) <= DeltaRMax) overlap=true;
154  }
155  if(overlap) continue;
156  Survivors.push_back(jetvec.at(itjet));
157  }
158  jetvec=Survivors;
159 
160  return;
161  }
162 
163  void LeptonJetOverlapRemoval(vector<const HEPUtils::Particle*> &lepvec, vector<const HEPUtils::Jet*> &jetvec, double DeltaRMax) {
164  //Routine to do lepton-jet check
165  //Discards leptons if they are within DeltaRMax of a jet
166 
167  vector<const HEPUtils::Particle*> Survivors;
168 
169  for(unsigned int itlep = 0; itlep < lepvec.size(); itlep++) {
170  bool overlap = false;
171  HEPUtils::P4 lepmom=lepvec.at(itlep)->mom();
172  for(unsigned int itjet= 0; itjet < jetvec.size(); itjet++) {
173  HEPUtils::P4 jetmom=jetvec.at(itjet)->mom();
174  double dR;
175 
176  dR=jetmom.deltaR_eta(lepmom);
177 
178  if(fabs(dR) <= DeltaRMax) overlap=true;
179  }
180  if(overlap) continue;
181  Survivors.push_back(lepvec.at(itlep));
182  }
183  lepvec=Survivors;
184 
185  return;
186  }
187 
188  void RemoveLeptonsMllLt12(vector<const HEPUtils::Particle*> &lepvec){
189 
190  ssize_t removeLep1=-1;
191  ssize_t removeLep2=-1;
192  vector<const HEPUtils::Particle*> Survivors;
193 
194  //Function removes SF lepton pairs with m_ll < 12 GeV
195  for(unsigned int itlep1 = 0; itlep1 < lepvec.size(); itlep1++) {
196  HEPUtils::P4 lepmom1=lepvec.at(itlep1)->mom();
197  for(unsigned int itlep2= itlep1; itlep2 < lepvec.size(); itlep2++) {
198  if(itlep1!=itlep2){
199  HEPUtils::P4 lepmom2=lepvec.at(itlep2)->mom();
200  double mass=(lepmom1+lepmom2).m();
201  if(mass<12.){
202  removeLep1=itlep1;
203  removeLep2=itlep2;
204  }
205  }
206  }
207  }
208  for(unsigned int itlep = 0; itlep < lepvec.size(); itlep++) {
209  if(itlep!=removeLep1 && itlep!=removeLep2) Survivors.push_back(lepvec.at(itlep));
210  }
211 
212  lepvec=Survivors;
213  }
214 
215  void run(const HEPUtils::Event* event) {
216 
217  // Missing energy
218  HEPUtils::P4 ptot = event->missingmom();
219  double met = event->met();
220 
221  // Now define vector of baseline electrons
222  vector<const HEPUtils::Particle*> signalElectrons;
223  for (const HEPUtils::Particle* electron : event->electrons()) {
224  if (electron->pT() > 10. && fabs(electron->eta()) < 2.47) signalElectrons.push_back(electron);
225  }
226 
227  // Apply electron efficiency
228  ATLAS::applyElectronEff(signalElectrons);
229 
230  // Now define vector of baseline muons
231  vector<const HEPUtils::Particle*> signalMuons;
232  for (const HEPUtils::Particle* muon : event->muons()) {
233  if (muon->pT() > 10. && fabs(muon->eta()) < 2.4) signalMuons.push_back(muon);
234  }
235 
236  // Apply muon efficiency
237  ATLAS::applyMuonEff(signalMuons);
238 
239  vector<const HEPUtils::Jet*> signalJets;
240  for (const HEPUtils::Jet* jet : event->jets()) {
241  if (jet->pT() > 20. && fabs(jet->eta()) < 4.5) signalJets.push_back(jet);
242  //if(jet->btag() && fabs(jet->eta()) < 2.5 && jet->pT() > 20.) bJets.push_back(jet);
243  }
244 
245  vector<const HEPUtils::Particle*> signalTaus;
246  for (const HEPUtils::Particle* tau : event->taus()) {
247  if (tau->pT() > 20. && tau->abseta() < 2.5) signalTaus.push_back(tau);
248  }
249  ATLAS::applyTauEfficiencyR1(signalTaus);
250 
251  // Overlap removal
252 
253  //Note that ATLAS use |eta|<10 for removing jets close to electrons
254  //Then 2.8 is used for the rest of the overlap process
255  //Then the signal cut is applied for signal jets
256 
257  EleEleOverlapRemoval(signalElectrons,signalElectrons,0.05);
258  JetLeptonOverlapRemoval(signalJets,signalElectrons,0.2);
259  LepLepOverlapRemoval(signalTaus,signalElectrons,0.2);
260  LepLepOverlapRemoval(signalTaus,signalMuons,0.2);
261  LeptonJetOverlapRemoval(signalElectrons,signalJets,0.4);
262  LeptonJetOverlapRemoval(signalMuons,signalJets,0.4);
263  //Note have not bothered with close-by electron and muon pairs (bremsstrahlung probably not significant in signal MC)
264 
265  RemoveLeptonsMllLt12(signalElectrons);
266  RemoveLeptonsMllLt12(signalMuons);
267  JetLeptonOverlapRemoval(signalJets,signalTaus,0.2);
268 
269  ATLAS::applyTightIDElectronSelection(signalElectrons);
270 
271  int numElectrons=signalElectrons.size();
272  int numMuons=signalMuons.size();
273  int numTaus=signalTaus.size();
274 
275  //Search for at least one SFOS pair
276  //m_SFOS must be > 12 GeV
277 
278 
279  //Classify jets into various categories
280  vector<const HEPUtils::Jet*> centralBJets;
281  vector<const HEPUtils::Jet*> centralNonBJets;
282  vector<const HEPUtils::Jet*> forwardJets;
283 
284  const std::vector<double> a = {0,10.};
285  const std::vector<double> b = {0,10000.};
286  const std::vector<double> c = {0.8};
287  HEPUtils::BinnedFn2D<double> _eff2d(a,b,c);
288 
289  for (const HEPUtils::Jet* jet : signalJets) {
290  bool hasTag=has_tag(_eff2d, jet->abseta(), jet->pT());
291  if(fabs(jet->eta()) < 2.4){
292  if(jet->btag() && hasTag){
293  centralBJets.push_back(jet);
294  }
295  else {
296  centralNonBJets.push_back(jet);
297  }
298  }
299  if(fabs(jet->eta()) > 2.4 && jet->pT()>30.)forwardJets.push_back(jet);
300  }
301 
302  //Common cuts for all signal regions
303 
304  bool leptonPTCut=false;
305  vector<const HEPUtils::Particle*> signalLeptons;
306  for (const HEPUtils::Particle* ele : signalElectrons) {
307  signalLeptons.push_back(ele);
308  }
309 
310  for (const HEPUtils::Particle* muo : signalMuons) {
311  signalLeptons.push_back(muo);
312  }
313 
314  std::sort(signalLeptons.begin(), signalLeptons.end(), sortByPT_2lep);
315 
316 
317  if(signalLeptons.size()==2 && signalLeptons[0]->pT()>35. && signalLeptons[1]->pT()>20.)leptonPTCut=true;
318 
319  bool mllCut=false;
320  if(signalLeptons.size()==2 && (signalLeptons[0]->mom()+signalLeptons[1]->mom()).m() > 20.)mllCut=true;
321 
322  bool isOS=false;
323  if(signalLeptons.size()==2 && (signalLeptons[0]->pid()*signalLeptons[1]->pid()<0))isOS=true;
324 
325  bool passZVeto=true;
326  double mLepLep=0.;
327  if(signalLeptons.size()==2)mLepLep=(signalLeptons[0]->mom()+signalLeptons[1]->mom()).m();
328  if(mLepLep>81.2 && mLepLep<101.2)passZVeto=false;
329 
330  bool cut_SRMT290=false;
331  bool cut_SRMT2120=false;
332  bool cut_SRMT2150=false;
333 
334  bool tauVeto=false;
335  if(numTaus==0)tauVeto=true;
336 
337  int numCentralNonBJets=centralNonBJets.size();
338  int numCentralBJets=centralBJets.size();
339  int numForwardJets=forwardJets.size();
340 
341  //Now do the MT2 signal regions
342 
343  if(leptonPTCut && mllCut && isOS && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0){
344 
345  //Calculate MT2
346  double pa[3] = { 0, signalLeptons[0]->mom().px(), signalLeptons[0]->mom().py() };
347  double pb[3] = { 0, signalLeptons[1]->mom().px(), signalLeptons[1]->mom().py() };
348  double pmiss[3] = { 0, ptot.px(), ptot.py() };
349  double mn = 0.;
350 
351  mt2_bisect::mt2 mt2_calc;
352 
353  mt2_calc.set_momenta(pa,pb,pmiss);
354  mt2_calc.set_mn(mn);
355  double mt2 = mt2_calc.get_mt2();
356 
357  double mll=(signalLeptons[0]->mom() + signalLeptons[1]->mom()).m();
358 
359 
360  if(mt2>90.)cut_SRMT290=true;
361  if(mt2>120.)cut_SRMT2120=true;
362  if(mt2>150.)cut_SRMT2150=true;
363 
364  //Signal region increments use the trigger efficiencies for ee, emu and mumu triggers
365  if(mt2 > 90. && (numElectrons==1 && numMuons==1)) _num_MT2_90_DF += event->weight() * 0.89;
366  if(passZVeto && mt2 > 90. && (numElectrons==2 && fabs(mll-91.)>10)) _num_MT2_90_SF += event->weight() * 0.97;
367  if(passZVeto && mt2 > 90. && (numMuons==2 && fabs(mll-91.)>10)) _num_MT2_90_SF += event->weight() * 0.75;
368 
369  if(mt2 > 120. && (numElectrons==1 && numMuons==1)) _num_MT2_120_DF += event->weight() * 0.89;
370  if(passZVeto && mt2 > 120. && (numElectrons==2 && fabs(mll-91.)>10)) _num_MT2_120_SF += event->weight() * 0.97;
371  if(passZVeto && mt2 > 120. && (numMuons==2 && fabs(mll-91.)>10)) _num_MT2_120_SF += event->weight() * 0.75;
372 
373  if(mt2 > 150. && (numElectrons==1 && numMuons==1)) _num_MT2_150_DF += event->weight() * 0.89;
374  if(passZVeto && mt2 > 150. && (numElectrons==2 && fabs(mll-91.)>10)) _num_MT2_150_SF += event->weight() * 0.97;
375  if(passZVeto && mt2 > 150. && (numMuons==2 && fabs(mll-91.)>10)) _num_MT2_150_SF += event->weight() * 0.75;
376 
377  }
378 
379  //Now do the WW channels
380  bool passZVeto_WWa=false;
381  bool passPTll_WWa=false;
382  bool passMetRel_WWa=false;
383  bool passMll_WWa=false;
384 
385  bool passMT2_WWb=false;
386  bool passMT2_WWc=false;
387  bool passMll_WWb=false;
388 
389  if(leptonPTCut && mllCut && isOS && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0){
390 
391  //Calculate ETmiss_rel
392  double dPhiMin=9999;
393  for(const HEPUtils::Jet* jet : centralBJets){
394  double dphi=fabs(jet->mom().deltaPhi(ptot));
395  if(fabs(dphi)<dPhiMin)dPhiMin=dphi;
396  }
397 
398  for(const HEPUtils::Jet* jet : centralNonBJets){
399  double dphi=fabs(jet->mom().deltaPhi(ptot));
400  if(fabs(dphi)<dPhiMin)dPhiMin=dphi;
401  }
402 
403  for(const HEPUtils::Particle* lep : signalLeptons){
404  double dphi=fabs(lep->mom().deltaPhi(ptot));
405  if(fabs(dphi)<dPhiMin)dPhiMin=dphi;
406  }
407 
408  double ETmiss_rel=0;
409  if(dPhiMin<(3.14/2)){
410  ETmiss_rel=met*sin(dPhiMin);
411  }
412  else {
413  ETmiss_rel=met;
414  }
415 
416  //Calculate MT2
417  double pa[3] = { 0, signalLeptons[0]->mom().px(), signalLeptons[0]->mom().py() };
418  double pb[3] = { 0, signalLeptons[1]->mom().px(), signalLeptons[1]->mom().py() };
419  double pmiss[3] = { 0, ptot.px(), ptot.py() };
420  double mn = 0.;
421 
422  mt2_bisect::mt2 mt2_calc;
423 
424  mt2_calc.set_momenta(pa,pb,pmiss);
425  mt2_calc.set_mn(mn);
426  double mt2 = mt2_calc.get_mt2();
427 
428  double mll=(signalLeptons[0]->mom() + signalLeptons[1]->mom()).m();
429 
430  //Variables for cutflow debugging
431 
432  if(fabs(mll-91.)>10.)passZVeto_WWa=true;
433  if((signalLeptons[0]->mom() + signalLeptons[1]->mom()).pT() > 80.)passPTll_WWa=true;
434  if(ETmiss_rel > 80.)passMetRel_WWa=true;
435  if(mll < 120.) passMll_WWa=true;
436 
437  if(mt2 > 90.)passMT2_WWb=true;
438  if(mll < 170.)passMll_WWb=true;
439 
440  if(mt2 > 100.)passMT2_WWc=true;
441 
442 
443  if((signalLeptons[0]->mom() + signalLeptons[1]->mom()).pT() > 80. &&
444  ETmiss_rel > 80. &&
445  mll < 120. &&
446  (numElectrons==1 && numMuons==1)) _num_WWa_DF += event->weight() * 0.89;
447 
448  if((signalLeptons[0]->mom() + signalLeptons[1]->mom()).pT() > 80. &&
449  ETmiss_rel > 80. &&
450  mll < 120. &&
451  (numElectrons==2 && fabs(mll-91.)>10.)) _num_WWa_SF += event->weight() * 0.97;
452 
453  if((signalLeptons[0]->mom() + signalLeptons[1]->mom()).pT() > 80. &&
454  ETmiss_rel > 80. &&
455  mll < 120. &&
456  (numMuons==2 && fabs(mll-91.)>10.)) _num_WWa_SF += event->weight() * 0.75;
457 
458  if(mt2 > 90. &&
459  mll < 170. &&
460  (numElectrons==1 && numMuons==1)) _num_WWb_DF += event->weight() * 0.89;
461 
462  if(mt2 > 90. &&
463  mll < 170. &&
464  (numElectrons==2 && fabs(mll-91.)>10.)) _num_WWb_SF += event->weight() * 0.97;
465 
466 
467  if(mt2 > 90. &&
468  mll < 170. &&
469  (numMuons==2 && fabs(mll-91.)>10.)) _num_WWb_SF += event->weight() * 0.75;
470 
471  if(mt2 > 100. && (numElectrons==1 && numMuons==1)) _num_WWc_DF += event->weight() * 0.89;
472 
473  if(mt2 > 100. && (numElectrons==2 && fabs(mll-91.)>10.)) _num_WWc_SF += event->weight() * 0.97;
474 
475  if(mt2 > 100. && (numMuons==2 && fabs(mll-91.)>10.)) _num_WWc_SF += event->weight() * 0.75;
476 
477  }
478 
479  //Finally, do the Z+jets signal region
480 
481  bool passZWindow=true;
482  bool passPTll=true;
483  bool passETmissRel=true;
484  bool passdRll=true;
485  bool passMjj=true;
486  bool passJetPT=true;
487 
488  if(leptonPTCut && mllCut && isOS && tauVeto && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0){
489 
490  double mll=(signalLeptons[0]->mom() + signalLeptons[1]->mom()).m();
491 
492  //Calculate ETmiss_rel
493  double dPhiMin=9999;
494  for(const HEPUtils::Jet* jet : centralBJets){
495  double dphi=jet->mom().deltaPhi(ptot);
496  if(dphi<dPhiMin)dPhiMin=dphi;
497  }
498 
499  for(const HEPUtils::Jet* jet : centralNonBJets){
500  double dphi=jet->mom().deltaPhi(ptot);
501  if(dphi<dPhiMin)dPhiMin=dphi;
502  }
503 
504  for(const HEPUtils::Particle* lep : signalLeptons){
505  double dphi=lep->mom().deltaPhi(ptot);
506  if(dphi<dPhiMin)dPhiMin=dphi;
507  }
508 
509  double ETmiss_rel=0;
510  if(dPhiMin<(3.14/2)){
511  ETmiss_rel=met*sin(dPhiMin);
512  }
513  else {
514  ETmiss_rel=met;
515  }
516 
517 
518  double dRll = signalLeptons[0]->mom().deltaR_eta(signalLeptons[1]->mom());
519 
520  double mjj = (centralNonBJets[0]->mom()+centralNonBJets[1]->mom()).m();
521 
522 
523  //Cuts for cutflow debugging
524  if(ETmiss_rel<=80.)passETmissRel=false;
525  if(fabs(mll-91.)>10.)passZWindow=false;
526  if((signalLeptons[0]->mom()+signalLeptons[1]->mom()).pT()<=80.)passPTll=false;
527  if(!(dRll > 0.3 && dRll < 1.5))passdRll=false;
528  if(!(mjj > 50. && mjj<100.))passMjj=false;
529  if(!(centralNonBJets[0]->pT()>45. && centralNonBJets[1]->pT()>45.))passJetPT=false;
530 
531  if(fabs(mll-91.)<10 && ETmiss_rel>80. && (signalLeptons[0]->mom()+signalLeptons[1]->mom()).pT()>80. && dRll > 0.3 && dRll < 1.5 && mjj > 50. && mjj<100. && passJetPT && (numElectrons==2 && numMuons==0)) _num_Zjets += event->weight() * 0.97;
532 
533  if(fabs(mll-91.)<10 && ETmiss_rel>80. && (signalLeptons[0]->mom()+signalLeptons[1]->mom()).pT()>80. && dRll > 0.3 && dRll < 1.5 && mjj > 50. && mjj<100. && passJetPT && (numElectrons==0 && numMuons==2)) _num_Zjets += event->weight() * 0.75;
534 
535  }
536 
537  cutFlowVector_str[0] = "No cuts ";
538  cutFlowVector_str[1] = "2 electrons ";
539  cutFlowVector_str[2] = "Lepton pT cuts (trigger) ";
540  cutFlowVector_str[3] = "mll cuts ";
541  cutFlowVector_str[4] = "OS leptons ";
542  cutFlowVector_str[5] = "tau veto ";
543  cutFlowVector_str[6] = "e+e-: Jet veto ";
544  cutFlowVector_str[7] = "e+e-: Z veto ";
545  cutFlowVector_str[8] = "e+e-: SR-MT290 ";
546  cutFlowVector_str[9] = "e+e-: SR-MT2120 ";
547  cutFlowVector_str[10] = "e+e-: SR-MT2150 ";
548  cutFlowVector_str[11] = "mu+mu-: 2 signal leptons ";
549  cutFlowVector_str[12] = "mu+mu-: Jet veto ";
550  cutFlowVector_str[13] = "mu+mu-: Z veto ";
551  cutFlowVector_str[14] = "mu+mu-: SR-MT290 ";
552  cutFlowVector_str[15] = "mu+mu-: SR-MT2120 ";
553  cutFlowVector_str[16] = "mu+mu-: SR-MT2150 ";
554  cutFlowVector_str[17] = "e+-mu-+: 2 signal leptons ";
555  cutFlowVector_str[18] = "e+-mu-+: Jet veto ";
556  cutFlowVector_str[19] = "e+-mu-+: Z veto ";
557  cutFlowVector_str[20] = "e+-mu-+: SR-MT290 ";
558  cutFlowVector_str[21] = "e+-mu-+: SR-MT2120 ";
559  cutFlowVector_str[22] = "e+-mu-+: SR-MT2150 ";
560  cutFlowVector_str[23] = "SRZjets e+e-: 2 signal leptons ";
561  cutFlowVector_str[24] = "SRZjets e+e-: >=2 light jets ";
562  cutFlowVector_str[25] = "SRZjets e+e-: No b and forward jets ";
563  cutFlowVector_str[26] = "SRZjets e+e-: Z window ";
564  cutFlowVector_str[27] = "SRZjets e+e-: pTll > 80 ";
565  cutFlowVector_str[28] = "SRZjets e+e-: ETmissrel ";
566  cutFlowVector_str[29] = "SRZjets e+e-: dRll ";
567  cutFlowVector_str[30] = "SRZjets e+e-: mjj ";
568  cutFlowVector_str[31] = "SRZjets e+e-: jet pT ";
569  cutFlowVector_str[32] = "SRZjets mu+mu-: 2 signal leptons ";
570  cutFlowVector_str[33] = "SRZjets mu+mu-: >=2 light jets ";
571  cutFlowVector_str[34] = "SRZjets mu+mu-: No b and forward jets ";
572  cutFlowVector_str[35] = "SRZjets mu+mu-: Z window ";
573  cutFlowVector_str[36] = "SRZjets mu+mu-: pTll > 80 ";
574  cutFlowVector_str[37] = "SRZjets mu+mu-: ETmissrel ";
575  cutFlowVector_str[38] = "SRZjets mu+mu-: dRll ";
576  cutFlowVector_str[39] = "SRZjets mu+mu-: mjj ";
577  cutFlowVector_str[40] = "SRZjets mu+mu-: jet pT ";
578  cutFlowVector_str[41] = "SRWWa e+e-: 2 leptons ";
579  cutFlowVector_str[42] = "SRWWa e+e-: Jet veto ";
580  cutFlowVector_str[43] = "SRWWa e+e-: Z veto ";
581  cutFlowVector_str[44] = "SRWWa e+e-: pTll ";
582  cutFlowVector_str[45] = "SRWWa e+e-: ETmissrel ";
583  cutFlowVector_str[46] = "SRWWa e+e-: mll ";
584  cutFlowVector_str[47] = "SRWWa mu+mu-: 2 leptons ";
585  cutFlowVector_str[48] = "SRWWa mu+mu-: Jet veto ";
586  cutFlowVector_str[49] = "SRWWa mu+mu-: Z veto ";
587  cutFlowVector_str[50] = "SRWWa mu+mu-: pTll ";
588  cutFlowVector_str[51] = "SRWWa mu+mu-: ETmissrel ";
589  cutFlowVector_str[52] = "SRWWa mu+mu-: mll ";
590  cutFlowVector_str[53] = "SRWWa e+mu-: 2 leptons ";
591  cutFlowVector_str[54] = "SRWWa e+mu-: Jet veto ";
592  cutFlowVector_str[55] = "SRWWa e+mu-: pTll ";
593  cutFlowVector_str[56] = "SRWWa e+mu-: ETmissrel ";
594  cutFlowVector_str[57] = "SRWWa e+mu-: mll ";
595  cutFlowVector_str[58] = "SRWWb e+e-: 2 leptons ";
596  cutFlowVector_str[59] = "SRWWb e+e-: Jet veto ";
597  cutFlowVector_str[60] = "SRWWb e+e-: Z veto ";
598  cutFlowVector_str[61] = "SRWWb e+e-: mT2 > 90 ";
599  cutFlowVector_str[62] = "SRWWb e+e-: mll < 170 ";
600  cutFlowVector_str[63] = "SRWWb mu+mu-: 2 leptons ";
601  cutFlowVector_str[64] = "SRWWb mu+mu-: Jet veto ";
602  cutFlowVector_str[65] = "SRWWb mu+mu-: Z veto ";
603  cutFlowVector_str[66] = "SRWWb mu+mu-: mT2 > 90 ";
604  cutFlowVector_str[67] = "SRWWb mu+mu-: mll < 170 ";
605  cutFlowVector_str[68] = "SRWWb e+mu-: 2 leptons ";
606  cutFlowVector_str[69] = "SRWWb e+mu-: Jet veto ";
607  cutFlowVector_str[70] = "SRWWb e+mu-: mT2 > 90 ";
608  cutFlowVector_str[71] = "SRWWb e+mu-: mll < 170 ";
609  cutFlowVector_str[72] = "SRWWc e+e-: 2 leptons ";
610  cutFlowVector_str[73] = "SRWWc e+e-: Jet veto ";
611  cutFlowVector_str[74] = "SRWWc e+e-: Z veto ";
612  cutFlowVector_str[75] = "SRWWc e+e-: mT2 > 100 ";
613  cutFlowVector_str[76] = "SRWWc mu+mu-: 2 leptons ";
614  cutFlowVector_str[77] = "SRWWc mu+mu-: Jet veto ";
615  cutFlowVector_str[78] = "SRWWc mu+mu-: Z veto ";
616  cutFlowVector_str[79] = "SRWWc mu+mu-: mT2 > 100 ";
617  cutFlowVector_str[80] = "SRWWc e+mu-: 2 leptons ";
618  cutFlowVector_str[81] = "SRWWc e+mu-: Jet veto ";
619  cutFlowVector_str[82] = "SRWWc e+mu-: mT2 > 100 ";
620 
621  for(int j=0;j<NCUTS;j++){
622  if(j>=0 && j<=10)cutFlowIncrements[j]=0.97;
623  if(j>=11 && j<=16)cutFlowIncrements[j]=0.75;
624  if(j>=17 && j<=22)cutFlowIncrements[j]=0.89;
625 
626  if(j>=23 && j<=31)cutFlowIncrements[j]=0.97;
627  if(j>=32 && j<=40)cutFlowIncrements[j]=0.75;
628 
629  if(j>=41 && j<=46)cutFlowIncrements[j]=0.97;
630  if(j>=47 && j<=52)cutFlowIncrements[j]=0.75;
631  if(j>=53 && j<=57)cutFlowIncrements[j]=0.89;
632 
633  if(j>=58 && j<=62)cutFlowIncrements[j]=0.97;
634  if(j>=63 && j<=67)cutFlowIncrements[j]=0.75;
635  if(j>=68 && j<=71)cutFlowIncrements[j]=0.89;
636 
637  if(j>=72 && j<=75)cutFlowIncrements[j]=0.97;
638  if(j>=76 && j<=79)cutFlowIncrements[j]=0.75;
639  if(j>=80 && j<=82)cutFlowIncrements[j]=0.89;
640  }
641 
642 
643  for(int j=0;j<NCUTS;j++){
644  if( (j==0) ||
645 
646  (j==1 && signalElectrons.size()==2) ||
647 
648  (j==2 && signalElectrons.size()==2 && leptonPTCut) ||
649 
650  (j==3 && signalElectrons.size()==2 && leptonPTCut && mllCut) ||
651 
652  (j==4 && signalElectrons.size()==2 && leptonPTCut && mllCut && isOS) ||
653 
654  (j==5 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto) ||
655 
656  (j==6 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
657 
658  (j==7 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto) ||
659 
660  (j==8 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT290) ||
661 
662  (j==9 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2120) ||
663 
664  (j==10 && leptonPTCut && mllCut && isOS && signalElectrons.size()==2 && tauVeto && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2150) ||
665 
666  //mumu MT2 regions
667 
668  (j==11 && leptonPTCut && mllCut && isOS && signalMuons.size()==2 && tauVeto) ||
669 
670  (j==12 && leptonPTCut && mllCut && isOS && tauVeto && signalMuons.size()==2 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
671 
672  (j==13 && leptonPTCut && mllCut && isOS && tauVeto && signalMuons.size()==2 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto) ||
673 
674  (j==14 && leptonPTCut && mllCut && isOS && tauVeto && signalMuons.size()==2 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT290) ||
675 
676  (j==15 && leptonPTCut && mllCut && isOS && tauVeto && signalMuons.size()==2 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2120) ||
677 
678  (j==16 && leptonPTCut && mllCut && isOS && tauVeto && signalMuons.size()==2 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2150) ||
679 
680  //emu MT2 regions
681 
682  (j==17 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && tauVeto) ||
683 
684  (j==18 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
685 
686  (j==19 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto) ||
687 
688  (j==20 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT290) ||
689 
690  (j==21 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2120) ||
691 
692  (j==22 && leptonPTCut && mllCut && isOS && tauVeto && signalElectrons.size()==1 && signalMuons.size()==1 && (signalElectrons[0]->pid()*signalMuons[0]->pid())<0 && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto && cut_SRMT2150) ||
693 
694  //Start SR Z jets e+e-
695  (j==23 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0)) ||
696 
697  (j==24 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2) ||
698 
699  (j==25 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0) ||
700 
701  (j==26 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow) ||
702 
703  (j==27 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll) ||
704 
705  (j==28 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel) ||
706 
707  (j==29 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll) ||
708 
709  (j==30 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll && passMjj) ||
710 
711  (j==31 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll && passMjj && passJetPT) ||
712 
713  //Start SR Z jets mu+mu-
714  (j==32 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2)) ||
715 
716  (j==33 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2) ||
717 
718  (j==34 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0) ||
719 
720  (j==35 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow) ||
721 
722  (j==36 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll) ||
723 
724  (j==37 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel) ||
725 
726  (j==38 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll) ||
727 
728  (j==39 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll && passMjj) ||
729 
730  (j==40 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets>=2 && numCentralBJets==0 && numForwardJets==0 && passZWindow && passPTll && passETmissRel && passdRll && passMjj && passJetPT) ||
731 
732  //Now start WWa e+e-
733  (j==41 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0)) ||
734 
735  (j==42 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
736 
737  (j==43 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
738 
739  (j==44 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa) ||
740 
741  (j==45 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa) ||
742 
743  (j==46 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa && passMll_WWa) ||
744 
745  //Now start WWa mu+mu-
746  (j==47 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2)) ||
747 
748  (j==48 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
749 
750  (j==49 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
751 
752  (j==50 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa) ||
753 
754  (j==51 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa) ||
755 
756  (j==52 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa && passMll_WWa) ||
757 
758  //Now start WWa e+mu-
759  (j==53 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1)) ||
760 
761  (j==54 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
762 
763  (j==55 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa) ||
764 
765  (j==56 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa) ||
766 
767  (j==57 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passPTll_WWa && passMetRel_WWa && passMll_WWa) ||
768 
769  //WWb e+ e-
770  (j==58 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0)) ||
771 
772  (j==59 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
773 
774  (j==60 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
775 
776  (j==61 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb) ||
777 
778  (j==62 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb && passMll_WWb) ||
779 
780  //WWb mu+ mu-
781  (j==63 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2)) ||
782 
783  (j==64 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
784 
785  (j==65 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
786 
787  (j==66 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb) ||
788 
789  (j==67 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb && passMll_WWb) ||
790 
791  //WWb e+mu-
792 
793  (j==68 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1)) ||
794 
795  (j==69 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
796 
797  (j==70 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb) ||
798 
799  (j==71 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWb && passMll_WWb) ||
800 
801  //WWc e+ e-
802  (j==72 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0)) ||
803 
804  (j==73 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
805 
806  (j==74 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
807 
808  (j==75 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==2 && numMuons==0) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWc) ||
809 
810  //WWc mu+ mu-
811  (j==76 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2)) ||
812 
813  (j==77 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
814 
815  (j==78 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa) ||
816 
817  (j==79 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==0 && numMuons==2) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWc) ||
818 
819  //WWc e+mu-
820 
821  (j==80 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1)) ||
822 
823  (j==81 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0) ||
824 
825  (j==82 && tauVeto && leptonPTCut && mllCut && isOS && (numElectrons==1 && numMuons==1) && numCentralNonBJets==0 && numCentralBJets==0 && numForwardJets==0 && passZVeto_WWa && passMT2_WWc)
826 
827  )cutFlowVector[j]=cutFlowVector[j]+cutFlowIncrements[j];
828 
829  }
830  return;
831  }
832 
834  void combine(const Analysis* other)
835  {
836  const Analysis_ATLAS_8TeV_2LEPEW_20invfb* specificOther
837  = dynamic_cast<const Analysis_ATLAS_8TeV_2LEPEW_20invfb*>(other);
838 
839  for (int j=0; j<NCUTS; j++)
840  {
841  cutFlowVector[j] += specificOther->cutFlowVector[j];
842  cutFlowIncrements[j] += specificOther->cutFlowIncrements[j];
843  cutFlowVector_str[j] = specificOther->cutFlowVector_str[j];
844  }
845 
846  _num_MT2_90_SF += specificOther->_num_MT2_90_SF;
847  _num_MT2_90_DF += specificOther->_num_MT2_90_DF;
848  _num_MT2_120_SF += specificOther->_num_MT2_120_SF;
849  _num_MT2_120_DF += specificOther->_num_MT2_120_DF;
850  _num_MT2_150_SF += specificOther->_num_MT2_150_SF;
851  _num_MT2_150_DF += specificOther->_num_MT2_150_DF;
852  _num_WWa_SF += specificOther->_num_WWa_SF;
853  _num_WWa_DF += specificOther->_num_WWa_DF;
854  _num_WWb_SF += specificOther->_num_WWb_SF;
855  _num_WWb_DF += specificOther->_num_WWb_DF;
856  _num_WWc_SF += specificOther->_num_WWc_SF;
857  _num_WWc_DF += specificOther->_num_WWc_DF;
858  _num_Zjets += specificOther->_num_Zjets;
859  }
860 
861 
862  void collect_results() {
863 
864  // add_result(SignalRegionData("SR label", n_obs, {n_sig_MC, n_sig_MC_sys}, {n_bkg, n_bkg_err}));
865 
866  add_result(SignalRegionData("MT2_90_SF", 33., {_num_MT2_90_SF, 0.}, {38.2, 5.1}));
867  add_result(SignalRegionData("MT2_90_DF", 21., {_num_MT2_90_DF, 0.}, {23.3, 3.7}));
868  add_result(SignalRegionData("MT2_120_SF", 5., {_num_MT2_120_SF, 0.}, {8.9, 2.1}));
869  add_result(SignalRegionData("MT2_120_DF", 5., {_num_MT2_120_DF, 0.}, {3.6, 1.2}));
870  add_result(SignalRegionData("MT2_150_SF", 3., {_num_MT2_150_SF, 0.}, {3.2, 0.7}));
871  add_result(SignalRegionData("MT2_150_DF", 2., {_num_MT2_150_DF, 0.}, {1.0, 0.5}));
872  add_result(SignalRegionData("WWa_SF", 73., {_num_WWa_SF, 0.}, {86.5, 7.4}));
873  add_result(SignalRegionData("WWa_DF", 70., {_num_WWa_DF, 0.}, {73.6, 7.9}));
874  add_result(SignalRegionData("WWb_SF", 26., {_num_WWb_SF, 0.}, {30.2, 3.5}));
875  add_result(SignalRegionData("WWb_DF", 17., {_num_WWb_DF, 0.}, {18.1, 2.6}));
876  add_result(SignalRegionData("WWc_SF", 10., {_num_WWc_SF, 0.}, {20.3, 3.5}));
877  add_result(SignalRegionData("WWc_DF", 11., {_num_WWc_DF, 0.}, {9.0, 2.2}));
878  add_result(SignalRegionData("Zjets", 1., {_num_Zjets, 0.}, {1.4, 0.6}));
879 
880  return;
881  }
882 
883 
884  protected:
885  void analysis_specific_reset() {
886  _num_MT2_90_SF=0;
887  _num_MT2_90_DF=0;
888  _num_MT2_120_SF=0;
889  _num_MT2_120_DF=0;
890  _num_MT2_150_SF=0;
891  _num_MT2_150_DF=0;
892  _num_WWa_SF=0;
893  _num_WWa_DF=0;
894  _num_WWb_SF=0;
895  _num_WWb_DF=0;
896  _num_WWc_SF=0;
897  _num_WWc_DF=0;
898  _num_Zjets=0;
899 
900  std::fill(cutFlowVector.begin(), cutFlowVector.end(), 0);
901  }
902 
903  };
904 
905 
906  DEFINE_ANALYSIS_FACTORY(ATLAS_8TeV_2LEPEW_20invfb)
907 
908 
909  }
910 }
bool has_tag(const HEPUtils::BinnedFn2D< double > &effmap, double eta, double pt)
Randomly get a tag result (can be anything) from a 2D |eta|-pT efficiency map.
Definition: Utils.hpp:207
void run(bool, bool, bool, int, double, double, int, int, int, int, int, double, char[], int, int[], bool, bool, bool, bool, double, int, double(*)(double *, int, int, void *), void(*)(int, int, int, double *, double *, double *, double, double, double, void *), void *)
STL namespace.
START_MODEL b
Definition: demo.hpp:270
void set_momenta(double *pa0, double *pb0, double *pmiss0)
Definition: mt2_bisect.cpp:62
A class for collider analyses within ColliderBit.
Definition: Analysis.hpp:41
bool sortByPT_2lep(const HEPUtils::Particle *lep1, const HEPUtils::Particle *lep2)
double get_mt2()
Definition: mt2_bisect.cpp:50
A simple container for the result of one signal region from one analysis.
void applyTightIDElectronSelection(std::vector< const HEPUtils::Particle *> &electrons)
Efficiency function for Tight ID electrons.
void applyMuonEff(std::vector< const HEPUtils::Particle *> &muons)
Randomly filter the supplied particle list by parameterised muon efficiency.
void set_mn(double mn)
Definition: mt2_bisect.cpp:130
def combine(region_file, pip_file)
Definition: colouring.py:169
DEFINE_ANALYSIS_FACTORY(ATLAS_13TeV_ZGammaGrav_CONFNOTE_80invfb)
void applyElectronEff(std::vector< const HEPUtils::Particle *> &electrons)
Randomly filter the supplied particle list by parameterised electron efficiency.
TODO: see if we can use this one:
Definition: Analysis.hpp:33
void applyTauEfficiencyR1(std::vector< const HEPUtils::Particle *> &taus)
Randomly filter the supplied particle list by parameterised Run 1 tau efficiency. ...
Class for ColliderBit analyses.
Functions that do super fast ATLAS detector simulation based on four-vector smearing.