Rivet analyses referenceATLAS_2020_I1803608Electroweak Zjj at 13 TeVExperiment: ATLAS (LHC) Inspire ID: 1803608 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
Differential cross-section measurements are presented for the electroweak production of two jets in association with a $ZZ$ boson. These measurements are sensitive to the vector-boson fusion production mechanism and provide a fundamental test of the gauge structure of the Standard Model. The analysis is performed using proton-proton collision data collected by ATLAS at $\sqrt{s} = 13$ TeV and with an integrated luminosity of 139 fb$^{-1}$. The differential cross-sections are measured in the $Z\rightarrow \ell^+\ell^-$ decay channel ($\ell=e,\mu$) as a function of four observables: the dijet invariant mass, the rapidity interval spanned by the two jets, the signed azimuthal angle between the two jets, and the transverse momentum of the dilepton pair. The data are corrected for the effects of detector inefficiency and resolution and are sufficiently precise to distinguish between different state-of-the-art theoretical predictions calculated using Powheg+Pythia8, Herwig7+Vbfnlo and Sherpa 2.2. The differential cross-sections are used to search for anomalous weak-boson self-interactions using a dimension-six effective field theory. The measurement of the signed azimuthal angle between the two jets is found to be particularly sensitive to the interference between the Standard Model and dimension-six scattering amplitudes and provides a direct test of charge-conjugation and parity invariance in the weak-boson self-interactions. Note that the default entry point is for the inclusive Z+2jet selections. For the EW-only measurement use the option TYPE=EW_ONLY. In both cases, electron and muon channels are to be summed. Source code: ATLAS_2020_I1803608.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/PromptFinalState.hh"
5#include "Rivet/Projections/LeptonFinder.hh"
6#include "Rivet/Projections/FastJets.hh"
7
8namespace Rivet {
9
10
11 /// VBFZ in pp at 13 TeV
12 class ATLAS_2020_I1803608 : public Analysis {
13 public:
14
15 /// Constructor
16 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2020_I1803608);
17
18
19 /// @name Analysis methods
20 /// @{
21
22 /// Book histograms and initialise projections before the run
23 void init() {
24 FinalState fs(Cuts::abseta < 5.0);
25
26 PromptFinalState photons(Cuts::abspid == PID::PHOTON);
27 PromptFinalState electrons(Cuts::abspid == PID::ELECTRON);
28 PromptFinalState muons(Cuts::abspid == PID::MUON);
29
30 Cut cuts_el = Cuts::pT > 25*GeV && ( Cuts::abseta < 1.37 || Cuts::absetaIn(1.52, 2.47) );
31 Cut cuts_mu = Cuts::pT > 25*GeV && Cuts::abseta < 2.4;
32
33 LeptonFinder dressed_electrons(electrons, photons, 0.1, cuts_el);
34 declare(dressed_electrons, "DressedElectrons");
35
36 LeptonFinder dressed_muons(muons, photons, 0.1, cuts_mu);
37 declare(dressed_muons, "DressedMuons");
38
39 FastJets jets(fs, JetAlg::ANTIKT, 0.4, JetMuons::NONE, JetInvisibles::NONE);
40 declare(jets, "Jets");
41
42 _doControl = bool(getOption("TYPE") != "EW_ONLY");
43 if (_doControl) {
44 initialisePlots(SRplots, "SR");
45 initialisePlots(CRAplots, "CRA");
46 initialisePlots(CRBplots, "CRB");
47 initialisePlots(CRCplots, "CRC");
48 } else {
49 initialisePlots(SRplots, "EW");
50 }
51 }
52
53
54 /// Perform the per-event analysis
55 void analyze(const Event& event) {
56
57 // Access fiducial electrons and muons
58 const Particle *l1 = nullptr, *l2 = nullptr;
59 Particles muons = apply<LeptonFinder>(event, "DressedMuons").particlesByPt();
60 Particles elecs = apply<LeptonFinder>(event, "DressedElectrons").particlesByPt();
61
62 // Dilepton selection 1: =2 leptons of the same kind
63 if (muons.size()+elecs.size() != 2) vetoEvent;
64 if (muons.size()==2) { l1=&muons[0]; l2=&muons[1]; }
65 else if (elecs.size()==2) { l1=&elecs[0]; l2=&elecs[1]; }
66 else vetoEvent;
67
68 // Dilepton selection 2: oppostie-charge and in mass range
69 if ( !oppCharge(*l1, *l2) ) vetoEvent;
70 if ( !inRange((l1->mom()+l2->mom()).mass()/GeV, 81.0, 101.0) ) vetoEvent;
71
72 // Electron-jet overlap removal (note: muons are not included in jet finding)
73 // make sure jets do not overlap with an electron within DR<0.2
74 Jets jets;
75 for (const Jet& j : apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 25*GeV && Cuts::absrap < 4.4)) {
76 if (elecs.size() == 2 && (deltaR(j, *l1, RAPIDITY) < 0.2 || deltaR(j, *l2, RAPIDITY) < 0.2 )) {
77 continue;
78 }
79 jets += j;
80 }
81
82 // Require 2 jets with pT > 85 and 80 GeV
83 if (jets.size() < 2) vetoEvent;
84
85 // Calculate the observables
86 Variables vars(jets, l1, l2);
87
88 // make sure neither lepton overlaps with a jet within 0.4
89 for (const Jet& j : jets) {
90 if (deltaR(j, *l1, RAPIDITY) < 0.4 || deltaR(j, *l2, RAPIDITY) < 0.4) vetoEvent;
91 }
92
93 if (jets[0].pt() < 85*GeV || jets[1].pt() < 80*GeV ) vetoEvent;
94
95 bool pass_VBFZtopo = (vars.mjj > 250*GeV && vars.Dyjj > 2.0 && vars.pTll > 20*GeV && vars.Zcent < 1.0 && vars.pTbal < 0.15);
96
97 if (pass_VBFZtopo) {
98 if (_doControl && vars.Ngj > 0 && vars.Zcent < 0.5) fillPlots(vars, CRAplots);
99 else if (_doControl && vars.Ngj > 0 && vars.Zcent >= 0.5) fillPlots(vars, CRBplots);
100 else if (_doControl && vars.Ngj == 0 && vars.Zcent >= 0.5) fillPlots(vars, CRCplots);
101
102 if ( vars.Ngj == 0 && vars.Zcent < 0.5 ) {
103 fillPlots(vars, SRplots);
104 }
105 }
106 }
107
108
109 void finalize() {
110 const double xsec = crossSectionPerEvent()/femtobarn;
111 scalePlots(SRplots, xsec);
112 scalePlots(CRAplots, xsec);
113 scalePlots(CRBplots, xsec);
114 scalePlots(CRCplots, xsec);
115 }
116
117 /// @}
118
119
120 /// @name Analysis helpers
121 /// @{
122
123 struct Variables {
124 Variables(const vector<Jet>& jets, const Particle* l1, const Particle* l2) {
125 // get the jets
126 assert(jets.size()>=2);
127 FourMomentum j1 = jets[0].mom(), j2 = jets[1].mom();
128 pTj1 = j1.pT()/GeV; pTj2 = j2.pT()/GeV;
129 assert(pTj1 >= pTj2);
130
131 // build dilepton system
132 FourMomentum ll = (l1->mom() + l2->mom());
133 pTll = ll.pT(); mll = ll.mass();
134
135 Nj = jets.size();
136 Dyjj = std::abs(j1.rap() - j2.rap());
137 mjj = (j1 + j2).mass();
138 Dphijj = ( j1.rap() > j2.rap() ) ? mapAngleMPiToPi(j1.phi() - j2.phi()) : mapAngleMPiToPi(j2.phi() - j1.phi());
139
140 Jets gjets = getGapJets(jets);
141 Ngj = gjets.size();
142 pTgj = Ngj? gjets[0].pT()/GeV : 0;
143
144 FourMomentum vecSum = (j1 + j2 + l1->mom() + l2->mom());
145 double HT = j1.pT() + j2.pT() + l1->pT() + l2->pT();
146 if (Ngj) {
147 vecSum += gjets[0].mom();
148 HT += pTgj;
149 }
150 pTbal = vecSum.pT() / HT;
151
152 Zcent = std::abs(ll.rap() - (j1.rap() + j2.rap())/2) / Dyjj;
153 }
154
155 double Zcent, pTj1, pTj2, pTgj, pTll, mll, Dyjj, mjj, Dphijj, pTbal;
156 size_t Nj, Ngj;
157
158 Jets getGapJets(const Jets& jets) {
159 Jets gjets;
160 if (jets.size() <= 2) return gjets;
161 FourMomentum j1 = jets[0].mom(), j2 = jets[1].mom();
162 double yFwd = j1.rap(), yBwd = j2.rap();
163 if (yFwd < yBwd) std::swap(yFwd,yBwd);
164 for (size_t i = 2; i < jets.size(); ++i)
165 if (inRange(jets[i].rap(), yBwd, yFwd)) gjets += jets[i];
166 return gjets;
167 }
168
169 }; // struct variables
170
171
172 struct Plots {
173 string label;
174 Histo1DPtr m_jj, Dphi_jj, Dy_jj, pT_ll;
175 };
176
177
178 void initialisePlots(Plots& plots, const string& phase_space) {
179 plots.label = phase_space;
180 size_t region = 0;
181 if (phase_space == "SR") region = 4;
182 if (phase_space == "CRA") region = 8;
183 if (phase_space == "CRB") region = 12;
184 if (phase_space == "CRC") region = 16;
185 book(plots.m_jj, 1 + region, 1, 1);
186 book(plots.Dy_jj, 2 + region, 1, 1);
187 book(plots.pT_ll, 3 + region, 1, 1);
188 book(plots.Dphi_jj, 4 + region, 1, 1);
189 }
190
191
192 void fillPlots(const Variables& vars, Plots& plots) {
193 // The mjj plot extends down to 250 GeV
194 plots.m_jj->fill(vars.mjj/GeV);
195 if (vars.mjj > 1000*GeV) {
196 plots.Dy_jj->fill(vars.Dyjj);
197 plots.Dphi_jj->fill(vars.Dphijj);
198 plots.pT_ll->fill(vars.pTll/GeV);
199 }
200 }
201
202
203 void scalePlots(Plots& plots, const double xsec) {
204 scale(plots.m_jj, xsec);
205 scale(plots.Dy_jj, xsec);
206 scale(plots.Dphi_jj, xsec);
207 scale(plots.pT_ll, xsec);
208 }
209
210 /// @}
211
212
213 private:
214
215 Plots SRplots, CRAplots, CRBplots, CRCplots;
216 bool _doControl;
217
218 };
219
220
221
222 RIVET_DECLARE_PLUGIN(ATLAS_2020_I1803608);
223
224}
|