Rivet Analyses Reference

ATLAS_2012_CONF_2012_103

High jet-multiplicity + MET squark and gluino search
Experiment: ATLAS (LHC)
Status: OBSOLETE
Authors:
  • Peter Richardson
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • BSM signal events at 8000 GeV.

Search for SUSY using events with 6 or more jets in association with missing transverse momentum produced in proton-proton collisions at a centre-of-mass energy of 8 TeV. The data sample has a total integrated luminosity of 5.8 fb$^{-1}$. Distributions in the $W$ and top control regions are not produced, while in addition to the plots from the paper the count of events in the different signal regions is included. The analysis is identical to the previous 7 TeV paper.

Source code: ATLAS_2012_CONF_2012_103.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/VisibleFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/RivetMT2.hh"

namespace Rivet {


  class ATLAS_2012_CONF_2012_103 : public Analysis {
  public:

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


    /// @name Analysis methods
    //@{

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

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

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

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

      VetoedFinalState vfs;
      vfs.addVetoPairId(PID::MUON);

      /// Jet finder
      declare(FastJets(vfs, FastJets::ANTIKT, 0.4), "AntiKtJets04");

      /// Book histograms
      book(_etmiss_HT_7j55 ,"etmiss_HT_7j55", 8, 0., 16.);
      book(_etmiss_HT_8j55 ,"etmiss_HT_8j55", 8, 0., 16.);
      book(_etmiss_HT_9j55 ,"etmiss_HT_9j55", 8, 0., 16.);
      book(_etmiss_HT_6j80 ,"etmiss_HT_6j80", 8, 0., 16.);
      book(_etmiss_HT_7j80 ,"etmiss_HT_7j80", 8, 0., 16.);
      book(_etmiss_HT_8j80 ,"etmiss_HT_8j80", 8, 0., 16.);

      book(_hist_njet55 ,"hist_njet55", 4, 5.5, 9.5);
      book(_hist_njet80 ,"hist_njet80", 4, 4.5, 8.5);

      book(_count_7j55 ,"count_7j55", 1, 0., 1.);
      book(_count_8j55 ,"count_8j55", 1, 0., 1.);
      book(_count_9j55 ,"count_9j55", 1, 0., 1.);
      book(_count_6j80 ,"count_6j80", 1, 0., 1.);
      book(_count_7j80 ,"count_7j80", 1, 0., 1.);
      book(_count_8j80 ,"count_8j80", 1, 0., 1.);

    }


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

      // get the jet candidates
      Jets cand_jets;
      for (const Jet& jet :
               apply<FastJets>(event, "AntiKtJets04").jetsByPt(20.0*GeV) ) {
        if ( fabs( jet.eta() ) < 2.8 ) {
          cand_jets.push_back(jet);
        }
      }

      // candidate muons
      Particles cand_mu =
        apply<IdentifiedFinalState>(event, "muons").particlesByPt();

      // candidate electrons
      Particles cand_e  =
        apply<IdentifiedFinalState>(event, "elecs").particlesByPt();

      // resolve jet/lepton ambiguity
      Jets recon_jets;
      for ( const Jet& jet : cand_jets ) {
        // candidates after |eta| < 2.8
        if ( fabs( jet.eta() ) >= 2.8 ) continue;
        bool away_from_e = true;
        for ( const Particle & e : cand_e ) {
          if ( deltaR(e.momentum(),jet.momentum()) <= 0.2 ) {
            away_from_e = false;
            break;
          }
        }
        if ( away_from_e ) recon_jets.push_back( jet );
      }

      // only keep electrons more than R=0.4 from jets
      Particles recon_e;
      for ( const Particle & e : cand_e ) {
        bool away = true;
        for ( const Jet& jet : recon_jets ) {
          if ( deltaR(e.momentum(),jet.momentum()) < 0.4 ) {
            away = false;
            break;
          }
        }
        if ( away )
          recon_e.push_back( e );
      }

