Rivet Analyses Reference

ATLAS_2012_I1180197

Search for supersymmetry at 7 TeV in final states with jets, missing transverse momentum and isolated leptons with the ATLAS detector.
Experiment: ATLAS (LHC)
Inspire ID: 1180197
Status: UNVALIDATED
Authors:
  • Peter Richardson
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • BSM signal events at 7000 GeV.

One and two lepton search for supersymmmetric particles by ATLAS at 7 TeV. Event counts in the signal regions are implemented as one-bin histograms. Histograms for effective mass are implemented for the two signal hard lepton signal regions and the ratio of missing transverse energy to effective mass for the soft lepton region. Only the one lepton plots are currently implemented as taken from the CONF note.

Source code: ATLAS_2012_I1180197.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
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/VisibleFinalState.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  class ATLAS_2012_I1180197 : public Analysis {
  public:

    /// Constructor
    ATLAS_2012_I1180197()
      : Analysis("ATLAS_2012_I1180197")
    {    }


    /// @name Analysis methods
    //@{

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

      // projection to find the electrons
      IdentifiedFinalState elecs(Cuts::abseta < 2.47 && Cuts::pT > 7*GeV);
      elecs.acceptIdPair(PID::ELECTRON);
      declare(elecs, "elecs");

      // projection to find the muons
      IdentifiedFinalState muons(Cuts::abseta < 2.4 && Cuts::pT > 6*GeV);
      muons.acceptIdPair(PID::MUON);
      declare(muons, "muons");

      // Jet finder
      VetoedFinalState vfs;
      vfs.addVetoPairId(PID::MUON);
      declare(FastJets(vfs, FastJets::ANTIKT, 0.4), "AntiKtJets04");

      // all tracks (to do deltaR with leptons)
      declare(ChargedFinalState(Cuts::abseta < 3 && Cuts::pT > 0.5*GeV), "cfs");

      // for pTmiss
      declare(VisibleFinalState(Cuts::abseta < 4.9), "vfs");

      // Book histograms
      book(_count_1l_3jet_all_channel  ,"count_1l_3jet_all_channel", 1, 0., 1.);
      book(_count_1l_3jet_e_channel    ,"count_1l_3jet_e_channel"  , 1, 0., 1.);
      book(_count_1l_3jet_mu_channel   ,"count_1l_3jet_mu_channel" , 1, 0., 1.);
      book(_count_1l_4jet_all_channel  ,"count_1l_4jet_all_channel", 1, 0., 1.);
      book(_count_1l_4jet_e_channel    ,"count_1l_4jet_e_channel"  , 1, 0., 1.);
      book(_count_1l_4jet_mu_channel   ,"count_1l_4jet_mu_channel" , 1, 0., 1.);
      book(_count_1l_soft_all_channel  ,"count_1l_soft_all_channel", 1, 0., 1.);
      book(_count_1l_soft_e_channel    ,"count_1l_soft_e_channel"  , 1, 0., 1.);
      book(_count_1l_soft_mu_channel   ,"count_1l_soft_mu_channel" , 1, 0., 1.);

      book(_count_2l_2jet_all_channel  ,"count_2l_2jet_all_channel" , 1, 0., 1.);
      book(_count_2l_2jet_ee_channel   ,"count_2l_2jet_ee_channel"  , 1, 0., 1.);
      book(_count_2l_2jet_emu_channel  ,"count_2l_2jet_emu_channel" , 1, 0., 1.);
      book(_count_2l_2jet_mumu_channel ,"count_2l_2jet_mumu_channel", 1, 0., 1.);
      book(_count_2l_4jet_all_channel  ,"count_2l_4jet_all_channel" , 1, 0., 1.);
      book(_count_2l_4jet_ee_channel   ,"count_2l_4jet_ee_channel"  , 1, 0., 1.);
      book(_count_2l_4jet_emu_channel  ,"count_2l_4jet_emu_channel" , 1, 0., 1.);
      book(_count_2l_4jet_mumu_channel ,"count_2l_4jet_mumu_channel", 1, 0., 1.);
      book(_hist_1l_m_eff_3jet        ,"hist_1l_m_eff_3jet"       ,  6, 400., 1600.);
      book(_hist_1l_m_eff_4jet        ,"hist_1l_m_eff_4jet"       ,  4, 800., 1600.);
      book(_hist_1l_eTmiss_m_eff_soft ,"hist_1l_eTmiss_m_eff_soft",  6, 0.1 , 0.7  );
      book(_hist_2l_m_eff_2jet        ,"hist_2l_m_eff_2jet"       ,  5, 700., 1700.);
      book(_hist_2l_m_eff_4jet        ,"hist_2l_m_eff_4jet"       ,  5, 600., 1600.);
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      const double weight = 1.0;

      // get the candiate jets
      Jets cand_jets;
      for ( const Jet& jet :
                apply<FastJets>(event, "AntiKtJets04").jetsByPt(20.0*GeV) ) {
        if ( fabs( jet.eta() ) < 4.5 ) {
          cand_jets.push_back(jet);
        }
      }
      // charged tracks for isolation
      Particles chg_tracks =
        apply<ChargedFinalState>(event, "cfs").particles();
      // find the electrons
      Particles cand_soft_e,cand_hard_e;
      for( const Particle & e :
               apply<IdentifiedFinalState>(event, "elecs").particlesByPt()) {
        double pT  = e.pT();
        double eta = e.eta();
        // remove any leptons within 0.4 of any candidate jets
        bool e_near_jet = false;
        for ( const Jet& jet : cand_jets ) {
          double dR = deltaR(e.momentum(),jet.momentum());
          if ( inRange(dR, 0.2, 0.4) ) {
            e_near_jet = true;
            break;
          }
        }
        if ( e_near_jet ) continue;
        // soft selection
        if(pT>7.&&!(fabs(eta)>1.37&&fabs(eta) < 1.52)) {
          cand_soft_e.push_back(e);
        }
        // hard selection
        if(pT>10.) cand_hard_e.push_back(e);
      }
      Particles cand_soft_mu,cand_hard_mu;
      for( const Particle & mu :
               apply<IdentifiedFinalState>(event, "muons").particlesByPt()) {
        double pT  = mu.pT();
        double eta = mu.eta();
        // remove any leptons within 0.4 of any candidate jets
        bool mu_near_jet = false;
        for ( const Jet& jet : cand_jets ) {
          if ( deltaR(mu.momentum(),jet.momentum()) < 0.4 ) {
            mu_near_jet = true;
            break;
          }
        }
        if ( mu_near_jet ) continue;
        // soft selection
        if (pT > 6*GeV && !inRange(fabs(eta), 1.37, 1.52)) {
          cand_soft_mu.push_back(mu);
        }
        // hard selection
        if (pT > 10*GeV) cand_hard_mu.push_back(mu);
      }
      // pTcone around muon track (hard)
      Particles recon_hard_mu;
      for ( const Particle & mu : cand_hard_mu ) {
        double pTinCone = -mu.pT();
        for ( const Particle & track : chg_tracks ) {
          if ( deltaR(mu.momentum(),track.momentum()) < 0.2 )
            pTinCone += track.pT();
        }
        if ( pTinCone < 1.8*GeV ) recon_hard_mu.push_back(mu);
      }
      // pTcone around muon track (soft)
      Particles recon_soft_mu;
      for ( const Particle & mu : cand_soft_mu ) {
        double pTinCone = -mu.pT();
        if(-pTinCone>20.) continue;
        for ( const Particle & track : chg_tracks ) {
          if ( deltaR(mu.momentum(),track.momentum()) < 0.2 )
            pTinCone += track.pT();
        }
        if ( pTinCone < 1.8*GeV ) recon_soft_mu.push_back(mu);
      }
      // pTcone around electron track (hard)
      Particles recon_hard_e;
      for ( const Particle & e : cand_hard_e ) {
        double pTinCone = -e.pT();
        for ( const Particle & track : chg_tracks ) {
          if ( deltaR(e.momentum(),track.momentum()) < 0.2 )
            pTinCone += track.pT();
        }
        if ( pTinCone < 0.1 * e.pT() ) recon_hard_e.push_back(e);
      }
      // pTcone around electron track (soft)
      Particles recon_soft_e;
      for ( const Particle & e : cand_soft_e ) {
        double pTinCone = -e.pT();
        if(-pTinCone>25.) continue;
        for ( const Particle & track : chg_tracks ) {
          if ( deltaR(e.momentum(),track.momentum()) < 0.2 )
            pTinCone += track.pT();
        }
        if ( pTinCone < 0.1 * e.pT() ) recon_soft_e.push_back(e);
      }

      // pTmiss
      FourMomentum pTmiss;
      for ( const Particle & p :
                apply<VisibleFinalState>(event, "vfs").particles() ) {
        pTmiss -= p.momentum();
      }
      double eTmiss = pTmiss.pT();

      // hard lepton selection
      if( ! recon_hard_e.empty() || !recon_hard_mu.empty() ) {
        // discard jets that overlap with electrons
        Jets recon_jets;
        for ( const Jet& jet : cand_jets ) {
          if(jet.abseta()>2.5||
             jet.pT() < 25*GeV) continue;
          bool away_from_e = true;
          for ( const Particle & e : cand_hard_e ) {
            if ( deltaR(e.momentum(),jet.momentum()) < 0.2 ) {
              away_from_e = false;
              break;
            }
          }
          if ( away_from_e ) recon_jets.push_back( jet );
        }
        // both selections require at least 2 jets
        // meff calculation
        double HT=0.;
        for( const Jet & jet : recon_jets) {
          HT += jet.pT();
        }
        double m_eff_inc  = HT+eTmiss;
        unsigned int njet = recon_jets.size();
        // 1 lepton only
        if( recon_hard_e.size() + recon_hard_mu.size() == 1 && njet >=3 )  {
          // get the lepton
          Particle lepton = recon_hard_e.empty() ?
            recon_hard_mu[0] : recon_hard_e[0];
          // lepton variables
          double pT = lepton.pT();
          double mT  = 2.*(pT*eTmiss -
                           lepton.px()*pTmiss.px() -
                           lepton.py()*pTmiss.py());
          mT = sqrt(mT);
          HT += pT;
          m_eff_inc += pT;
          // apply the cuts on the leptons and min no. of jets
          if( ( ( lepton.abspid() == PID::ELECTRON && pT > 25. ) ||
                ( lepton.abspid() == PID::MUON     && pT > 20. ) ) &&
              mT > 100. && eTmiss > 250. ) {
            double m_eff = pT+eTmiss;
            for (size_t ix = 0; ix < 3; ++ix)
              m_eff += recon_jets[ix].pT();
            // 3 jet channel
            if ( (njet == 3 || recon_jets[3].pT() < 80*GeV ) &&
                recon_jets[0].pT() > 100*GeV ) {
              if (eTmiss/m_eff > 0.3) {
                if (m_eff_inc > 1200*GeV) {
                  _count_1l_3jet_all_channel->fill(0.5,weight);
                  if (lepton.abspid() == PID::ELECTRON )
                    _count_1l_3jet_e_channel->fill(0.5, weight);
                  else
                    _count_1l_3jet_mu_channel->fill(0.5, weight);
                }
                _hist_1l_m_eff_3jet->fill(min(1599., m_eff_inc), weight);
              }
            }
            // 4 jet channel
            else if (njet >=4 && recon_jets[3].pT() > 80*GeV) {
              m_eff += recon_jets[3].pT();
              if (eTmiss/m_eff>0.2) {
                if (m_eff_inc > 800*GeV) {
                  _count_1l_4jet_all_channel->fill(0.5, weight);
                  if(lepton.abspid() == PID::ELECTRON )
                    _count_1l_4jet_e_channel->fill(0.5, weight);
                  else
                    _count_1l_4jet_mu_channel->fill(0.5, weight);
                }
                _hist_1l_m_eff_4jet->fill(min(1599., m_eff_inc), weight);
              }
            }
          }
        }
        // multi lepton
        else if( recon_hard_e.size() + recon_hard_mu.size() >= 2 && njet >=2 ) {
          // get all the leptons and sort them by pT
          Particles leptons(recon_hard_e.begin(),recon_hard_e.end());
          leptons.insert(leptons.begin(),recon_hard_mu.begin(),recon_hard_mu.end());
          std::sort(leptons.begin(), leptons.end(), cmpMomByPt);
          double m_eff(0.0);
          for (size_t ix = 0; ix < leptons.size(); ++ix)
            m_eff += leptons[ix].pT();
          m_eff_inc += m_eff;
          m_eff += eTmiss;
          for (size_t ix = 0; ix < (size_t) min(4, int(recon_jets.size())); ++ix)
            m_eff += recon_jets[ix].pT();
          // require opposite sign leptons
          if (leptons[0].pid()*leptons[1].pid()<0) {
            // 2 jet
            if (recon_jets[1].pT()>200 &&
               ( njet<4 || (njet>=4 && recon_jets[3].pT() < 50*GeV)) && eTmiss > 300*GeV) {
              _count_2l_2jet_all_channel->fill(0.5, weight);
              if (leptons[0].abspid() == PID::ELECTRON && leptons[1].abspid() == PID::ELECTRON )
                _count_2l_2jet_ee_channel->fill(0.5, weight);
              else if (leptons[0].abspid() == PID::MUON && leptons[1].abspid() == PID::MUON )
                _count_2l_2jet_mumu_channel->fill(0.5, weight);
              else
                _count_2l_2jet_emu_channel->fill(0.5, weight);
              _hist_2l_m_eff_2jet->fill(min(1699., m_eff_inc), weight);
            }
            // 4 jet
            else if (njet >= 4 && recon_jets[3].pT() > 50*GeV &&
                     eTmiss > 100*GeV && eTmiss/m_eff > 0.2) {
              if ( m_eff_inc > 650*GeV ) {
                _count_2l_4jet_all_channel->fill(0.5, weight);
                if (leptons[0].abspid() == PID::ELECTRON && leptons[1].abspid() == PID::ELECTRON )
                  _count_2l_4jet_ee_channel->fill(0.5, weight);
                else if (leptons[0].abspid() == PID::MUON && leptons[1].abspid() == PID::MUON )
                  _count_2l_4jet_mumu_channel->fill(0.5, weight);
                else
                  _count_2l_4jet_emu_channel->fill(0.5, weight);
              }
              _hist_2l_m_eff_4jet->fill(min(1599., m_eff_inc), weight);
            }
          }
        }
      }
      // soft lepton selection
      if ( recon_soft_e.size() + recon_soft_mu.size() == 1 ) {
        // discard jets that overlap with electrons
        Jets recon_jets;
        for ( const Jet& jet : cand_jets ) {
          if (jet.abseta() > 2.5 || jet.pT() < 25*GeV) continue;
          bool away_from_e = true;
          for ( const Particle & e : cand_soft_e ) {
            if ( deltaR(e.momentum(), jet.momentum()) < 0.2 ) {
              away_from_e = false;
              break;
            }
          }
          if ( away_from_e ) recon_jets.push_back( jet );
        }
        // meff calculation
        double HT=0.;
        for (const Jet & jet : recon_jets) {
          HT += jet.pT();
        }
        double m_eff_inc  = HT+eTmiss;
        // get the lepton
        Particle lepton = recon_soft_e.empty() ?
          recon_soft_mu[0] : recon_soft_e[0];
        // lepton variables
        double pT = lepton.pT();
        double mT  = 2.*(pT*eTmiss -
                         lepton.px()*pTmiss.px() -
                         lepton.py()*pTmiss.py());
        mT = sqrt(mT);
        m_eff_inc += pT;
        double m_eff = pT+eTmiss;
        // apply final cuts
        if (recon_jets.size() >= 2 && recon_jets[0].pT()>130*GeV && mT > 100*GeV && eTmiss > 250*GeV) {
          for (size_t ix = 0; ix < 2; ++ix)
            m_eff += recon_jets[0].pT();
          if (eTmiss/m_eff > 0.3) {
            _count_1l_soft_all_channel->fill(0.5, weight);
            if (lepton.abspid() == PID::ELECTRON )
              _count_1l_soft_e_channel->fill(0.5, weight);
            else
              _count_1l_soft_mu_channel->fill(0.5, weight);
          }
          _hist_1l_eTmiss_m_eff_soft->fill( eTmiss/m_eff_inc, weight);
        }
      }
    }


    void finalize() {

      double norm = 4.7* crossSection()/sumOfWeights()/femtobarn;
      scale(_count_1l_3jet_all_channel  ,norm);
      scale(_count_1l_3jet_e_channel    ,norm);
      scale(_count_1l_3jet_mu_channel   ,norm);
      scale(_count_1l_4jet_all_channel  ,norm);
      scale(_count_1l_4jet_e_channel    ,norm);
      scale(_count_1l_4jet_mu_channel   ,norm);
      scale(_count_1l_soft_all_channel  ,norm);
      scale(_count_1l_soft_e_channel    ,norm);
      scale(_count_1l_soft_mu_channel   ,norm);
      scale(_count_2l_2jet_all_channel  ,norm);
      scale(_count_2l_2jet_ee_channel   ,norm);
      scale(_count_2l_2jet_emu_channel  ,norm);
      scale(_count_2l_2jet_mumu_channel ,norm);
      scale(_count_2l_4jet_all_channel  ,norm);
      scale(_count_2l_4jet_ee_channel   ,norm);
      scale(_count_2l_4jet_emu_channel  ,norm);
      scale(_count_2l_4jet_mumu_channel ,norm);

      scale(_hist_1l_m_eff_3jet         ,200.*norm);
      scale(_hist_1l_m_eff_4jet         ,200.*norm);
      scale(_hist_1l_eTmiss_m_eff_soft  ,0.1*norm);
      scale(_hist_2l_m_eff_2jet         ,200.*norm);
      scale(_hist_2l_m_eff_4jet         ,200.*norm);

    }

  private:

    /// @name Histos
    //@{
    Histo1DPtr _count_1l_3jet_all_channel;
    Histo1DPtr _count_1l_3jet_e_channel;
    Histo1DPtr _count_1l_3jet_mu_channel;
    Histo1DPtr _count_1l_4jet_all_channel;
    Histo1DPtr _count_1l_4jet_e_channel;
    Histo1DPtr _count_1l_4jet_mu_channel;
    Histo1DPtr _count_1l_soft_all_channel;
    Histo1DPtr _count_1l_soft_e_channel;
    Histo1DPtr _count_1l_soft_mu_channel;
    Histo1DPtr _count_2l_2jet_all_channel;
    Histo1DPtr _count_2l_2jet_ee_channel;
    Histo1DPtr _count_2l_2jet_emu_channel;
    Histo1DPtr _count_2l_2jet_mumu_channel;
    Histo1DPtr _count_2l_4jet_all_channel;
    Histo1DPtr _count_2l_4jet_ee_channel;
    Histo1DPtr _count_2l_4jet_emu_channel;
    Histo1DPtr _count_2l_4jet_mumu_channel;
    Histo1DPtr _hist_1l_m_eff_3jet;
    Histo1DPtr _hist_1l_m_eff_4jet;
    Histo1DPtr _hist_1l_eTmiss_m_eff_soft;
    Histo1DPtr _hist_2l_m_eff_2jet;
    Histo1DPtr _hist_2l_m_eff_4jet;
    //@}

  };

  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2012_I1180197);

}