ATLAS Offline Software
Loading...
Searching...
No Matches
Rivet::HiggsTemplateCrossSections Class Reference

Rivet routine for classifying MC events according to the Higgs template cross section categories. More...

#include <HiggsTemplateCrossSections.h>

Inheritance diagram for Rivet::HiggsTemplateCrossSections:
Collaboration diagram for Rivet::HiggsTemplateCrossSections:

Public Member Functions

 HiggsTemplateCrossSections ()
HiggsClassification classifyEvent (const Event &event, const HTXS::HiggsProdMode prodMode) const
 Main classificaion method.
Utility methods

Methods to identify the Higgs boson and associated vector boson and to build jets

Particle getLastInstance (const Particle &ptcl) const
 follow a "propagating" particle and return its last instance
bool originateFrom (const Particle &p, const Particles &ptcls) const
 Whether particle p originate from any of the ptcls.
bool originateFrom (const Particle &p, const Particle &p2) const
 Whether particle p originates from p2.
bool hasChild (HepMC::ConstGenParticlePtr ptcl, int pdgID) const
 Checks whether the input particle has a child with a given PDGID.
bool hasParent (HepMC::ConstGenParticlePtr ptcl, int pdgID) const
 Checks whether the input particle has a parent with a given PDGID.
bool quarkDecay (const Particle &p) const
 Return true is particle decays to quarks.
bool ChLeptonDecay (const Particle &p) const
 Return true if particle decays to charged leptons.
HiggsClassification error (HiggsClassification &cat, HTXS::ErrorCode err, std::string msg="", int NmaxWarnings=20) const
 Returns the classification object with the error code set.
Categorization methods

Methods to assign the truth category based on the identified Higgs boson and associated vector bosons and/or reconstructed jets

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.
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, but fail additional cut pT(Hjj)<25.
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 pT(Hjj)<25.
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, but fail additional cut pT(Hjj)<25.
bool isVH (HTXS::HiggsProdMode p) const
 Whether the Higgs is produced in association with a vector boson (VH)
HTXS::Stage0::Category getStage0Category (const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Particle &V) const
 Stage-0 HTXS categorization.
HTXS::Stage1::Category getStage1Category (const HTXS::HiggsProdMode prodMode, const Particle &higgs, const Jets &jets, const Particle &V) const
 Stage-1 categorization.
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.
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.
Default Rivet analysis methods and steering methods
void setHiggsProdMode (HTXS::HiggsProdMode prodMode)
 Sets the Higgs production mode.
void init ()
 default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Higgs production mode from shell variable if not set manually using setHiggsProdMode
void analyze (const Event &event)
void printClassificationSummary ()
void finalize ()
void initializeHistos ()

Private Attributes

double m_sumw =0.0
HTXS::HiggsProdMode m_HiggsProdMode
std::array< std::atomic< size_t >, HTXS::NUM_ERRORCODES > m_errorCount ATLAS_THREAD_SAFE {}
Histo1DPtr m_hist_stage0
Histo1DPtr m_hist_stage1_pTjet25
Histo1DPtr m_hist_stage1_pTjet30
Histo1DPtr m_hist_stage1_2_pTjet25
Histo1DPtr m_hist_stage1_2_pTjet30
Histo1DPtr m_hist_stage1_2_fine_pTjet25
Histo1DPtr m_hist_stage1_2_fine_pTjet30
Histo1DPtr m_hist_pT_Higgs
Histo1DPtr m_hist_y_Higgs
Histo1DPtr m_hist_pT_V
Histo1DPtr m_hist_pT_jet1
Histo1DPtr m_hist_deltay_jj
Histo1DPtr m_hist_dijet_mass
Histo1DPtr m_hist_pT_Hjj
Histo1DPtr m_hist_Njets25
Histo1DPtr m_hist_Njets30
Histo1DPtr m_hist_isZ2vv

Detailed Description

Rivet routine for classifying MC events according to the Higgs template cross section categories.

Author
Jim Lacey (DESY) <james.nosp@m..lac.nosp@m.ey@ce.nosp@m.rn.c.nosp@m.h,jlace.nosp@m.y@de.nosp@m.sy.de>
Dag Gillberg (Carleton University) dag.g.nosp@m.illb.nosp@m.erg@c.nosp@m.ern..nosp@m.ch

Definition at line 33 of file HiggsTemplateCrossSections.h.

Constructor & Destructor Documentation

◆ HiggsTemplateCrossSections()

Rivet::HiggsTemplateCrossSections::HiggsTemplateCrossSections ( )
inline

Definition at line 36 of file HiggsTemplateCrossSections.h.

37 : Analysis("HiggsTemplateCrossSections"),

Member Function Documentation

◆ analyze()

void Rivet::HiggsTemplateCrossSections::analyze ( const Event & event)
inline

Definition at line 675 of file HiggsTemplateCrossSections.h.

