Rivet Analyses Reference

ATLAS_2014_I1326641

3-jet cross section with 7 TeV data
Experiment: ATLAS (LHC)
Inspire ID: 1326641
Status: VALIDATED
Authors:
  • Aliaksei Hrynevich
  • Pavel Starovoitov
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • QCD jet production with at least three leading jets

Double-differential three-jet production cross-sections have been measured in proton-proton collisions at a centre-of-mass energy of $\sqrt{s} = 7\,$TeV using the ATLAS detector at the Large Hadron Collider. The measurements are presented as a function of the three-jet mass ($m_{jjj}$), in bins of the sum of the absolute rapidity separations between the three leading jets ($\left|Y^\ast\right|$). Invariant masses extending up to 5\,TeV are reached for $8<\left|Y^\ast\right|<10$. These measurements use a sample of data recorded using the ATLAS detector in 2011, which corresponds to an integrated luminosity of $4.51\,\text{fb}^{-1}$. Jets are identified using the anti-$k_\text{t}$ algorithm with two different jet radius parameters, $R=0.4$ and $R=0.6$.

Source code: ATLAS_2014_I1326641.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


    class ATLAS_2014_I1326641 : public Analysis {
        public:

            /// @name Constructors etc.
            //@{

            /// Constructor
            ATLAS_2014_I1326641()
                : Analysis("ATLAS_2014_I1326641")
            {

            }

            //@}


        public:

            /// @name Analysis methods
            //@{

            /// Book histograms and initialise projections before the run
            void init() {
                //std::cout << " HELLO ANALYSIS : init " << std::'\n';
                const FinalState fs;

                FastJets fj04(fs, FastJets::ANTIKT, 0.4);
                fj04.useInvisibles();
                declare(fj04, "AntiKT04");

                FastJets fj06(fs, FastJets::ANTIKT, 0.6);
                fj06.useInvisibles();
                declare(fj06, "AntiKT06");

                double ystarBins[] = { 0.0, 2.0, 4.0, 6.0, 8.0, 10.0 };

                size_t massDsOffset(0);
                for (size_t alg = 0; alg < 2; ++alg) {
                    for (size_t i = 0; i < 5; ++i) {
                        Histo1DPtr tmp;
                        h_trijet_Mass[alg].add(ystarBins[i], ystarBins[i+1], book(tmp, 1 + massDsOffset, 1, 1));
                        massDsOffset += 1;
                    }
                }
            }

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

                Jets jetAr[2];
                jetAr[AKT4] = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 50*GeV);
                jetAr[AKT6] = apply<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT > 50*GeV);

                const size_t nJets = 3;
                double ptCut[nJets] = { 150., 100., 50.};

                // Loop over jet "radii" used in analysis
                for (size_t alg = 0; alg < 2; ++alg) {
                    // Identify 3jets
                    vector<FourMomentum> leadJets;
                    for (const Jet& jet : jetAr[alg]) {
                        if (jet.absrap() < 3.0 && leadJets.size() < nJets){
                            int filledJets = leadJets.size();
                            if (jet.pT() < ptCut[filledJets])  continue;
                            leadJets.push_back(jet.momentum());
                        }
                    }

                    if (leadJets.size() < nJets) {
                        MSG_DEBUG("Could not find three suitable leading jets");
                        continue;
                    }

                    const double y1 = leadJets[0].rapidity();
                    const double y2 = leadJets[1].rapidity();
                    const double y3 = leadJets[2].rapidity();

                    const double yStar = fabs(y1-y2) + fabs(y2-y3) + fabs(y1-y3);
                    const double m = (leadJets[0] + leadJets[1] + leadJets[2]).mass();
                    h_trijet_Mass[alg].fill(yStar, m, 1.0);
                }
            }


            /// Normalise histograms etc., after the run
            void finalize() {

                const double sf( 2.0 * crossSection() / sumOfWeights() );
                for (size_t alg = 0; alg < 2; ++alg) {
                    h_trijet_Mass[alg].scale(sf, this);
                }

            }

            //@}


        private:

            // Data members like post-cuts event weight counters go here
            enum Alg { AKT4=0, AKT6=1 };

        private:

            // The 3 jets mass spectrum for anti-kt 4 and anti-kt 6 jets (array index is jet type from enum above)
            BinnedHistogram  h_trijet_Mass[2];
    };

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