Rivet Analyses Reference

ATLAS_2016_CONF_2016_078

ATLAS ICHEP16 0-lepton SUSY search at 13 \text{TeV} with 13.2/fb
Experiment: ATLAS (LHC)
Status: VALIDATED
Authors:
  • Andy Buckley
No references listed
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • BSM signal events.

ATLAS search for SUSY in 13 TeV $pp$ collisions at LHC Run 2, using 13.2/fb of integrated luminosity and events containing missing transverse momentum and no isolated high-energy leptons. Detailed info: https://atlas.web.cern.ch/Atlas/GROUPS/PHYSICS/CONFNOTES/ATLAS-CONF-2016-078/

Source code: ATLAS_2016_CONF_2016_078.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/Sphericity.hh"
#include "Rivet/Projections/SmearedParticles.hh"
#include "Rivet/Projections/SmearedJets.hh"
#include "Rivet/Projections/SmearedMET.hh"
#include "Rivet/Tools/Cutflow.hh"

namespace Rivet {


  /// @brief ATLAS 2016 0-lepton SUSY search, from 13/fb ICHEP'16 CONF note
  class ATLAS_2016_CONF_2016_078 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_CONF_2016_078);


    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      // Initialise and register projections
      FinalState calofs(Cuts::abseta < 3.2);
      FastJets fj(calofs, FastJets::ANTIKT, 0.4);
      declare(fj, "TruthJets");
      declare(SmearedJets(fj, JET_SMEAR_ATLAS_RUN2, //JET_BTAG_ATLAS_RUN2_MV2C10
                          [](const Jet& j) {
                            if (j.abseta() > 2.5) return 0.;
                            return j.bTagged(Cuts::pT > 5*GeV) ? 0.77 : j.cTagged(Cuts::pT > 5*GeV) ? 1/6. : 1/134.;
                          }), "RecoJets");

      MissingMomentum mm(calofs);
      declare(mm, "TruthMET");
      declare(SmearedMET(mm, MET_SMEAR_ATLAS_RUN2), "RecoMET");

      PromptFinalState es(Cuts::abseta < 2.47 && Cuts::abspid == PID::ELECTRON, true, true);
      declare(es, "TruthElectrons");
      declare(SmearedParticles(es, ELECTRON_RECOEFF_ATLAS_RUN2, ELECTRON_SMEAR_ATLAS_RUN2), "RecoElectrons");

      PromptFinalState mus(Cuts::abseta < 2.7 && Cuts::abspid == PID::MUON, true);
      declare(mus, "TruthMuons");
      declare(SmearedParticles(mus, MUON_EFF_ATLAS_RUN2, MUON_SMEAR_ATLAS_RUN2), "RecoMuons");


      // Book histograms/counters
      book(_h_2j_0800,"2j-0800");
      book(_h_2j_1200,"2j-1200");
      book(_h_2j_1600,"2j-1600");
      book(_h_2j_2000,"2j-2000");
      book(_h_3j_1200,"2j-2000");
      book(_h_4j_1000,"4j-1000");
      book(_h_4j_1400,"4j-1400");
      book(_h_4j_1800,"4j-1800");
      book(_h_4j_2200,"4j-2200");
      book(_h_4j_2600,"4j-2600");
      book(_h_5j_1400,"5j-1400");
      book(_h_6j_1800,"6j-1800");
      book(_h_6j_2200,"6j-2200");


      // Book cut-flows
      const vector<string> cuts23j = {"Pre-sel+MET+pT1+meff", "Njet", "Dphi_min(j123,MET)", "Dphi_min(j4+,MET)", "pT2", "eta_j12", "MET/sqrtHT", "m_eff(incl)"};
      _flows.addCutflow("2j-0800", cuts23j);
      _flows.addCutflow("2j-1200", cuts23j);
      _flows.addCutflow("2j-1600", cuts23j);
      _flows.addCutflow("2j-2000", cuts23j);
      _flows.addCutflow("3j-1200", cuts23j);
      const vector<string> cuts456j = {"Pre-sel+MET+pT1+meff", "Njet", "Dphi_min(j123,MET)", "Dphi_min(j4+,MET)", "pT4", "eta_j1234", "Aplanarity", "MET/m_eff(Nj)", "m_eff(incl)"};
      _flows.addCutflow("4j-1000", cuts456j);
      _flows.addCutflow("4j-1400", cuts456j);
      _flows.addCutflow("4j-1800", cuts456j);
      _flows.addCutflow("4j-2200", cuts456j);
      _flows.addCutflow("4j-2600", cuts456j);
      _flows.addCutflow("5j-1400", cuts456j);
      _flows.addCutflow("6j-1800", cuts456j);
      _flows.addCutflow("6j-2200", cuts456j);

    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      _flows.fillinit();

      // Same MET cut for all signal regions
      const Vector3 vmet = -apply<SmearedMET>(event, "RecoMET").vectorEt();
      const double met = vmet.mod();
      if (met < 250*GeV) vetoEvent;

      // Get baseline electrons, muons, and jets
      Particles elecs = apply<ParticleFinder>(event, "RecoElectrons").particles(Cuts::pT > 10*GeV);
      Particles muons = apply<ParticleFinder>(event, "RecoMuons").particles(Cuts::pT > 10*GeV);
      Jets jets = apply<JetAlg>(event, "RecoJets").jetsByPt(Cuts::pT > 20*GeV && Cuts::abseta < 2.8); ///< @todo Pile-up subtraction

      // Jet/electron/muons overlap removal and selection
      // Remove electrons within dR = 0.2 of a b-tagged jet
      for (const Jet& j : jets)
        if (j.abseta() < 2.5 && j.pT() > 50*GeV && j.bTagged(Cuts::pT > 5*GeV))
          ifilter_discard(elecs, deltaRLess(j, 0.2, RAPIDITY));
      // Remove any |eta| < 2.8 jet within dR = 0.2 of a remaining electron
      for (const Particle& e : elecs)
        ifilter_discard(jets, deltaRLess(e, 0.2, RAPIDITY));
      // Remove any electron with dR in [0.2, 0.4] of a remaining jet
      for (const Jet& j : jets)
        ifilter_discard(elecs, [&](const Particle& e) { return inRange(deltaR(e,j, RAPIDITY), 0.2, 0.4); });
      // Remove any muon with dR close to a remaining jet, via a functional form
      for (const Jet& j : jets)
        ifilter_discard(muons, [&](const Particle& m) { return deltaR(m,j, RAPIDITY) < min(0.4, 0.04 + 10*GeV/m.pT()); });
      // Remove any |eta| < 2.8 jet within dR = 0.2 of a remaining muon if track conditions are met
      for (const Particle& m : muons)
        /// @todo Add track efficiency random filtering
        ifilter_discard(jets, [&](const Jet& j) {
            if (deltaR(j,m, RAPIDITY) > 0.2) return false;
            const Particles trks = j.particles(Cuts::abscharge > 0 && Cuts::pT > 0.5*GeV);
            return trks.size() < 3 || (m.pT() > 2*j.pT() && m.pT() > 0.7*sum(trks, pT, 0.0));
          });
      // Loose electron selection
      ifilter_select(elecs, ParticleEffFilter(ELECTRON_EFF_ATLAS_RUN2_LOOSE));

      // Veto the event if there are any remaining baseline leptons
      if (!elecs.empty()) vetoEvent;
      if (!muons.empty()) vetoEvent;

      // Passed presel & MET
      _flows.fill(1);

      // Get jets and their pTs
      const Jets jets20 = jets;
      const Jets jets50 = filterBy(jets, Cuts::pT > 50*GeV);
      const size_t njets50 = jets50.size(), njets20 = jets20.size();
      if (jets50.size() < 2) vetoEvent;
      vector<double> jetpts20, jetpts50;
      transform(jets20, jetpts20, pT);
      transform(jets50, jetpts50, pT);

      // Construct multi-jet observables
      const double ht = sum(jetpts20, 0.0);
      const double met_sqrtHT = met / sqrt(ht);
      const double meff_incl = sum(jetpts50, met);
      const double meff_4 = (njets50 >= 4) ? sum(head(jetpts50, 4), met) : -1;
      const double meff_5 = (njets50 >= 5) ? sum(head(jetpts50, 5), met) : -1;
      const double meff_6 = (njets50 >= 6) ? sum(head(jetpts50, 6), met) : -1;
      const double met_meff_4 = met / meff_4;
      const double met_meff_5 = met / meff_5;
      const double met_meff_6 = met / meff_6;

      // Jet |eta|s
      vector<double> jetetas20; transform(jets20, jetetas20, abseta);
      const double etamax_2 = (njets20 >= 2) ? max(head(jetetas20, 2)) : -1;
      const double etamax_4 = (njets20 >= 4) ? max(head(jetetas20, 4)) : -1;
      const double etamax_6 = (njets20 >= 6) ? max(head(jetetas20, 6)) : -1;

      // Get dphis between MET and jets
      vector<double> dphimets50; transform(jets50, dphimets50, deltaPhiWRT(vmet));
      const vector<double> dphimets50_123 = head(dphimets50, 3);
      const vector<double> dphimets50_more = tail(dphimets50, -3);
      const double dphimin_123 = !dphimets50_123.empty() ? min(dphimets50_123) : -1;
      const double dphimin_more = !dphimets50_more.empty() ? min(dphimets50_more) : -1;

      // Jet aplanarity
      Sphericity sph; sph.calc(jets50);
      const double aplanarity = sph.aplanarity();


      //////////////////


      // 2 jet regions
      if (dphimin_123 > 0.8 && dphimin_more > 0.4) {
        if (jetpts50[1] > 200*GeV && etamax_2 < 0.8) { //< implicit pT[0] cut
          if (met_sqrtHT > 14*sqrt(GeV) && meff_incl > 800*GeV) _h_2j_0800->fill();
        }
        if (jetpts50[1] > 250*GeV && etamax_2 < 1.2) { //< implicit pT[0] cut
          if (met_sqrtHT > 16*sqrt(GeV) && meff_incl > 1200*GeV) _h_2j_1200->fill();
          if (met_sqrtHT > 18*sqrt(GeV) && meff_incl > 1600*GeV) _h_2j_1600->fill();
          if (met_sqrtHT > 20*sqrt(GeV) && meff_incl > 2000*GeV) _h_2j_2000->fill();
        }
      }

      // 3 jet region
      if (njets50 >= 3 && dphimin_123 > 0.4 && dphimin_more > 0.2) {
        if (jetpts50[0] > 600*GeV && jetpts50[2] > 50*GeV) { //< implicit pT[1] cut
          if (met_sqrtHT > 16*sqrt(GeV) && meff_incl > 1200*GeV) _h_3j_1200->fill();
        }
      }

      // 4 jet regions (note implicit pT[1,2] cuts)
      if (njets50 >= 4 && dphimin_123 > 0.4 && dphimin_more > 0.4 && jetpts50[0] > 200*GeV && aplanarity > 0.04) {
        if (jetpts50[3] > 100*GeV && etamax_4 < 1.2 && met_meff_4 > 0.25*sqrt(GeV) && meff_incl > 1000*GeV) _h_4j_1000->fill();
        if (jetpts50[3] > 100*GeV && etamax_4 < 2.0 && met_meff_4 > 0.25*sqrt(GeV) && meff_incl > 1400*GeV) _h_4j_1400->fill();
        if (jetpts50[3] > 100*GeV && etamax_4 < 2.0 && met_meff_4 > 0.20*sqrt(GeV) && meff_incl > 1800*GeV) _h_4j_1800->fill();
        if (jetpts50[3] > 150*GeV && etamax_4 < 2.0 && met_meff_4 > 0.20*sqrt(GeV) && meff_incl > 2200*GeV) _h_4j_2200->fill();
        if (jetpts50[3] > 150*GeV &&                   met_meff_4 > 0.20*sqrt(GeV) && meff_incl > 2600*GeV) _h_4j_2600->fill();
      }

      // 5 jet region (note implicit pT[1,2,3] cuts)
      if (njets50 >= 5 && dphimin_123 > 0.4 && dphimin_more > 0.2 && jetpts50[0] > 500*GeV) {
        if (jetpts50[4] > 50*GeV && met_meff_5 > 0.3*sqrt(GeV) && meff_incl > 1400*GeV) _h_5j_1400->fill();
      }

      // 6 jet regions (note implicit pT[1,2,3,4] cuts)
      if (njets50 >= 6 && dphimin_123 > 0.4 && dphimin_more > 0.2 && jetpts50[0] > 200*GeV && aplanarity > 0.08) {
        if (jetpts50[5] >  50*GeV && etamax_6 < 2.0 && met_meff_6*sqrt(GeV) > 0.20 && meff_incl > 1800*GeV) _h_6j_1800->fill();
        if (jetpts50[5] > 100*GeV &&                   met_meff_6*sqrt(GeV) > 0.15 && meff_incl > 2200*GeV) _h_6j_2200->fill();
      }

      // Cutflows
      _flows["2j-0800"].filltail({true, dphimin_123 > 0.8, dphimin_more > 0.4, jetpts50[1] > 200*GeV, etamax_2 < 0.8, met_sqrtHT > 14*sqrt(GeV), meff_incl >  800*GeV});
      _flows["2j-1200"].filltail({true, dphimin_123 > 0.8, dphimin_more > 0.4, jetpts50[1] > 250*GeV, etamax_2 < 1.2, met_sqrtHT > 16*sqrt(GeV), meff_incl > 1200*GeV});
      _flows["2j-1600"].filltail({true, dphimin_123 > 0.8, dphimin_more > 0.4, jetpts50[1] > 250*GeV, etamax_2 < 1.2, met_sqrtHT > 18*sqrt(GeV), meff_incl > 1600*GeV});
      _flows["2j-2000"].filltail({true, dphimin_123 > 0.8, dphimin_more > 0.4, jetpts50[1] > 250*GeV, etamax_2 < 1.2, met_sqrtHT > 20*sqrt(GeV), meff_incl > 2000*GeV});
      _flows["3j-1200"].filltail({njets50 >= 3, dphimin_123 > 0.4, dphimin_more > 0.2, jetpts50[0] > 600*GeV && jetpts50[2] > 50*GeV, true, met_sqrtHT > 16*sqrt(GeV), meff_incl > 1200*GeV});
      _flows["4j-1000"].filltail({njets50 >= 4, dphimin_123 > 0.4, dphimin_more > 0.4, jetpts50[0] > 200*GeV && jetpts50[3] > 100*GeV, etamax_4 < 1.2, aplanarity > 0.04, met_meff_4 > 0.25*sqrt(GeV), meff_incl > 1000*GeV});
      _flows["4j-1400"].filltail({njets50 >= 4, dphimin_123 > 0.4, dphimin_more > 0.4, jetpts50[0] > 200*GeV && jetpts50[3] > 100*GeV, etamax_4 < 2.0, aplanarity > 0.04, met_meff_4 > 0.25*sqrt(GeV), meff_incl > 1400*GeV});
      _flows["4j-1800"].filltail({njets50 >= 4, dphimin_123 > 0.4, dphimin_more > 0.4, jetpts50[0] > 200*GeV && jetpts50[3] > 100*GeV, etamax_4 < 2.0, aplanarity > 0.04, met_meff_4 > 0.20*sqrt(GeV), meff_incl > 1800*GeV});
      _flows["4j-2200"].filltail({njets50 >= 4, dphimin_123 > 0.4, dphimin_more > 0.4, jetpts50[0] > 200*GeV && jetpts50[3] > 150*GeV, etamax_4 < 2.0, aplanarity > 0.04, met_meff_4 > 0.20*sqrt(GeV), meff_incl > 2200*GeV});
      _flows["4j-2600"].filltail({njets50 >= 4, dphimin_123 > 0.4, dphimin_more > 0.4, jetpts50[0] > 200*GeV && jetpts50[3] > 150*GeV, true,           aplanarity > 0.04, met_meff_4 > 0.20*sqrt(GeV), meff_incl > 2600*GeV});
      _flows["5j-1400"].filltail({njets50 >= 5, dphimin_123 > 0.4, dphimin_more > 0.2, jetpts50[0] > 500*GeV && jetpts50[4] > 50*GeV, true, true, met_meff_5 > 0.3*sqrt(GeV), meff_incl > 1400*GeV});
      _flows["6j-1800"].filltail({njets50 >= 6, dphimin_123 > 0.4, dphimin_more > 0.2, jetpts50[0] > 200*GeV && jetpts50[5] >  50*GeV, etamax_6 < 2.0, aplanarity > 0.08, met_meff_6 > 0.20*sqrt(GeV), meff_incl > 1800*GeV});
      _flows["6j-2200"].filltail({njets50 >= 6, dphimin_123 > 0.4, dphimin_more > 0.2, jetpts50[0] > 200*GeV && jetpts50[5] > 100*GeV, true,           aplanarity > 0.08, met_meff_6 > 0.15*sqrt(GeV), meff_incl > 2200*GeV});

    }


    /// Normalise counters after the run
    void finalize() {

      const double sf = 13.3*crossSection()/femtobarn/sumOfWeights();
      scale(_h_2j_0800, sf); scale(_h_2j_1200, sf); scale(_h_2j_1600, sf);
      scale(_h_2j_2000, sf); scale(_h_3j_1200, sf); scale(_h_4j_1000, sf);
      scale(_h_4j_1400, sf); scale(_h_4j_1800, sf); scale(_h_4j_2200, sf);
      scale(_h_4j_2600, sf); scale(_h_5j_1400, sf); scale(_h_6j_1800, sf);
      scale(_h_6j_2200, sf);

      _flows.scale(sf);
      MSG_INFO("CUTFLOWS:\n\n" << _flows);

    }

    //@}


  private:

    /// @name Histograms
    //@{
    CounterPtr _h_2j_0800, _h_2j_1200, _h_2j_1600, _h_2j_2000, _h_3j_1200;
    CounterPtr _h_4j_1000, _h_4j_1400, _h_4j_1800, _h_4j_2200, _h_4j_2600;
    CounterPtr _h_5j_1400, _h_6j_1800, _h_6j_2200;
    //@}

    /// Cut-flows
    Cutflows _flows;

  };

  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2016_CONF_2016_078);
}