675 {
676
677 // get the classification
678 HiggsClassification cat = classifyEvent(event,m_HiggsProdMode);
679
680 // Fill histograms: categorization --> linerize the categories
681 const double weight = 1.; // Event weights are now all 1 in Rivet
682 m_sumw += weight;
683
684 int F=cat.stage0_cat%10, P=cat.stage1_cat_pTjet30GeV/100;
685 m_hist_stage0->fill( cat.stage0_cat/10*2+F, weight );
686
687 // Stage 1 enum offsets for each production mode: GGF=12, VBF=6, WH= 5, QQ2ZH=5, GG2ZH=4, TTH=2, BBH=2, TH=2
688 static const vector<int> offset({0,1,13,19,24,29,33,35,37,39});
689 int off = offset[P];
690 // 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
691 static const vector<int> offset1_2({0,1,18,29,35,41,47,53,55,57});
692 int off1_2 = offset1_2[P];
693 // 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
694 static const vector<int> offset1_2_Fine({0,1,29,54,70,86,102,109,111,113});
695 int off1_2_Fine = offset1_2_Fine[P];
696
697 m_hist_stage1_pTjet25->fill(cat.stage1_cat_pTjet25GeV%100 + off, weight);
698 m_hist_stage1_pTjet30->fill(cat.stage1_cat_pTjet30GeV%100 + off, weight);
699 m_hist_stage1_2_pTjet25->fill(cat.stage1_2_cat_pTjet25GeV%100 + off1_2, weight);
700 m_hist_stage1_2_pTjet30->fill(cat.stage1_2_cat_pTjet30GeV%100 + off1_2, weight);
701 m_hist_stage1_2_fine_pTjet25->fill(cat.stage1_2_fine_cat_pTjet25GeV%100 + off1_2_Fine, weight);
702 m_hist_stage1_2_fine_pTjet30->fill(cat.stage1_2_fine_cat_pTjet30GeV%100 + off1_2_Fine, weight);
703
704 // Fill histograms: variables used in the categorization
705 m_hist_pT_Higgs->fill(cat.higgs.pT(),weight);
706 m_hist_y_Higgs->fill(cat.higgs.rapidity(),weight);
707 m_hist_pT_V->fill(cat.V.pT(),weight);
708
709 m_hist_Njets25->fill(cat.jets25.size(),weight);
710 m_hist_Njets30->fill(cat.jets30.size(),weight);
711
712 m_hist_isZ2vv->fill(cat.isZ2vvDecay, weight);
713
714 // Jet variables. Use jet collection with pT threshold at 30 GeV
715 if (cat.jets30.size()) m_hist_pT_jet1->fill(cat.jets30[0].pt(),weight);
716 if (cat.jets30.size()>=2) {
717 const FourMomentum &j1 = cat.jets30[0].momentum(), &j2 = cat.jets30[1].momentum();
718 m_hist_deltay_jj->fill(std::abs(j1.rapidity()-j2.rapidity()),weight);
719 m_hist_dijet_mass->fill((j1+j2).mass(),weight);
720 m_hist_pT_Hjj->fill((j1+j2+cat.higgs.momentum()).pt(),weight);
721 }
722 }
static Double_t P(Double_t *tt, Double_t *par)
#define F(x, y, z)
Definition MD5.cxx:112
HiggsClassification classifyEvent(const Event &event, const HTXS::HiggsProdMode prodMode) const
Main classificaion method.

◆ ChLeptonDecay()

bool Rivet::HiggsTemplateCrossSections::ChLeptonDecay ( const Particle & p) const
inline

Return true if particle decays to charged leptons.

Definition at line 96 of file HiggsTemplateCrossSections.h.

96 {
97 for (const Particle& child:p.children())
98 if (
99#if RIVET_VERSION_CODE >= 40000
100 PID::isChargedLepton(child.pid())
101#else
102 PID::isChLepton(child.pid())
103#endif // RIVET_VERSION_CODE
104 ) return true;
105 return false;
106 }

◆ classifyEvent()

HiggsClassification Rivet::HiggsTemplateCrossSections::classifyEvent ( const Event & event,
const HTXS::HiggsProdMode prodMode ) const
inline

Main classificaion method.

Definition at line 126 of file HiggsTemplateCrossSections.h.

