Rivet analyses referenceATLAS_2017_I1637587Soft-Drop Jet Mass at 13 TeVExperiment: ATLAS (LHC) Inspire ID: 1637587 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
Jet substructure observables have significantly extended the search program for physics beyond the standard model at the Large Hadron Collider. The state-of-the-art tools have been motivated by theoretical calculations, but there has never been a direct comparison between data and calculations of jet substructure observables that are accurate beyond leading-logarithm approximation. Such observables are significant not only for probing the collinear regime of QCD that is largely unexplored at a hadron collider, but also for improving the understanding of jet substructure properties that are used in many studies at the Large Hadron Collider. This Letter documents a measurement of the first jet substructure quantity at a hadron collider to be calculated at next-to-next-to-leading-logarithm accuracy. The normalized, differential cross section is measured as a function of $\text{log}_{10}\rho^2$, where $\rho$ is the ratio of the soft-drop mass to the ungroomed jet transverse momentum. This quantity is measured in dijet events from 32.9fb$^{-1}$ of $\sqrt{s} = 13$ TeV proton-proton collisions recorded by the ATLAS detector. The data are unfolded to correct for detector effects and compared to precise QCD calculations and leading-logarithm particle-level Monte Carlo simulations. Source code: ATLAS_2017_I1637587.cc 1#include "Rivet/Analysis.hh"
2#include "Rivet/Projections/FinalState.hh"
3#include "Rivet/Projections/FastJets.hh"
4
5#include "fastjet/contrib/SoftDrop.hh"
6
7namespace Rivet {
8
9 /// @brief Soft drop mass at 13 TeV
10 class ATLAS_2017_I1637587: public Analysis {
11 public:
12
13 /// Constructor
14 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2017_I1637587);
15
16 /// Book cuts and projections
17 void init() {
18 // All final state particles
19 const FinalState fs(Cuts::abseta < 5.0);
20
21 FastJets jets(fs, JetAlg::ANTIKT, 0.8, JetMuons::NONE, JetInvisibles::NONE);
22 declare(jets, "jets");
23
24 book(_h_Table1, 1,1,1);
25 book(_h_Table2, 2,1,1);
26 book(_h_Table3, 3,1,1);
27
28 book(_h_Table4, 4,1,1);
29 book(_h_Table5, 5,1,1);
30 book(_h_Table6, 6,1,1);
31
32 betas = { 0., 1., 2. };
33 ptBins = { 600, 650, 700, 750, 800, 850, 900, 950, 1000, 2000 };
34 rhoBins = { -4.5, -4.1, -3.7, -3.3, -2.9, -2.5, -2.1, -1.7, -1.3, -0.9, -0.5 };
35 }
36
37 void analyze(const Event& event) {
38
39 const Jets& myJets = apply<FastJets>(event, "jets").jetsByPt(Cuts::pT > 400*GeV);
40
41 if (myJets.size() < 2) vetoEvent;
42 if (myJets[0].pT() > 1.5*myJets[1].pT()) vetoEvent;
43 if (myJets[0].abseta() > 1.5 || myJets[1].abseta() > 1.5) vetoEvent;
44
45 for (size_t i = 0; i < 2; ++i) {
46 if (myJets[i].pT() < 600*GeV) continue;
47 ClusterSequence cs_ca(myJets[i].constituents(), JetDefinition(fastjet::cambridge_algorithm, 0.8));
48 PseudoJets myJet_ca = sorted_by_pt(cs_ca.inclusive_jets(400.0));
49 if(myJet_ca.size()==0) continue;
50 for (size_t ibeta = 0; ibeta < 3; ++ibeta) {
51 fastjet::contrib::SoftDrop sd(betas[ibeta], 0.1); //beta, zcut
52 PseudoJet sdJet = sd(myJet_ca[0]);
53 double rho2 = pow(sdJet.m()/myJets[i].pT(),2);
54 double log10rho2 = log(rho2)/log(10.);
55
56 if (log10rho2 < -4.5) continue;
57 if (ibeta==0) _h_Table1->fill(log10rho2);
58 if (ibeta==1) _h_Table2->fill(log10rho2);
59 if (ibeta==2) _h_Table3->fill(log10rho2);
60
61 if (ibeta==0) _h_Table4->fill(return_bin(rho2, myJets[i].pT()));
62 if (ibeta==1) _h_Table5->fill(return_bin(rho2, myJets[i].pT()));
63 if (ibeta==2) _h_Table6->fill(return_bin(rho2, myJets[i].pT()));
64 }
65 }
66 }
67
68 void finalize() {
69 //Normalization comes here.
70 double norm0 = 0.;
71 double norm1 = 0.;
72 double norm2 = 0.;
73 for (size_t i = 4; i <= 7; ++i) { //only normalize in the resummation region.
74 norm0 += _h_Table1->bin(i+1).sumW()/_h_Table1->bin(i+1).dVol();
75 norm1 += _h_Table2->bin(i+1).sumW()/_h_Table2->bin(i+1).dVol();
76 norm2 += _h_Table3->bin(i+1).sumW()/_h_Table3->bin(i+1).dVol();
77 }
78
79 if (norm0 != 0) {
80 _h_Table1->scaleW(1.0/norm0);
81 } else {
82 MSG_WARNING("Zero entries, cannot normalise Table 1");
83 }
84
85 if (norm1 != 0) {
86 _h_Table2->scaleW(1.0/norm1);
87 } else {
88 MSG_WARNING("Zero entries, cannot normalise Table 2");
89 }
90
91 if (norm2 != 0) {
92 _h_Table3->scaleW(1.0/norm2);
93 } else {
94 MSG_WARNING("Zero entries, cannot normalise Table 3");
95 }
96
97 ptNorm( _h_Table4 );
98 ptNorm( _h_Table5 );
99 ptNorm( _h_Table6 );
100 }
101
102 void ptNorm(Histo1DPtr ptBinnedHist) {
103 for (size_t k = 0; k < 9; ++k){
104 double normalization = 0;
105 for (size_t j = 4; j <= 7; ++j) {
106 double height = ptBinnedHist->bin(k*10 + j+1).sumW();
107 height /= ptBinnedHist->bin(k*10 + j+1).dVol();
108 normalization += height;
109 }
110 if( normalization == 0 ) continue;
111
112 for (size_t j = 0; j < 10; ++j) {
113 ptBinnedHist->bin(k*10 + j+1).scaleW(1. / normalization);
114 }
115 }
116
117 return;
118 }
119
120 size_t return_bin(double rho, double pT){
121 if (pT < 600.) return -1;
122 if (rho < pow(10,-4.5)) return -1;
123
124 size_t pTbin = 1;
125 if (pT < 600) pTbin = 0; //should not happen
126 else if (pT < 650) pTbin = 1;
127 else if (pT < 700) pTbin = 2;
128 else if (pT < 750) pTbin = 3;
129 else if (pT < 800) pTbin = 4;
130 else if (pT < 850) pTbin = 5;
131 else if (pT < 900) pTbin = 6;
132 else if (pT < 950) pTbin = 7;
133 else if (pT < 1000) pTbin = 8;
134 else pTbin = 9;
135
136 size_t rhobin = 1;
137 if (rho < pow(10,-4.5)) rhobin = 0; //this should not happen.
138 else if (rho < pow(10,-4.1)) rhobin = 1;
139 else if (rho < pow(10,-3.7)) rhobin = 2;
140 else if (rho < pow(10,-3.3)) rhobin = 3;
141 else if (rho < pow(10,-2.9)) rhobin = 4;
142 else if (rho < pow(10,-2.5)) rhobin = 5;
143 else if (rho < pow(10,-2.1)) rhobin = 6;
144 else if (rho < pow(10,-1.7)) rhobin = 7;
145 else if (rho < pow(10,-1.3)) rhobin = 8;
146 else if (rho < pow(10,-0.9)) rhobin = 9;
147 else if (rho < pow(10,-0.5)) rhobin = 10;
148 else rhobin = 10;
149
150 return (rhobin-1) + (pTbin-1)*10;
151 }
152
153
154 private:
155 /// Histograms
156 Histo1DPtr _h_Table1, _h_Table2, _h_Table3, _h_Table4, _h_Table5, _h_Table6;
157 vector<double> betas, ptBins, rhoBins;
158
159 };
160
161 RIVET_DECLARE_PLUGIN(ATLAS_2017_I1637587);
162}
|