      // only keep muons more than R=0.4 from jets
      Particles recon_mu;
      for ( const Particle & mu : cand_mu ) {
        bool away = true;
        for ( const Jet& jet : recon_jets ) {
          if ( deltaR(mu.momentum(),jet.momentum()) < 0.4 ) {
            away = false;
            break;
          }
        }
        if ( away )
          recon_mu.push_back( mu );
      }

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

      // now only use recon_jets, recon_mu, recon_e

      // reject events with electrons and muons
      if ( ! ( recon_mu.empty() && recon_e.empty() ) ) {
        MSG_DEBUG("Charged leptons left after selection");
        vetoEvent;
      }

      // calculate H_T
      double HT=0;
      for ( const Jet& jet : recon_jets ) {
        if ( jet.pT() > 40 * GeV )
          HT += jet.pT() ;
      }

      // number of jets
      unsigned int njet55=0, njet80=0;
      for (unsigned int ix=0;ix<recon_jets.size();++ix) {
        if(recon_jets[ix].pT()>80.*GeV) ++njet80;
        if(recon_jets[ix].pT()>55.*GeV) ++njet55;
      }

      double ratio = eTmiss/sqrt(HT);

      if(ratio>4.) {
        if(njet55>9) njet55 = 9;
        if(njet80>8) njet80 = 8;
        _hist_njet55->fill(njet55,weight);
        _hist_njet80->fill(njet80,weight);
        // 7j55
        if(njet55>=7)
          _count_7j55->fill( 0.5, weight);
        // 8j55
        if(njet55>=8)
          _count_8j55->fill( 0.5, weight) ;
        // 8j55
        if(njet55==9)
          _count_9j55->fill( 0.5, weight) ;
        // 6j80
        if(njet80>=6)
          _count_6j80->fill( 0.5, weight) ;
        // 7j80
        if(njet80>=7)
          _count_7j80->fill( 0.5, weight) ;
        // 8j80
        if(njet80==8)
          _count_8j80->fill( 0.5, weight) ;
      }

      if(njet55>=7)
        _etmiss_HT_7j55->fill( ratio, weight);
      // 8j55
      if(njet55>=8)
        _etmiss_HT_8j55->fill( ratio, weight) ;
      // 8j55
      if(njet55>=9)
        _etmiss_HT_9j55->fill( ratio, weight) ;
      // 6j80
      if(njet80>=6)
        _etmiss_HT_6j80->fill( ratio, weight) ;
      // 7j80
      if(njet80>=7)
        _etmiss_HT_7j80->fill( ratio, weight) ;
      // 8j80
      if(njet80>=8)
        _etmiss_HT_8j80->fill( ratio, weight) ;

    }

    //@}

    void finalize() {
      double norm = crossSection()/femtobarn*5.8/sumOfWeights();

      scale(_etmiss_HT_7j55,2.*norm);
      scale(_etmiss_HT_8j55,2.*norm);
      scale(_etmiss_HT_9j55,2.*norm);
      scale(_etmiss_HT_6j80,2.*norm);
      scale(_etmiss_HT_7j80,2.*norm);
      scale(_etmiss_HT_8j80,2.*norm);

      scale(_hist_njet55,norm);
      scale(_hist_njet80,norm);

      scale(_count_7j55,norm);
      scale(_count_8j55,norm);
      scale(_count_9j55,norm);
      scale(_count_6j80,norm);
      scale(_count_7j80,norm);
      scale(_count_8j80,norm);
    }

  private:

    /// @name Histograms
    //@{
    Histo1DPtr _etmiss_HT_7j55;
    Histo1DPtr _etmiss_HT_8j55;
    Histo1DPtr _etmiss_HT_9j55;
    Histo1DPtr _etmiss_HT_6j80;
    Histo1DPtr _etmiss_HT_7j80;
    Histo1DPtr _etmiss_HT_8j80;

    Histo1DPtr _hist_njet55;
    Histo1DPtr _hist_njet80;

    Histo1DPtr _count_7j55;
    Histo1DPtr _count_8j55;
    Histo1DPtr _count_9j55;
    Histo1DPtr _count_6j80;
    Histo1DPtr _count_7j80;
    Histo1DPtr _count_8j80;
    //@}

  };

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

}