126 {
127
128 // the classification object
129 HiggsClassification cat;
130 cat.prodMode = prodMode;
131 cat.errorCode = HTXS::UNDEFINED;
132 cat.stage0_cat = HTXS::Stage0::UNKNOWN;
133 cat.stage1_cat_pTjet25GeV = HTXS::Stage1::UNKNOWN;
134 cat.stage1_cat_pTjet30GeV = HTXS::Stage1::UNKNOWN;
135 cat.stage1_2_cat_pTjet25GeV = HTXS::Stage1_2::UNKNOWN;
136 cat.stage1_2_cat_pTjet30GeV = HTXS::Stage1_2::UNKNOWN;
137 cat.stage1_2_fine_cat_pTjet25GeV = HTXS::Stage1_2_Fine::UNKNOWN;
138 cat.stage1_2_fine_cat_pTjet30GeV = HTXS::Stage1_2_Fine::UNKNOWN;
139
140 if (prodMode == HTXS::UNKNOWN)
141 return error(cat,HTXS::PRODMODE_DEFINED,
142 "Unkown Higgs production mechanism. Cannot classify event."
143 " Classification for all events will most likely fail.");
144
145 /*****
146 * Step 1.
147 * Idenfify the Higgs boson and the hard scatter vertex
148 * There should be only one of each.
149 */
150
151 auto HSvtx = HepMC::signal_process_vertex(event.genEvent());
152 int Nhiggs=0;
153 for (auto ptcl : Rivet::HepMCUtils::particles(event.genEvent()) ) {
154
155 // a) Reject all non-Higgs particles
156 if ( !PID::isHiggs(ptcl->pdg_id()) ) continue;
157 // b) select only the final Higgs boson copy, prior to decay
158 if ( ptcl->end_vertex() && !hasChild(ptcl,PID::HIGGS) ) {
159 cat.higgs = Particle(ptcl); ++Nhiggs;
160 }
161 // c) if HepMC::signal_proces_vertex is missing
162 // set hard-scatter vertex based on first Higgs boson
163 if ( HSvtx==nullptr && ptcl->production_vertex() && !hasParent(ptcl,PID::HIGGS) )
164 HSvtx = ptcl->production_vertex();
165 }
166
167 // Make sure things are in order so far
168 if (Nhiggs!=1)
170 "Current event has "+std::to_string(Nhiggs)+" Higgs bosons. There must be only one.");
171 if (cat.higgs.children().size()<2)
173 "Could not identify Higgs boson decay products.");
174
175 if (HSvtx == nullptr)
176 return error(cat,HTXS::HS_VTX_IDENTIFICATION,"Cannot find hard-scatter vertex of current event.");
177
178 /*****
179 * Step 2.
180 * Identify associated vector bosons
181 */
182
183 // Find associated vector bosons
184 bool is_uncatdV = false;
185 Particles uncatV_decays;
186 FourMomentum uncatV_p4(0,0,0,0);
187 FourVector uncatV_v4(0,0,0,0);
188 int nWs=0, nZs=0;
189 if ( isVH(prodMode) ) {
190 for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
191 if (PID::isW(ptcl->pdg_id())) { ++nWs; cat.V=Particle(ptcl); }
192 if (PID::isZ(ptcl->pdg_id())) { ++nZs; cat.V=Particle(std::move(ptcl)); }
193 }
194 if(nWs+nZs>0) cat.V = getLastInstance(cat.V);
195 else {
196 for (auto ptcl:Rivet::HepMCUtils::particles(HSvtx,Relatives::CHILDREN)) {
197 if (!PID::isHiggs(ptcl->pdg_id())) {
198 uncatV_decays += Particle(ptcl);
199 uncatV_p4 += Particle(ptcl).momentum();
200 uncatV_v4 += Particle(std::move(ptcl)).origin();
201 }
202 }
203 is_uncatdV = true; cat.V = Particle(24,uncatV_p4,uncatV_v4);
204 }
205 }
206
207 if ( !is_uncatdV ){
208
209 if ( isVH(prodMode) && !cat.V.genParticle()->end_vertex() )
210 return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
211
212 if ( isVH(prodMode) && cat.V.children().size()<2 )
213 return error(cat,HTXS::VH_DECAY_IDENTIFICATION,"Vector boson does not decay!");
214
215 if ( ( prodMode==HTXS::WH && (nZs>0||nWs!=1) ) ||
216 ( (prodMode==HTXS::QQ2ZH||prodMode==HTXS::GG2ZH) && (nZs!=1||nWs>0) ) )
217 return error(cat,HTXS::VH_IDENTIFICATION,"Found "+std::to_string(nWs)+" W-bosons and "+
218 std::to_string(nZs)+" Z-bosons. Inconsitent with VH expectation.");
219 }
220
221 // Find and store the W-bosons from ttH->WbWbH
222 Particles Ws;
223 if ( prodMode==HTXS::TTH || prodMode==HTXS::TH ){
224 // loop over particles produced in hard-scatter vertex
225 for ( auto ptcl : Rivet::HepMCUtils::particles(std::move(HSvtx),Relatives::CHILDREN) ) {
226 if ( !PID::isTop(ptcl->pdg_id()) ) continue;
227 Particle top = getLastInstance(Particle(std::move(ptcl)));
228 if ( top.genParticle()->end_vertex() )
229 for (const auto &child:top.children())
230 if ( PID::isW(child.pid()) ) Ws += getLastInstance(child);
231 }
232 }
233
234 // Make sure result make sense
235 if ( (prodMode==HTXS::TTH && Ws.size()<2) || (prodMode==HTXS::TH && Ws.size()<1 ) )
236 return error(cat,HTXS::TOP_W_IDENTIFICATION,"Failed to identify W-boson(s) from t-decay!");
237
238 /*****
239 * Step 3.
240 * Build jets
241 * Make sure all stable particles are present
242 */
243
244 // Create a list of the vector bosons that decay leptonically
245 // Either the vector boson produced in association with the Higgs boson,
246 // or the ones produced from decays of top quarks produced with the Higgs
247 Particles leptonicVs;
248 if ( !is_uncatdV ){
249 if ( isVH(prodMode) && !quarkDecay(cat.V) ) leptonicVs += cat.V;
250 }else leptonicVs = std::move(uncatV_decays);
251 for ( const auto & W:Ws ) if ( W.genParticle()->end_vertex() && !quarkDecay(W) ) leptonicVs += W;
252
253 // Obtain all stable, final-state particles
254 const Particles FS = apply<FinalState>(event, "FS").particles();
256 FourMomentum sum(0,0,0,0), vSum(0,0,0,0), hSum(0,0,0,0);
257 for ( const Particle &p : FS ) {
258 // Add up the four momenta of all stable particles as a cross check
259 sum += p.momentum();
260 // ignore particles from the Higgs boson
261 if ( originateFrom(p,cat.higgs) ) { hSum += p.momentum(); continue; }
262 // Cross-check the V decay products for VH
263 if ( isVH(prodMode) && !is_uncatdV && originateFrom(p,Ws) ) vSum += p.momentum();
264 // ignore final state particles from leptonic V decays
265 if ( leptonicVs.size() && originateFrom(p,leptonicVs) ) continue;
266 // All particles reaching here are considered hadrons and will be used to build jets
267 hadrons += p;
268 }
269
270 cat.p4decay_higgs = hSum;
271 cat.p4decay_V = vSum;
272
273 FinalState fps_temp;
274 FastJets jets(fps_temp,
275#if RIVET_VERSION_CODE >= 40000
276 JetAlg::ANTIKT,
277#else
278 FastJets::ANTIKT,
279#endif // RIVET_VERSION_CODE
280 0.4 );
281 jets.calc(hadrons);
282
283 cat.jets25 = jets.jetsByPt( Cuts::pT > 25.0 );
284 cat.jets30 = jets.jetsByPt( Cuts::pT > 30.0 );
285
286 // check that four mometum sum of all stable particles satisfies momentum consevation
287 if ( sum.pt()>0.1 )
288 return error(cat,HTXS::MOMENTUM_CONSERVATION,"Four vector sum does not amount to pT=0, m=E=sqrt(s), but pT="+
289 std::to_string(sum.pt())+" GeV and m = "+std::to_string(sum.mass())+" GeV");
290
291 // check if V-boson was not included in the event record but decay particles were
292 // EFT contact interaction: return UNKNOWN for category but set all event/particle kinematics
293 if(is_uncatdV)
294 return error(cat,HTXS::VH_IDENTIFICATION,"Failed to identify associated V-boson!");
295
296 /*****
297 * Step 4.
298 * Classify and save output
299 */
300
301 // Apply the categorization categorization
302 cat.isZ2vvDecay = false;
303 if( (prodMode==HTXS::GG2ZH || prodMode==HTXS::QQ2ZH) && !quarkDecay(cat.V) && !ChLeptonDecay(cat.V) ) cat.isZ2vvDecay = true;
304 cat.stage0_cat = getStage0Category(prodMode,cat.higgs,cat.V);
305 cat.stage1_cat_pTjet25GeV = getStage1Category(prodMode,cat.higgs,cat.jets25,cat.V);
306 cat.stage1_cat_pTjet30GeV = getStage1Category(prodMode,cat.higgs,cat.jets30,cat.V);
307 cat.stage1_2_cat_pTjet25GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets25,cat.V);
308 cat.stage1_2_cat_pTjet30GeV = getStage1_2_Category(prodMode,cat.higgs,cat.jets30,cat.V);
309 cat.stage1_2_fine_cat_pTjet25GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets25,cat.V);
310 cat.stage1_2_fine_cat_pTjet30GeV = getStage1_2_Fine_Category(prodMode,cat.higgs,cat.jets30,cat.V);
311 cat.errorCode = HTXS::SUCCESS; ++m_errorCount[HTXS::SUCCESS];
312
313 return cat;
314 }
@ top
bool originateFrom(const Particle &p, const Particles &ptcls) const
Whether particle p originate from any of the ptcls.
HiggsClassification error(HiggsClassification &cat, HTXS::ErrorCode err, std::string msg="", int NmaxWarnings=20) const
Returns the classification object with the error code set.
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.
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.
@ 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
@ VH_IDENTIFICATION
failed to identify associated vector boson
GenVertex * signal_process_vertex(const GenEvent *e)
Definition GenEvent.h:626
list Particles
Definition makePlot.py:25
Particle_v1 Particle
Define the latest version of the particle class.

