ATLAS Offline Software
Loading...
Searching...
No Matches
HiggsTemplateCrossSections.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_H
6#define TRUTHRIVETTOOLS_HIGGSTEMPLATECROSSSECTIONS_H
7
8// -*- C++ -*-
9#include "Rivet/Analysis.hh"
10#include "Rivet/Particle.hh"
11#include "Rivet/Projections/FastJets.hh"
12
13// Definition of the StatusCode and Category enums
14// Note: the Template XSec Defs *depends* on having included
15// the TLorentzVector header *before* it is included -- it
16// uses the include guard from TLorentzVector to decide
17// what is available
18#include "TLorentzVector.h"
21#include "AtlasHepMC/GenEvent.h"
26#include <string>
27#include <string_view>
28
29
30namespace Rivet {
31
37 public:
38 // Constructor
40 : Analysis("HiggsTemplateCrossSections"),
41 m_HiggsProdMode(HTXS::UNKNOWN) {}
42
43 public:
44
49
51 Particle getLastInstance(const Particle & ptcl) const {
52 if ( ptcl.genParticle()->end_vertex() ) {
53 if ( !hasChild(ptcl.genParticle(),ptcl.pid()) ) return ptcl;
54 else return getLastInstance(ptcl.children()[0]);
55 }
56 return ptcl;
57 }
58
60 bool originateFrom(const Particle& p, const Particles& ptcls ) const {
61 auto prodVtx = p.genParticle()->production_vertex();
62 if (prodVtx == nullptr) return false;
63 // for each ancestor, check if it matches any of the input particles
64 for (auto ancestor:Rivet::HepMCUtils::particles(std::move(prodVtx),Relatives::ANCESTORS)){
65 for ( const auto & part:ptcls )
66 if ( ancestor==part.genParticle() ) return true;
67 }
68 // if we get here, no ancestor matched any input particle
69 return false;
70 }
71
73 bool originateFrom(const Particle& p, const Particle& p2 ) const {
74 Particles ptcls = {p2}; return originateFrom(p,ptcls);
75 }
76
78 bool hasChild(HepMC::ConstGenParticlePtr ptcl, int pdgID) const {
79 for (const Particle& child:Particle(*ptcl).children())
80 if (child.pid()==pdgID) return true;
81 return false;
82 }
83
85 bool hasParent(HepMC::ConstGenParticlePtr ptcl, int pdgID) const {
86 for (auto parent:Rivet::HepMCUtils::particles(ptcl->production_vertex(),Relatives::PARENTS))
87 if (parent->pdg_id()==pdgID) return true;
88 return false;
89 }
90
92 bool quarkDecay(const Particle &p) const {
93 for (const Particle& child:p.children())
94 if (PID::isQuark(child.pid())) return true;
95 return false;
96 }
97
99 bool ChLeptonDecay(const Particle &p) const {
100 for (const Particle& child:p.children())
101 if (
102#if RIVET_VERSION_CODE >= 40000
103 PID::isChargedLepton(child.pid())
104#else
105 PID::isChLepton(child.pid())
106#endif // RIVET_VERSION_CODE
107 ) return true;
108 return false;
109 }
110
113 HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err,
114 std::string_view msg={}, int NmaxWarnings=20) const {
115 // Set the error, and keep statistics
116 cat.errorCode = err;
117 const auto errIndex = static_cast<std::size_t>(err);
118 if (errIndex < std::size(m_errorCount)) {
119 ++m_errorCount[errIndex];
120 } else {
121 MSG_WARNING("Invalid HTXS error code: " << errIndex);
122 }
123 // Print warning message to the screen/log
124 static std::atomic<int> Nwarnings = 0;
125 if ( !msg.empty() && ++Nwarnings < NmaxWarnings )
127
128 return cat;
129 }
130
131
133 HiggsClassification classifyEvent(const Event& event, const HTXS::HiggsProdMode prodMode ) const {
134
135 // the classification object
136 HiggsClassification cat;
137 cat.prodMode = prodMode;
138 cat.errorCode = HTXS::UNDEFINED;
139 cat.stage0_cat = HTXS::Stage0::UNKNOWN;
140 cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
141 cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
142 cat.stage1_2_cat_pTjet25GeV = HTXS::Stage1_2::UNKNOWN;
143 cat.stage1_2_cat_pTjet30GeV = HTXS::Stage1_2::UNKNOWN;
144 cat.stage1_2_fine_cat_pTjet25GeV = HTXS::Stage1_2_Fine::UNKNOWN;
145 cat.stage1_2_fine_cat_pTjet30GeV = HTXS::Stage1_2_Fine::UNKNOWN;
146
147 if (prodMode == HTXS::UNKNOWN)
148 return error(cat,HTXS::PRODMODE_DEFINED,
149 "Unkown Higgs production mechanism. Cannot classify event."
150 " Classification for all events will most likely fail.");
151
152 /*****
153 * Step 1.
154 * Idenfify the Higgs boson and the hard scatter vertex
155 * There should be only one of each.
156 */
157
158 auto HSvtx = HepMC::signal_process_vertex(event.genEvent());
159 int Nhiggs=0;
160 for (auto ptcl : Rivet::HepMCUtils::particles(event.genEvent()) ) {
161
162 // a) Reject all non-Higgs particles
163 if ( !PID::isHiggs(ptcl->pdg_id()) ) continue;
164 // b) select only the final Higgs boson copy, prior to decay
165 if ( ptcl->end_vertex() && !hasChild(ptcl,PID::HIGGS) ) {
166 cat.higgs = Particle(ptcl); ++Nhiggs;
167 }
168 // c) if HepMC::signal_proces_vertex is missing
169 // set hard-scatter vertex based on first Higgs boson
170 if ( HSvtx==nullptr && ptcl->production_vertex() && !hasParent(ptcl,PID::HIGGS) )
171 HSvtx = ptcl->production_vertex();
172 }
173
174 // Make sure things are in order so far
175 if (Nhiggs!=1)
177 "Current event has "+std::to_string(Nhiggs)+" Higgs bosons. There must be only one.");
178 if (cat.higgs.children().size()<2)
180 "Could not identify Higgs boson decay products.");
181
182 if (HSvtx == nullptr)
183 return error(cat,HTXS::HS_VTX_IDENTIFICATION,"Cannot find hard-scatter vertex of current event.");
184
185 /*****
186 * Step 2.
187 * Identify associated vector bosons
188 */
189
190 // Find associated vector bosons
191 bool is_uncatdV = false;
192 Particles uncatV_decays;
193 FourMomentum uncatV_p4(0,0,0,0);
194 FourVector uncatV_v4(0,0,0,0);
195 int nWs=0, nZs=0;
196 if ( isVH(prodMode) ) {
197 for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
198 if (PID::isW(ptcl->pdg_id())) { ++nWs; cat.V=Particle(ptcl); }
199 if (PID::isZ(ptcl->pdg_id())) { ++nZs; cat.V=Particle(std::move(ptcl)); }
200 }
201 if(nWs+nZs>0) cat.V = getLastInstance(cat.V);
202 else {
203 for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
204 if (!PID::isHiggs(ptcl->pdg_id())) {
205 uncatV_decays += Particle(ptcl);
206 uncatV_p4 += Particle(ptcl).momentum();
207 uncatV_v4 += Particle(std::move(ptcl)).origin();
208 }
209 }
210 is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
211 }
212 }
213
214 if ( !is_uncatdV ){
215
216 if ( isVH(prodMode) && !cat.V.genParticle()->end_vertex() )
217 return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
218
219 if ( isVH(prodMode) && cat.V.children().size()<2 )
220 return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
221
222 if ( ( prodMode==HTXS::WH && (nZs>0||nWs!=1) ) ||
223 ( (prodMode==HTXS::QQ2ZH||prodMode==HTXS::GG2ZH) && (nZs!=1||nWs>0) ) )
224 return error(cat,HTXS::VH_IDENTIFICATION,"Found "+std::to_string(nWs)+" W-bosons and "+
225 std::to_string(nZs)+" Z-bosons. Inconsitent with VH expectation.");
226 }
227
228 // Find and store the W-bosons from ttH->WbWbH
229 Particles Ws;
230 if ( prodMode==HTXS::TTH || prodMode==HTXS::TH ){
231 // loop over particles produced in hard-scatter vertex
232 for ( auto ptcl : Rivet::HepMCUtils::particles(std::move(HSvtx),Relatives::CHILDREN) ) {
233 if ( !PID::isTop(ptcl->pdg_id()) ) continue;
234 Particle top = getLastInstance(Particle(std::move(ptcl)));
235 if ( top.genParticle()->end_vertex() )
236 for (const auto &child:top.children())
237 if ( PID::isW(child.pid()) ) Ws += getLastInstance(child);
238 }
239 }
240
241 // Make sure result make sense
242 if ( (prodMode==HTXS::TTH && Ws.size()<2) || (prodMode==HTXS::TH && Ws.size()<1 ) )
243 return error(cat,HTXS::TOP_W_IDENTIFICATION,"Failed to identify W-boson(s) from t-decay!");
244
245 /*****
246 * Step 3.
247 * Build jets
248 * Make sure all stable particles are present
249 */
250
251 // Create a list of the vector bosons that decay leptonically
252 // Either the vector boson produced in association with the Higgs boson,
253 // or the ones produced from decays of top quarks produced with the Higgs
254 Particles leptonicVs;
255 if ( !is_uncatdV ){
256 if ( isVH(prodMode) && !quarkDecay(cat.V) ) leptonicVs += cat.V;
257 }else leptonicVs = std::move(uncatV_decays);
258 for ( const auto & W:Ws ) if ( W.genParticle()->end_vertex() && !quarkDecay(W) ) leptonicVs += W;
259
260 // Obtain all stable, final-state particles
261 const Particles FS = apply<FinalState>(event, "FS").particles();
262 Particles hadrons;
263 FourMomentum sum(0,0,0,0), vSum(0,0,0,0), hSum(0,0,0,0);
264 for ( const Particle &p : FS ) {
265 // Add up the four momenta of all stable particles as a cross check
266 sum += p.momentum();
267 // ignore particles from the Higgs boson
268 if ( originateFrom(p,cat.higgs) ) { hSum += p.momentum(); continue; }
269 // Cross-check the V decay products for VH
270 if ( isVH(prodMode) && !is_uncatdV && originateFrom(p,Ws) ) vSum += p.momentum();
271 // ignore final state particles from leptonic V decays
272 if ( leptonicVs.size() && originateFrom(p,leptonicVs) ) continue;
273 // All particles reaching here are considered hadrons and will be used to build jets
274 hadrons += p;
275 }
276
277 cat.p4decay_higgs = hSum;
278 cat.p4decay_V = vSum;
279
280 FinalState fps_temp;
281 FastJets jets(fps_temp,
282#if RIVET_VERSION_CODE >= 40000
283 JetAlg::ANTIKT,
284#else
285 FastJets::ANTIKT,
286#endif // RIVET_VERSION_CODE
287 0.4 );
288 jets.calc(hadrons);
289
290 cat.jets25 = jets.jetsByPt( Cuts::pT > 25.0 );
291 cat.jets30 = jets.jetsByPt( Cuts::pT > 30.0 );
292
293 // check that four mometum sum of all stable particles satisfies momentum consevation
294 if ( sum.pt()>0.1 )
295 return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
296 std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
297
298 // check if V-boson was not included in the event record but decay particles were
299 // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
300 if(is_uncatdV)
301 return error(cat,HTXS::VH_IDENTIFICATION,"Failed to identify associated V-boson!");
302
303 /*****
304 * Step 4.
305 * Classify and save output
306 */
307
308 // Apply the categorization categorization
309 cat.isZ2vvDecay = false;
310 if( (prodMode==HTXS::GG2ZH || prodMode==HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V) ) cat.isZ2vvDecay = true;
311 cat.stage0_cat = getStage0Category(prodMode,cat.higgs,cat.V);
312 cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode,cat.higgs,cat.jets25,cat.V);
313 cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode,cat.higgs,cat.jets30,cat.V);
314 cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets25,cat.V);
315 cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets30,cat.V);
316 cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets25,cat.V);
317 cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets30,cat.V);
318 cat.errorCode = HTXS::SUCCESS; ++m_errorCount[HTXS::SUCCESS];
319
320 return cat;
321 }
322
328
330int getBin(double x, const std::vector<double>& bins) const {
331 if (bins.empty() || x < bins.front()) {
332 throw std::invalid_argument("Input value is out of bin range or bins vector is empty.");
333 }
334
335 for (size_t i = 1; i < bins.size(); ++i) {
336 if (x < bins[i]) {
337 return static_cast<int>(i - 1);
338 }
339 }
340
341 return static_cast<int>(bins.size() - 1);
342}
343
347 int vbfTopology(const Jets &jets, const Particle &higgs) const {
348 if (jets.size()<2) return 0;
349 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
350 bool VBFtopo = (j1+j2).mass() > 400.0 && std::abs(j1.rapidity()-j2.rapidity()) > 2.8;
351 return VBFtopo ? (j1+j2+higgs.momentum()).pt()<25 ? 2 : 1 : 0;
352 }
353
357 int vbfTopology_Stage1_2(const Jets &jets, const Particle &higgs) const {
358 if (jets.size()<2) return 0;
359 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
360 double mjj = (j1+j2).mass();
361 if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
362 else if(mjj>700) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
363 else return 0;
364 }
365
371 int vbfTopology_Stage1_2_Fine(const Jets &jets, const Particle &higgs) const {
372 if (jets.size()<2) return 0;
373 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
374 double mjj = (j1+j2).mass();
375 if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
376 else if(mjj>700 && mjj<=1000) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
377 else if(mjj>1000 && mjj<=1500) return (j1+j2+higgs.momentum()).pt()<25 ? 5 : 6;
378 else if(mjj>1500) return (j1+j2+higgs.momentum()).pt()<25 ? 7 : 8;
379 else return 0;
380 }
381
383 bool isVH(HTXS::HiggsProdMode p) const { return p==HTXS::WH || p==HTXS::QQ2ZH || p==HTXS::GG2ZH; }
384
387 const Particle &higgs,
388 const Particle &V) const {
389 using namespace HTXS::Stage0;
390 int ctrlHiggs = std::abs(higgs.rapidity())<2.5;
391 // Special cases first, qq→Hqq
392 if ( (prodMode==HTXS::WH||prodMode==HTXS::QQ2ZH) && quarkDecay(V) ) {
393 return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
394 } else if ( prodMode==HTXS::GG2ZH && quarkDecay(V) ) {
395 return Category(HTXS::GGF*10 + ctrlHiggs);
396 }
397 // General case after
398 return Category(prodMode*10 + ctrlHiggs);
399 }
400
403 const Particle &higgs,
404 const Jets &jets,
405 const Particle &V) const {
406 using namespace HTXS::Stage1;
407 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
408 double pTj1 = jets.size() ? jets[0].momentum().pt() : 0;
409 int vbfTopo = vbfTopology(jets,higgs);
410
411 // 1. GGF Stage 1 categories
412 // Following YR4 write-up: XXXXX
413 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
414 if (fwdHiggs) return GG2H_FWDH;
415 if (Njets==0) return GG2H_0J;
416 else if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
417 else if (Njets>=2) {
418 // events with pT_H>200 get priority over VBF cuts
419 if(higgs.pt()<=200){
420 if (vbfTopo==2) return GG2H_VBFTOPO_JET3VETO;
421 else if (vbfTopo==1) return GG2H_VBFTOPO_JET3;
422 }
423 // Njets >= 2jets without VBF topology
424 return Category(GG2H_GE2J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
425 }
426 }
427 // 2. Electroweak qq->Hqq Stage 1 categories
428 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
429 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
430 if (pTj1>200) return QQ2HQQ_PTJET1_GT200;
431 if (vbfTopo==2) return QQ2HQQ_VBFTOPO_JET3VETO;
432 if (vbfTopo==1) return QQ2HQQ_VBFTOPO_JET3;
433 double mjj = jets.size()>1 ? (jets[0].mom()+jets[1].mom()).mass():0;
434 if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_VH2JET;
435 return QQ2HQQ_REST;
436 }
437 // 3. WH->Hlv categories
438 else if (prodMode==HTXS::WH) {
439 if (fwdHiggs) return QQ2HLNU_FWDH;
440 else if (V.pt()<150) return QQ2HLNU_PTV_0_150;
441 else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
442 // 150 < pTV/GeV < 250
443 return jets.size()==0 ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
444 }
445 // 4. qq->ZH->llH categories
446 else if (prodMode==HTXS::QQ2ZH) {
447 if (fwdHiggs) return QQ2HLL_FWDH;
448 else if (V.pt()<150) return QQ2HLL_PTV_0_150;
449 else if (V.pt()>250) return QQ2HLL_PTV_GT250;
450 // 150 < pTV/GeV < 250
451 return jets.size()==0 ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
452 }
453 // 5. gg->ZH->llH categories
454 else if (prodMode==HTXS::GG2ZH ) {
455 if (fwdHiggs) return GG2HLL_FWDH;
456 if (V.pt()<150) return GG2HLL_PTV_0_150;
457 else if (jets.size()==0) return GG2HLL_PTV_GT150_0J;
459 }
460 // 6.ttH,bbH,tH categories
461 else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
462 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
463 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
464 return UNKNOWN;
465 }
466
469 const Particle &higgs,
470 const Jets &jets,
471 const Particle &V) const {
472 using namespace HTXS::Stage1_2;
473 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
474 int vbfTopo = vbfTopology_Stage1_2(jets,higgs);
475
476 // 1. GGF Stage 1 categories
477 // Following YR4 write-up: XXXXX
478 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
479 if (fwdHiggs) return GG2H_FWDH;
480 if ( higgs.pt()>200 ) return Category(GG2H_PTH_200_300+getBin(higgs.pt(),{200,300,450,650}));
481 if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
482 if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
483 if (Njets>1){
484 //VBF topology
485 if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
486 //Njets >= 2jets without VBF topology (mjj<350)
487 return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
488 }
489 }
490
491 // 2. Electroweak qq->Hqq Stage 1.2 categories
492 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
493 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
494 int Njets=jets.size();
495 if (Njets==0) return QQ2HQQ_0J;
496 else if (Njets==1) return QQ2HQQ_1J;
497 else if (Njets>=2) {
498 double mjj = (jets[0].mom()+jets[1].mom()).mass();
499 if ( mjj < 60 ) return QQ2HQQ_GE2J_MJJ_0_60;
500 else if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_GE2J_MJJ_60_120;
501 else if ( 120 < mjj && mjj < 350 ) return QQ2HQQ_GE2J_MJJ_120_350;
502 else if ( mjj > 350 ) {
503 if (higgs.pt()>200) return QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200;
504 if(vbfTopo) return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200+vbfTopo);
505 }
506 }
507 }
508 // 3. WH->Hlv categories
509 else if (prodMode==HTXS::WH) {
510 if (fwdHiggs) return QQ2HLNU_FWDH;
511 else if (V.pt()<75) return QQ2HLNU_PTV_0_75;
512 else if (V.pt()<150) return QQ2HLNU_PTV_75_150;
513 else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
514 // 150 < pTV/GeV < 250
515 return jets.size()==0 ? QQ2HLNU_PTV_150_250_0J : QQ2HLNU_PTV_150_250_GE1J;
516 }
517 // 4. qq->ZH->llH categories
518 else if (prodMode==HTXS::QQ2ZH) {
519 if (fwdHiggs) return QQ2HLL_FWDH;
520 else if (V.pt()<75) return QQ2HLL_PTV_0_75;
521 else if (V.pt()<150) return QQ2HLL_PTV_75_150;
522 else if (V.pt()>250) return QQ2HLL_PTV_GT250;
523 // 150 < pTV/GeV < 250
524 return jets.size()==0 ? QQ2HLL_PTV_150_250_0J : QQ2HLL_PTV_150_250_GE1J;
525 }
526 // 5. gg->ZH->llH categories
527 else if (prodMode==HTXS::GG2ZH ) {
528 if (fwdHiggs) return GG2HLL_FWDH;
529 else if (V.pt()<75) return GG2HLL_PTV_0_75;
530 else if (V.pt()<150) return GG2HLL_PTV_75_150;
531 else if (V.pt()>250) return GG2HLL_PTV_GT250;
532 return jets.size()==0 ? GG2HLL_PTV_150_250_0J : GG2HLL_PTV_150_250_GE1J;
533 }
534 // 6.ttH,bbH,tH categories
535 else if (prodMode==HTXS::TTH) {
536 if (fwdHiggs) return TTH_FWDH;
537 else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300}));
538 }
539 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
540 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
541 return UNKNOWN;
542 }
543
546 const Particle &higgs,
547 const Jets &jets,
548 const Particle &V) const {
549 using namespace HTXS::Stage1_2_Fine;
550 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
551 int vbfTopo = vbfTopology_Stage1_2_Fine(jets,higgs);
552
553 // 1. GGF Stage 1.2 categories
554 // Following YR4 write-up: XXXXX
555 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
556 if (fwdHiggs) return GG2H_FWDH;
557 if ( higgs.pt()>200 ){
558 if (Njets>0){
559 double pTHj = (jets[0].momentum()+higgs.momentum()).pt();
560 if( pTHj/higgs.pt()>0.15 ) return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15+getBin(higgs.pt(),{200,300,450,650}));
561 else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
562 }
563 else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
564 }
565 if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
566 if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
567 if (Njets>1){
568 //double mjj = (jets[0].mom()+jets[1].mom()).mass();
569 double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
570 //VBF topology
571 if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
572 //Njets >= 2jets without VBF topology (mjj<350)
573 if (pTHjj<25) return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25+getBin(higgs.pt(),{0,60,120,200}));
574 else return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25+getBin(higgs.pt(),{0,60,120,200}));
575 }
576 }
577
578 // 2. Electroweak qq->Hqq Stage 1.2 categories
579 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
580 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
581 int Njets=jets.size();
582 if (Njets==0) return QQ2HQQ_0J;
583 else if (Njets==1) return QQ2HQQ_1J;
584 else if (Njets>=2) {
585 double mjj = (jets[0].mom()+jets[1].mom()).mass();
586 double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
587 if (mjj<350){
588 if (pTHjj<25) return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25+getBin(mjj,{0,60,120,350}));
589 else return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25+getBin(mjj,{0,60,120,350}));
590 } else { //mjj>350 GeV
591 if (higgs.pt()<200){
593 } else {
595 }
596 }
597 }
598 }
599 // 3. WH->Hlv categories
600 else if (prodMode==HTXS::WH) {
601 if (fwdHiggs) return QQ2HLNU_FWDH;
602 int Njets=jets.size();
603 if (Njets==0) return Category(QQ2HLNU_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
604 if (Njets==1) return Category(QQ2HLNU_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
605 return Category(QQ2HLNU_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
606 }
607 // 4. qq->ZH->llH categories
608 else if (prodMode==HTXS::QQ2ZH) {
609 if (fwdHiggs) return QQ2HLL_FWDH;
610 int Njets=jets.size();
611 if (Njets==0) return Category(QQ2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
612 if (Njets==1) return Category(QQ2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
613 return Category(QQ2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
614 }
615 // 5. gg->ZH->llH categories
616 else if (prodMode==HTXS::GG2ZH ) {
617 if (fwdHiggs) return GG2HLL_FWDH;
618 int Njets=jets.size();
619 if (Njets==0) return Category(GG2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
620 if (Njets==1) return Category(GG2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
621 return Category(GG2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
622 }
623 // 6.ttH,bbH,tH categories
624 else if (prodMode==HTXS::TTH) {
625 if (fwdHiggs) return TTH_FWDH;
626 else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300,450}));
627 }
628 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
629 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
630 return UNKNOWN;
631 }
632
634
635
638
641
645 void init() {
646 printf("==============================================================\n");
647 printf("======== HiggsTemplateCrossSections Initialization =========\n");
648 printf("==============================================================\n");
649 // check that the production mode has been set
650 // if running in standalone Rivet the production mode is set through an env variable
652 char *pm_env = getenv("HIGGSPRODMODE");
653 string pm(pm_env==nullptr?"":pm_env);
654 if ( pm == "GGF" ) m_HiggsProdMode = HTXS::GGF;
655 else if ( pm == "VBF" ) m_HiggsProdMode = HTXS::VBF;
656 else if ( pm == "WH" ) m_HiggsProdMode = HTXS::WH;
657 else if ( pm == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
658 else if ( pm == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
659 else if ( pm == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
660 else if ( pm == "TTH" ) m_HiggsProdMode = HTXS::TTH;
661 else if ( pm == "BBH" ) m_HiggsProdMode = HTXS::BBH;
662 else if ( pm == "TH" ) m_HiggsProdMode = HTXS::TH;
663 else {
664 MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
665 }
666 }
667
668 // Projections for final state particles
669 const FinalState FS;
670 declare(FS,"FS");
671
672 // initialize the histograms with for each of the stages
674 m_sumw = 0.0;
675 printf("==============================================================\n");
676 printf("======== Higgs prod mode %d =========\n",m_HiggsProdMode);
677 printf("======== Sucessful Initialization =========\n");
678 printf("==============================================================\n");
679 }
680
681 // Perform the per-event analysis
682 void analyze(const Event& event) {
683
684 // get the classification
685 HiggsClassification cat = classifyEvent(event,m_HiggsProdMode);
686
687 // Fill histograms: categorization --> linerize the categories
688 const double weight = 1.; // Event weights are now all 1 in Rivet
689 m_sumw += weight;
690
691 int F=cat.stage0_cat%10, P=cat.stage1_cat_pTjet30GeV/100;
692 m_hist_stage0->fill( cat.stage0_cat/10*2+F, weight );
693
694 // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
695 static const vector<int> offset({0,1,13,19,24,29,33,35,37,39});
696 int off = offset[P];
697 // Stage 1.2 enum offsets for each production mode: GGF=17, VBF=11, WH= 6, QQ2ZH=6, GG2ZH=6, TTH=6, BBH=2, TH=2
698 static const vector<int> offset1_2({0,1,18,29,35,41,47,53,55,57});
699 int off1_2 = offset1_2[P];
700 // Stage 1.2-Fine enum offsets for each production mode: GGF=28, VBF=25, WH= 16, QQ2ZH=16, GG2ZH=16, TTH=7, BBH=2, TH=2
701 static const vector<int> offset1_2_Fine({0,1,29,54,70,86,102,109,111,113});
702 int off1_2_Fine = offset1_2_Fine[P];
703
704 m_hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV%100 + off, weight);
705 m_hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV%100 + off, weight);
706 m_hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV%100 + off1_2, weight);
707 m_hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV%100 + off1_2, weight);
708 m_hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV%100 + off1_2_Fine, weight);
709 m_hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV%100 + off1_2_Fine, weight);
710
711 // Fill histograms: variables used in the categorization
712 m_hist_pT_Higgs->fill(cat.higgs.pT(),weight);
713 m_hist_y_Higgs->fill(cat.higgs.rapidity(),weight);
714 m_hist_pT_V->fill(cat.V.pT(),weight);
715
716 m_hist_Njets25->fill(cat.jets25.size(),weight);
717 m_hist_Njets30->fill(cat.jets30.size(),weight);
718
719 m_hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
720
721 // Jet variables. Use jet collection with pT threshold at 30 GeV
722 if (cat.jets30.size()) m_hist_pT_jet1->fill(cat.jets30[0].pt(),weight);
723 if (cat.jets30.size()>=2) {
724 const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
725 m_hist_deltay_jj->fill(std::abs(j1.rapidity()-j2.rapidity()),weight);
726 m_hist_dijet_mass->fill((j1+j2).mass(),weight);
727 m_hist_pT_Hjj->fill((j1+j2+cat.higgs.momentum()).pt(),weight);
728 }
729 }
730
732 MSG_INFO (" ====================================================== ");
733 MSG_INFO (" Higgs Template X-Sec Categorization Tool ");
734 MSG_INFO (" Status Code Summary ");
735 MSG_INFO (" ====================================================== ");
736 bool allSuccess = (numEvents()==m_errorCount[HTXS::SUCCESS]);
737 if ( allSuccess ) MSG_INFO (" >>>> All "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized!");
738 else{
739 MSG_INFO (" >>>> "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized");
740 MSG_INFO (" >>>> --> the following errors occured:");
741 MSG_INFO (" >>>> "<< m_errorCount[HTXS::PRODMODE_DEFINED] <<" had an undefined Higgs production mode.");
742 MSG_INFO (" >>>> "<< m_errorCount[HTXS::MOMENTUM_CONSERVATION] <<" failed momentum conservation.");
743 MSG_INFO (" >>>> "<< m_errorCount[HTXS::HIGGS_IDENTIFICATION] <<" failed to identify a valid Higgs boson.");
744 MSG_INFO (" >>>> "<< m_errorCount[HTXS::HS_VTX_IDENTIFICATION] <<" failed to identify the hard scatter vertex.");
745 MSG_INFO (" >>>> "<< m_errorCount[HTXS::VH_IDENTIFICATION] <<" VH: to identify a valid V-boson.");
746 MSG_INFO (" >>>> "<< m_errorCount[HTXS::TOP_W_IDENTIFICATION] <<" failed to identify valid Ws from top decay.");
747 }
748 MSG_INFO (" ====================================================== ");
749 MSG_INFO (" ====================================================== ");
750 }
751
752
760
761 /*
762 * initialize histograms
763 */
764
766 book(m_hist_stage0,"HTXS_stage0",20,0,20);
767 book(m_hist_stage1_pTjet25,"HTXS_stage1_pTjet25",40,0,40);
768 book(m_hist_stage1_pTjet30,"HTXS_stage1_pTjet30",40,0,40);
769 book(m_hist_stage1_2_pTjet25,"HTXS_stage1_2_pTjet25",57,0,57);
770 book(m_hist_stage1_2_pTjet30,"HTXS_stage1_2_pTjet30",57,0,57);
771 book(m_hist_stage1_2_fine_pTjet25,"HTXS_stage1_2_fine_pTjet25",113,0,113);
772 book(m_hist_stage1_2_fine_pTjet30,"HTXS_stage1_2_fine_pTjet30",113,0,113);
773 book(m_hist_pT_Higgs,"pT_Higgs",80,0,400);
774 book(m_hist_y_Higgs,"y_Higgs",80,-4,4);
775 book(m_hist_pT_V,"pT_V",80,0,400);
776 book(m_hist_pT_jet1,"pT_jet1",80,0,400);
777 book(m_hist_deltay_jj ,"deltay_jj",50,0,10);
778 book(m_hist_dijet_mass,"m_jj",50,0,2000);
779 book(m_hist_pT_Hjj,"pT_Hjj",50,0,250);
780 book(m_hist_Njets25,"Njets25",10,0,10);
781 book(m_hist_Njets30,"Njets30",10,0,10);
782 book(m_hist_isZ2vv,"isZ2vv",2,0,2);
783 }
784
785
786 /*
787 * initialize private members used in the classification procedure
788 */
789
790 private:
791 double m_sumw=0.0;
793 mutable std::array<std::atomic<size_t>, HTXS::NUM_ERRORCODES> m_errorCount ATLAS_THREAD_SAFE {};
794 Histo1DPtr m_hist_stage0;
802 Histo1DPtr m_hist_isZ2vv;
803 };
804
805 // the PLUGIN only needs to be decleared when running standalone Rivet
806 // and causes compilation / linking issues if included in Athena / RootCore
807 //check for Rivet environment variable RIVET_ANALYSIS_PATH
808#ifdef RIVET_ANALYSIS_PATH
809 // The hook for the plugin system
810 DECLARE_RIVET_PLUGIN(HiggsTemplateCrossSections);
811#endif
812
813}
814
815#endif
static Double_t P(Double_t *tt, Double_t *par)
#define F(x, y, z)
Definition MD5.cxx:112
static const std::vector< std::string > bins
@ top
#define x
Define macros for attributes used to control the static checker.
Rivet routine for classifying MC events according to the Higgs template cross section categories.
int vbfTopology_Stage1_2_Fine(const Jets &jets, const Particle &higgs) const
VBF topology selection for Stage1_2 0 = fail loose selection: m_jj > 350 GeV 1 pass loose,...
HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err, std::string_view msg={}, int NmaxWarnings=20) const
Returns the classification object with the error code set.
bool originateFrom(const Particle &p, const Particle &p2) const
Whether particle p originates from p2.
void setHiggsProdMode(HTXS::HiggsProdMode prodMode)
Sets the Higgs production mode.
void init()
default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Hig...
bool originateFrom(const Particle &p, const Particles &ptcls) const
Whether particle p originate from any of the ptcls.
int vbfTopology_Stage1_2(const Jets &jets, const Particle &higgs) const
VBF topology selection 0 = fail loose selection: m_jj > 350 GeV 1 pass loose, but fail additional cut...
bool isVH(HTXS::HiggsProdMode p) const
Whether the Higgs is produced in association with a vector boson (VH).
Particle getLastInstance(const Particle &ptcl) const
follow a "propagating" particle and return its last instance
HTXS::Stage1_2_Fine::Category getStage1_2_Fine_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V) const
Stage-1.2_Fine categorization.
bool quarkDecay(const Particle &p) const
Return true is particle decays to quarks.
bool hasParent(HepMC::ConstGenParticlePtr ptcl, int pdgID) const
Checks whether the input particle has a parent with a given PDGID.
int getBin(double x, const std::vector< double > &bins) const
Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode) const
Main classificaion method.
bool hasChild(HepMC::ConstGenParticlePtr ptcl, int pdgID) const
Checks whether the input particle has a child with a given PDGID.
HTXS::Stage1::Category getStage1Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V) const
Stage-1 categorization.
HTXS::Stage0::Category getStage0Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V) const
Stage-0 HTXS categorization.
bool ChLeptonDecay(const Particle &p) const
Return true if particle decays to charged leptons.
HTXS::Stage1_2::Category getStage1_2_Category(const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V) const
Stage-1.2 categorization.
std::array< std::atomic< size_t >, HTXS::NUM_ERRORCODES > m_errorCount ATLAS_THREAD_SAFE
int vbfTopology(const Jets &jets, const Particle &higgs) const
VBF topolog selection 0 = fail loose selction: m_jj > 400 GeV and Dy_jj > 2.8 1 pass loose,...
The namespace of all packages in PhysicsAnalysis/JetTagging.
Two digit number of format PF P is digit for the physics process and F is 0 for |yH|>2....
Categorization Stage 1.2: Three digit integer of format PF Where P is a digit representing the proces...
Categorization Stage 1: Three digit integer of format PF Where P is a digit representing the process ...
Higgs Template Cross Section namespace.
ErrorCode
Error code: whether the classification was successful or failed.
@ HS_VTX_IDENTIFICATION
failed to identify hard scatter vertex
@ PRODMODE_DEFINED
production mode not defined
@ SUCCESS
successful classification
@ VH_DECAY_IDENTIFICATION
failed to identify associated vector boson decay products
@ HIGGS_IDENTIFICATION
failed to identify Higgs boson
@ TOP_W_IDENTIFICATION
failed to identify top decay
@ HIGGS_DECAY_IDENTIFICATION
failed to identify Higgs boson decay products
@ MOMENTUM_CONSERVATION
failed momentum conservation
@ NUM_ERRORCODES
number of error codes (keep this unnumbered and last)
@ VH_IDENTIFICATION
failed to identify associated vector boson
HiggsProdMode
Higgs production modes, corresponding to input sample.
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
GenVertex * signal_process_vertex(const GenEvent *e)
Definition GenEvent.h:645
Definition Jets.py:1
#define MSG_WARNING(ARG)
#define MSG_INFO(ARG)
MsgStream & msg
Definition testRead.cxx:32