◆ error()

HiggsClassification Rivet::HiggsTemplateCrossSections::error ( HiggsClassification & cat,
HTXS::ErrorCode err,
std::string msg = "",
int NmaxWarnings = 20 ) const
inline

Returns the classification object with the error code set.

Prints an warning message, and keeps track of number of errors

Definition at line 110 of file HiggsTemplateCrossSections.h.

111 {
112 // Set the error, and keep statistics
113 cat.errorCode = err;
114 ++m_errorCount[err];
115
116 // Print warning message to the screen/log
117 static std::atomic<int> Nwarnings = 0;
118 if ( !msg.empty() && ++Nwarnings < NmaxWarnings )
120
121 return cat;
122 }
#define MSG_WARNING(ARG)
MsgStream & msg
Definition testRead.cxx:32

◆ finalize()

◆ getBin()

int Rivet::HiggsTemplateCrossSections::getBin ( double x,
const std::vector< double > & bins ) const
inline

Return bin index of x given the provided bin edges. 0=first bin, -1=underflow bin.

Definition at line 323 of file HiggsTemplateCrossSections.h.

323 {
324 if (bins.empty() || x < bins.front()) {
325 throw std::invalid_argument("Input value is out of bin range or bins vector is empty.");
326 }
327
328 for (size_t i = 1; i < bins.size(); ++i) {
329 if (x < bins[i]) {
330 return static_cast<int>(i - 1);
331 }
332 }
333
334 return static_cast<int>(bins.size() - 1);
335}
static const std::vector< std::string > bins
#define x

◆ getLastInstance()

Particle Rivet::HiggsTemplateCrossSections::getLastInstance ( const Particle & ptcl) const
inline

follow a "propagating" particle and return its last instance

Definition at line 48 of file HiggsTemplateCrossSections.h.

48 {
49 if ( ptcl.genParticle()->end_vertex() ) {
50 if ( !hasChild(ptcl.genParticle(),ptcl.pid()) ) return ptcl;
51 else return getLastInstance(ptcl.children()[0]);
52 }
53 return ptcl;
54 }

◆ getStage0Category()

HTXS::Stage0::Category Rivet::HiggsTemplateCrossSections::getStage0Category ( const HTXS::HiggsProdMode prodMode,
const Particle & higgs,
const Particle & V ) const
inline

Stage-0 HTXS categorization.

Definition at line 379 of file HiggsTemplateCrossSections.h.

381 {
382 using namespace HTXS::Stage0;
383 int ctrlHiggs = std::abs(higgs.rapidity())<2.5;
384 // Special cases first, qq→Hqq
385 if ( (prodMode==HTXS::WH||prodMode==HTXS::QQ2ZH) && quarkDecay(V) ) {
386 return ctrlHiggs ? VH2HQQ : VH2HQQ_FWDH;
387 } else if ( prodMode==HTXS::GG2ZH && quarkDecay(V) ) {
388 return Category(HTXS::GGF*10 + ctrlHiggs);
389 }
390 // General case after
391 return Category(prodMode*10 + ctrlHiggs);
392 }

◆ getStage1_2_Category()

HTXS::Stage1_2::Category Rivet::HiggsTemplateCrossSections::getStage1_2_Category ( const HTXS::HiggsProdMode prodMode,
const Particle & higgs,
const Jets & jets,
const Particle & V ) const
inline

Stage-1.2 categorization.

Definition at line 461 of file HiggsTemplateCrossSections.h.

464 {
465 using namespace HTXS::Stage1_2;
466 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
467 int vbfTopo = vbfTopology_Stage1_2(jets,higgs);
468
469 // 1. GGF Stage 1 categories
470 // Following YR4 write-up: XXXXX
471 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
472 if (fwdHiggs) return GG2H_FWDH;
473 if ( higgs.pt()>200 ) return Category(GG2H_PTH_200_300+getBin(higgs.pt(),{200,300,450,650}));
474 if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
475 if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
476 if (Njets>1){
477 //VBF topology
478 if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
479 //Njets >= 2jets without VBF topology (mjj<350)
480 return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
481 }
482 }
483
484 // 2. Electroweak qq->Hqq Stage 1.2 categories
485 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
486 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
487 int Njets=jets.size();
488 if (Njets==0) return QQ2HQQ_0J;
489 else if (Njets==1) return QQ2HQQ_1J;
490 else if (Njets>=2) {
491 double mjj = (jets[0].mom()+jets[1].mom()).mass();
492 if ( mjj < 60 ) return QQ2HQQ_GE2J_MJJ_0_60;
493 else if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_GE2J_MJJ_60_120;
494 else if ( 120 < mjj && mjj < 350 ) return QQ2HQQ_GE2J_MJJ_120_350;
495 else if ( mjj > 350 ) {
496 if (higgs.pt()>200) return QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200;
497 if(vbfTopo) return Category(QQ2HQQ_GE2J_MJJ_GT350_PTH_GT200+vbfTopo);
498 }
499 }
500 }
501 // 3. WH->Hlv categories
502 else if (prodMode==HTXS::WH) {
503 if (fwdHiggs) return QQ2HLNU_FWDH;
504 else if (V.pt()<75) return QQ2HLNU_PTV_0_75;
505 else if (V.pt()<150) return QQ2HLNU_PTV_75_150;
506 else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
507 // 150 < pTV/GeV < 250
509 }
510 // 4. qq->ZH->llH categories
511 else if (prodMode==HTXS::QQ2ZH) {
512 if (fwdHiggs) return QQ2HLL_FWDH;
513 else if (V.pt()<75) return QQ2HLL_PTV_0_75;
514 else if (V.pt()<150) return QQ2HLL_PTV_75_150;
515 else if (V.pt()>250) return QQ2HLL_PTV_GT250;
516 // 150 < pTV/GeV < 250
518 }
519 // 5. gg->ZH->llH categories
520 else if (prodMode==HTXS::GG2ZH ) {
521 if (fwdHiggs) return GG2HLL_FWDH;
522 else if (V.pt()<75) return GG2HLL_PTV_0_75;
523 else if (V.pt()<150) return GG2HLL_PTV_75_150;
524 else if (V.pt()>250) return GG2HLL_PTV_GT250;
526 }
527 // 6.ttH,bbH,tH categories
528 else if (prodMode==HTXS::TTH) {
529 if (fwdHiggs) return TTH_FWDH;
530 else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300}));
531 }
532 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
533 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
534 return UNKNOWN;
535 }
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...
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.

◆ getStage1_2_Fine_Category()

HTXS::Stage1_2_Fine::Category Rivet::HiggsTemplateCrossSections::getStage1_2_Fine_Category ( const HTXS::HiggsProdMode prodMode,
const Particle & higgs,
const Jets & jets,
const Particle & V ) const
inline

Stage-1.2_Fine categorization.

Definition at line 538 of file HiggsTemplateCrossSections.h.

541 {
542 using namespace HTXS::Stage1_2_Fine;
543 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
544 int vbfTopo = vbfTopology_Stage1_2_Fine(jets,higgs);
545
546 // 1. GGF Stage 1.2 categories
547 // Following YR4 write-up: XXXXX
548 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
549 if (fwdHiggs) return GG2H_FWDH;
550 if ( higgs.pt()>200 ){
551 if (Njets>0){
552 double pTHj = (jets[0].momentum()+higgs.momentum()).pt();
553 if( pTHj/higgs.pt()>0.15 ) return Category(GG2H_PTH_200_300_PTHJoverPTH_GT15+getBin(higgs.pt(),{200,300,450,650}));
554 else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
555 }
556 else return Category(GG2H_PTH_200_300_PTHJoverPTH_0_15+getBin(higgs.pt(),{200,300,450,650}));
557 }
558 if (Njets==0) return higgs.pt()<10 ? GG2H_0J_PTH_0_10 : GG2H_0J_PTH_GT10;
559 if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
560 if (Njets>1){
561 //double mjj = (jets[0].mom()+jets[1].mom()).mass();
562 double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
563 //VBF topology
564 if(vbfTopo) return Category(GG2H_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
565 //Njets >= 2jets without VBF topology (mjj<350)
566 if (pTHjj<25) return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_0_25+getBin(higgs.pt(),{0,60,120,200}));
567 else return Category(GG2H_GE2J_MJJ_0_350_PTH_0_60_PTHJJ_GT25+getBin(higgs.pt(),{0,60,120,200}));
568 }
569 }
570
571 // 2. Electroweak qq->Hqq Stage 1.2 categories
572 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
573 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
574 int Njets=jets.size();
575 if (Njets==0) return QQ2HQQ_0J;
576 else if (Njets==1) return QQ2HQQ_1J;
577 else if (Njets>=2) {
578 double mjj = (jets[0].mom()+jets[1].mom()).mass();
579 double pTHjj = (jets[0].momentum()+jets[1].momentum()+higgs.momentum()).pt();
580 if (mjj<350){
581 if (pTHjj<25) return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_0_25+getBin(mjj,{0,60,120,350}));
582 else return Category(QQ2HQQ_GE2J_MJJ_0_60_PTHJJ_GT25+getBin(mjj,{0,60,120,350}));
583 } else { //mjj>350 GeV
584 if (higgs.pt()<200){
585 return Category(QQ2HQQ_GE2J_MJJ_350_700_PTH_0_200_PTHJJ_0_25+vbfTopo-1);
586 } else {
587 return Category(QQ2HQQ_GE2J_MJJ_350_700_PTH_GT200_PTHJJ_0_25+vbfTopo-1);
588 }
589 }
590 }
591 }
592 // 3. WH->Hlv categories
593 else if (prodMode==HTXS::WH) {
594 if (fwdHiggs) return QQ2HLNU_FWDH;
595 int Njets=jets.size();
596 if (Njets==0) return Category(QQ2HLNU_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
597 if (Njets==1) return Category(QQ2HLNU_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
598 return Category(QQ2HLNU_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
599 }
600 // 4. qq->ZH->llH categories
601 else if (prodMode==HTXS::QQ2ZH) {
602 if (fwdHiggs) return QQ2HLL_FWDH;
603 int Njets=jets.size();
604 if (Njets==0) return Category(QQ2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
605 if (Njets==1) return Category(QQ2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
606 return Category(QQ2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
607 }
608 // 5. gg->ZH->llH categories
609 else if (prodMode==HTXS::GG2ZH ) {
610 if (fwdHiggs) return GG2HLL_FWDH;
611 int Njets=jets.size();
612 if (Njets==0) return Category(GG2HLL_PTV_0_75_0J+getBin(V.pt(),{0,75,150,250,400}));
613 if (Njets==1) return Category(GG2HLL_PTV_0_75_1J+getBin(V.pt(),{0,75,150,250,400}));
614 return Category(GG2HLL_PTV_0_75_GE2J+getBin(V.pt(),{0,75,150,250,400}));
615 }
616 // 6.ttH,bbH,tH categories
617 else if (prodMode==HTXS::TTH) {
618 if (fwdHiggs) return TTH_FWDH;
619 else return Category(TTH_PTH_0_60+getBin(higgs.pt(),{0,60,120,200,300,450}));
620 }
621 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
622 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
623 return UNKNOWN;
624 }
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,...

◆ getStage1Category()

HTXS::Stage1::Category Rivet::HiggsTemplateCrossSections::getStage1Category ( const HTXS::HiggsProdMode prodMode,
const Particle & higgs,
const Jets & jets,
const Particle & V ) const
inline

Stage-1 categorization.

Definition at line 395 of file HiggsTemplateCrossSections.h.

398 {
399 using namespace HTXS::Stage1;
400 int Njets=jets.size(), ctrlHiggs = std::abs(higgs.rapidity())<2.5, fwdHiggs = !ctrlHiggs;
401 double pTj1 = jets.size() ? jets[0].momentum().pt() : 0;
402 int vbfTopo = vbfTopology(jets,higgs);
403
404 // 1. GGF Stage 1 categories
405 // Following YR4 write-up: XXXXX
406 if (prodMode==HTXS::GGF || (prodMode==HTXS::GG2ZH && quarkDecay(V)) ) {
407 if (fwdHiggs) return GG2H_FWDH;
408 if (Njets==0) return GG2H_0J;
409 else if (Njets==1) return Category(GG2H_1J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
410 else if (Njets>=2) {
411 // events with pT_H>200 get priority over VBF cuts
412 if(higgs.pt()<=200){
413 if (vbfTopo==2) return GG2H_VBFTOPO_JET3VETO;
414 else if (vbfTopo==1) return GG2H_VBFTOPO_JET3;
415 }
416 // Njets >= 2jets without VBF topology
417 return Category(GG2H_GE2J_PTH_0_60+getBin(higgs.pt(),{0,60,120,200}));
418 }
419 }
420 // 2. Electroweak qq->Hqq Stage 1 categories
421 else if (prodMode==HTXS::VBF || ( isVH(prodMode) && quarkDecay(V)) ) {
422 if (std::abs(higgs.rapidity())>2.5) return QQ2HQQ_FWDH;
423 if (pTj1>200) return QQ2HQQ_PTJET1_GT200;
424 if (vbfTopo==2) return QQ2HQQ_VBFTOPO_JET3VETO;
425 if (vbfTopo==1) return QQ2HQQ_VBFTOPO_JET3;
426 double mjj = jets.size()>1 ? (jets[0].mom()+jets[1].mom()).mass():0;
427 if ( 60 < mjj && mjj < 120 ) return QQ2HQQ_VH2JET;
428 return QQ2HQQ_REST;
429 }
430 // 3. WH->Hlv categories
431 else if (prodMode==HTXS::WH) {
432 if (fwdHiggs) return QQ2HLNU_FWDH;
433 else if (V.pt()<150) return QQ2HLNU_PTV_0_150;
434 else if (V.pt()>250) return QQ2HLNU_PTV_GT250;
435 // 150 < pTV/GeV < 250
437 }
438 // 4. qq->ZH->llH categories
439 else if (prodMode==HTXS::QQ2ZH) {
440 if (fwdHiggs) return QQ2HLL_FWDH;
441 else if (V.pt()<150) return QQ2HLL_PTV_0_150;
442 else if (V.pt()>250) return QQ2HLL_PTV_GT250;
443 // 150 < pTV/GeV < 250
445 }
446 // 5. gg->ZH->llH categories
447 else if (prodMode==HTXS::GG2ZH ) {
448 if (fwdHiggs) return GG2HLL_FWDH;
449 if (V.pt()<150) return GG2HLL_PTV_0_150;
450 else if (jets.size()==0) return GG2HLL_PTV_GT150_0J;
452 }
453 // 6.ttH,bbH,tH categories
454 else if (prodMode==HTXS::TTH) return Category(TTH_FWDH+ctrlHiggs);
455 else if (prodMode==HTXS::BBH) return Category(BBH_FWDH+ctrlHiggs);
456 else if (prodMode==HTXS::TH ) return Category(TH_FWDH+ctrlHiggs);
457 return UNKNOWN;
458 }
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,...

◆ hasChild()

bool Rivet::HiggsTemplateCrossSections::hasChild ( HepMC::ConstGenParticlePtr ptcl,
int pdgID ) const
inline

Checks whether the input particle has a child with a given PDGID.

Definition at line 75 of file HiggsTemplateCrossSections.h.

75 {
76 for (const Particle& child:Particle(*ptcl).children())
77 if (child.pid()==pdgID) return true;
78 return false;
79 }

◆ hasParent()

bool Rivet::HiggsTemplateCrossSections::hasParent ( HepMC::ConstGenParticlePtr ptcl,
int pdgID ) const
inline

Checks whether the input particle has a parent with a given PDGID.

Definition at line 82 of file HiggsTemplateCrossSections.h.

82 {
83 for (auto parent:Rivet::HepMCUtils::particles(ptcl->production_vertex(),Relatives::PARENTS))
84 if (parent->pdg_id()==pdgID) return true;
85 return false;
86 }

◆ init()

void Rivet::HiggsTemplateCrossSections::init ( )
inline

default Rivet Analysis::init method Booking of histograms, initializing Rivet projection Extracts Higgs production mode from shell variable if not set manually using setHiggsProdMode

Definition at line 638 of file HiggsTemplateCrossSections.h.

638 {
639 printf("==============================================================\n");
640 printf("======== HiggsTemplateCrossSections Initialization =========\n");
641 printf("==============================================================\n");
642 // check that the production mode has been set
643 // if running in standalone Rivet the production mode is set through an env variable
645 char *pm_env = getenv("HIGGSPRODMODE");
646 string pm(pm_env==nullptr?"":pm_env);
647 if ( pm == "GGF" ) m_HiggsProdMode = HTXS::GGF;
648 else if ( pm == "VBF" ) m_HiggsProdMode = HTXS::VBF;
649 else if ( pm == "WH" ) m_HiggsProdMode = HTXS::WH;
650 else if ( pm == "ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
651 else if ( pm == "QQ2ZH" ) m_HiggsProdMode = HTXS::QQ2ZH;
652 else if ( pm == "GG2ZH" ) m_HiggsProdMode = HTXS::GG2ZH;
653 else if ( pm == "TTH" ) m_HiggsProdMode = HTXS::TTH;
654 else if ( pm == "BBH" ) m_HiggsProdMode = HTXS::BBH;
655 else if ( pm == "TH" ) m_HiggsProdMode = HTXS::TH;
656 else {
657 MSG_WARNING("No HIGGSPRODMODE shell variable found. Needed when running Rivet stand-alone.");
658 }
659 }
660
661 // Projections for final state particles
662 const FinalState FS;
663 declare(FS,"FS");
664
665 // initialize the histograms with for each of the stages
667 m_sumw = 0.0;
668 printf("==============================================================\n");
669 printf("======== Higgs prod mode %d =========\n",m_HiggsProdMode);
670 printf("======== Sucessful Initialization =========\n");
671 printf("==============================================================\n");
672 }
std::string getenv(const std::string &variableName)
get an environment variable

◆ initializeHistos()

void Rivet::HiggsTemplateCrossSections::initializeHistos ( )
inline

Definition at line 758 of file HiggsTemplateCrossSections.h.

758 {
759 book(m_hist_stage0,"HTXS_stage0",20,0,20);
760 book(m_hist_stage1_pTjet25,"HTXS_stage1_pTjet25",40,0,40);
761 book(m_hist_stage1_pTjet30,"HTXS_stage1_pTjet30",40,0,40);
762 book(m_hist_stage1_2_pTjet25,"HTXS_stage1_2_pTjet25",57,0,57);
763 book(m_hist_stage1_2_pTjet30,"HTXS_stage1_2_pTjet30",57,0,57);
764 book(m_hist_stage1_2_fine_pTjet25,"HTXS_stage1_2_fine_pTjet25",113,0,113);
765 book(m_hist_stage1_2_fine_pTjet30,"HTXS_stage1_2_fine_pTjet30",113,0,113);
766 book(m_hist_pT_Higgs,"pT_Higgs",80,0,400);
767 book(m_hist_y_Higgs,"y_Higgs",80,-4,4);
768 book(m_hist_pT_V,"pT_V",80,0,400);
769 book(m_hist_pT_jet1,"pT_jet1",80,0,400);
770 book(m_hist_deltay_jj ,"deltay_jj",50,0,10);
771 book(m_hist_dijet_mass,"m_jj",50,0,2000);
772 book(m_hist_pT_Hjj,"pT_Hjj",50,0,250);
773 book(m_hist_Njets25,"Njets25",10,0,10);
774 book(m_hist_Njets30,"Njets30",10,0,10);
775 book(m_hist_isZ2vv,"isZ2vv",2,0,2);
776 }

◆ isVH()

bool Rivet::HiggsTemplateCrossSections::isVH ( HTXS::HiggsProdMode p) const
inline

Whether the Higgs is produced in association with a vector boson (VH)

Definition at line 376 of file HiggsTemplateCrossSections.h.

376{ return p==HTXS::WH || p==HTXS::QQ2ZH || p==HTXS::GG2ZH; }

◆ originateFrom() [1/2]

bool Rivet::HiggsTemplateCrossSections::originateFrom ( const Particle & p,
const Particle & p2 ) const
inline

Whether particle p originates from p2.

Definition at line 70 of file HiggsTemplateCrossSections.h.

70 {
71 Particles ptcls = {p2}; return originateFrom(p,ptcls);
72 }

◆ originateFrom() [2/2]

bool Rivet::HiggsTemplateCrossSections::originateFrom ( const Particle & p,
const Particles & ptcls ) const
inline

Whether particle p originate from any of the ptcls.

Definition at line 57 of file HiggsTemplateCrossSections.h.

57 {
58 auto prodVtx = p.genParticle()->production_vertex();
59 if (prodVtx == nullptr) return false;
60 // for each ancestor, check if it matches any of the input particles
61 for (auto ancestor:Rivet::HepMCUtils::particles(std::move(prodVtx),Relatives::ANCESTORS)){
62 for ( const auto & part:ptcls )
63 if ( ancestor==part.genParticle() ) return true;
64 }
65 // if we get here, no ancestor matched any input particle
66 return false;
67 }

◆ printClassificationSummary()

void Rivet::HiggsTemplateCrossSections::printClassificationSummary ( )
inline

Definition at line 724 of file HiggsTemplateCrossSections.h.

724 {
725 MSG_INFO (" ====================================================== ");
726 MSG_INFO (" Higgs Template X-Sec Categorization Tool ");
727 MSG_INFO (" Status Code Summary ");
728 MSG_INFO (" ====================================================== ");
729 bool allSuccess = (numEvents()==m_errorCount[HTXS::SUCCESS]);
730 if ( allSuccess ) MSG_INFO (" >>>> All "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized!");
731 else{
732 MSG_INFO (" >>>> "<< m_errorCount[HTXS::SUCCESS] <<" events successfully categorized");
733 MSG_INFO (" >>>> --> the following errors occured:");
734 MSG_INFO (" >>>> "<< m_errorCount[HTXS::PRODMODE_DEFINED] <<" had an undefined Higgs production mode.");
735 MSG_INFO (" >>>> "<< m_errorCount[HTXS::MOMENTUM_CONSERVATION] <<" failed momentum conservation.");
736 MSG_INFO (" >>>> "<< m_errorCount[HTXS::HIGGS_IDENTIFICATION] <<" failed to identify a valid Higgs boson.");
737 MSG_INFO (" >>>> "<< m_errorCount[HTXS::HS_VTX_IDENTIFICATION] <<" failed to identify the hard scatter vertex.");
738 MSG_INFO (" >>>> "<< m_errorCount[HTXS::VH_IDENTIFICATION] <<" VH: to identify a valid V-boson.");
739 MSG_INFO (" >>>> "<< m_errorCount[HTXS::TOP_W_IDENTIFICATION] <<" failed to identify valid Ws from top decay.");
740 }
741 MSG_INFO (" ====================================================== ");
742 MSG_INFO (" ====================================================== ");
743 }
#define MSG_INFO(ARG)

◆ quarkDecay()

bool Rivet::HiggsTemplateCrossSections::quarkDecay ( const Particle & p) const
inline

Return true is particle decays to quarks.

Definition at line 89 of file HiggsTemplateCrossSections.h.

89 {
90 for (const Particle& child:p.children())
91 if (PID::isQuark(child.pid())) return true;
92 return false;
93 }

◆ setHiggsProdMode()

void Rivet::HiggsTemplateCrossSections::setHiggsProdMode ( HTXS::HiggsProdMode prodMode)
inline

Sets the Higgs production mode.

Definition at line 633 of file HiggsTemplateCrossSections.h.

633{ m_HiggsProdMode = prodMode; }

◆ vbfTopology()

int Rivet::HiggsTemplateCrossSections::vbfTopology ( const Jets & jets,
const Particle & higgs ) const
inline

VBF topolog selection 0 = fail loose selction: m_jj > 400 GeV and Dy_jj > 2.8 1 pass loose, but fail additional cut pT(Hjj)<25.

2 pass tight selection

Definition at line 340 of file HiggsTemplateCrossSections.h.

340 {
341 if (jets.size()<2) return 0;
342 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
343 bool VBFtopo = (j1+j2).mass() > 400.0 && std::abs(j1.rapidity()-j2.rapidity()) > 2.8;
344 return VBFtopo ? (j1+j2+higgs.momentum()).pt()<25 ? 2 : 1 : 0;
345 }

◆ vbfTopology_Stage1_2()

int Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_2 ( const Jets & jets,
const Particle & higgs ) const
inline

VBF topology selection 0 = fail loose selection: m_jj > 350 GeV 1 pass loose, but fail additional cut pT(Hjj)<25.

2 pass pT(Hjj)>25 selection 3 pass tight (m_jj>700 GeV), but fail additional cut pT(Hjj)<25. 4 pass pT(Hjj)>25 selection

Definition at line 350 of file HiggsTemplateCrossSections.h.

350 {
351 if (jets.size()<2) return 0;
352 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
353 double mjj = (j1+j2).mass();
354 if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
355 else if(mjj>700) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
356 else return 0;
357 }

◆ vbfTopology_Stage1_2_Fine()

int Rivet::HiggsTemplateCrossSections::vbfTopology_Stage1_2_Fine ( const Jets & jets,
const Particle & higgs ) const
inline

VBF topology selection for Stage1_2 0 = fail loose selection: m_jj > 350 GeV 1 pass loose, but fail additional cut pT(Hjj)<25.

2 pass pT(Hjj)>25 selection 3 pass 700<m_jj<1000 GeV, but fail additional cut pT(Hjj)<25. 4 pass pT(Hjj)>25 selection 5 pass 1000<m_jj<1500 GeV, but fail additional cut pT(Hjj)<25. 6 pass pT(Hjj)>25 selection 7 pass m_jj>1500 GeV, but fail additional cut pT(Hjj)<25. 8 pass pT(Hjj)>25 selection

Definition at line 364 of file HiggsTemplateCrossSections.h.

364 {
365 if (jets.size()<2) return 0;
366 const FourMomentum &j1=jets[0].momentum(), &j2=jets[1].momentum();
367 double mjj = (j1+j2).mass();
368 if(mjj>350 && mjj<=700) return (j1+j2+higgs.momentum()).pt()<25 ? 1 : 2;
369 else if(mjj>700 && mjj<=1000) return (j1+j2+higgs.momentum()).pt()<25 ? 3 : 4;
370 else if(mjj>1000 && mjj<=1500) return (j1+j2+higgs.momentum()).pt()<25 ? 5 : 6;
371 else if(mjj>1500) return (j1+j2+higgs.momentum()).pt()<25 ? 7 : 8;
372 else return 0;
373 }

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::array<std::atomic<size_t>, HTXS::NUM_ERRORCODES> m_errorCount Rivet::HiggsTemplateCrossSections::ATLAS_THREAD_SAFE {}
mutableprivate

Definition at line 786 of file HiggsTemplateCrossSections.h.

786{};

◆ m_HiggsProdMode

HTXS::HiggsProdMode Rivet::HiggsTemplateCrossSections::m_HiggsProdMode
private

Definition at line 785 of file HiggsTemplateCrossSections.h.

◆ m_hist_deltay_jj

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_deltay_jj
private

Definition at line 793 of file HiggsTemplateCrossSections.h.

◆ m_hist_dijet_mass

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_dijet_mass
private

Definition at line 793 of file HiggsTemplateCrossSections.h.

◆ m_hist_isZ2vv

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_isZ2vv
private

Definition at line 795 of file HiggsTemplateCrossSections.h.

◆ m_hist_Njets25

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_Njets25
private

Definition at line 794 of file HiggsTemplateCrossSections.h.

◆ m_hist_Njets30

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_Njets30
private

Definition at line 794 of file HiggsTemplateCrossSections.h.

◆ m_hist_pT_Higgs

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_pT_Higgs
private

Definition at line 791 of file HiggsTemplateCrossSections.h.

◆ m_hist_pT_Hjj

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_pT_Hjj
private

Definition at line 793 of file HiggsTemplateCrossSections.h.

◆ m_hist_pT_jet1

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_pT_jet1
private

Definition at line 792 of file HiggsTemplateCrossSections.h.

◆ m_hist_pT_V

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_pT_V
private

Definition at line 792 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage0

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage0
private

Definition at line 787 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_2_fine_pTjet25

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_fine_pTjet25
private

Definition at line 790 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_2_fine_pTjet30

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_fine_pTjet30
private

Definition at line 790 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_2_pTjet25

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_pTjet25
private

Definition at line 789 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_2_pTjet30

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_2_pTjet30
private

Definition at line 789 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_pTjet25

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_pTjet25
private

Definition at line 788 of file HiggsTemplateCrossSections.h.

◆ m_hist_stage1_pTjet30

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_stage1_pTjet30
private

Definition at line 788 of file HiggsTemplateCrossSections.h.

◆ m_hist_y_Higgs

Histo1DPtr Rivet::HiggsTemplateCrossSections::m_hist_y_Higgs
private

Definition at line 791 of file HiggsTemplateCrossSections.h.

◆ m_sumw

double Rivet::HiggsTemplateCrossSections::m_sumw =0.0
private

Definition at line 784 of file HiggsTemplateCrossSections.h.


The documentation for this class was generated from the following file: