ATLAS Offline Software
Loading...
Searching...
No Matches
SkimmingToolHIGG5VBF.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// SkimmingToolHIGG5VBF.cxx, (c) ATLAS Detector software
8// Author: Yasu Okumura (yasuyuki.okumura@cern.ch)
9// Based on DerivationFramework::SkimmingToolExample and DerivationFramework::SkimmingToolHIGG2
10
13#include <vector>
14#include <string>
15
16#include "CLHEP/Units/SystemOfUnits.h"
17
18// Constructor
20 const std::string& n,
21 const IInterface* p) :
22 base_class(t, n, p),
23 m_trigDecisionTool("Trig::TrigDecisionTool/TrigDecisionTool"),
24 m_ntot(0),
25 m_npass(0)
26{
27
28 // for jet multiplicity requirement
29 declareProperty("JetContainerKey", m_jetSGKey="AntiKt4EMTopoJets");
30 declareProperty("CalibedJetMomentKey", m_calibedJetMomKey="DFCommonJets_Calib");
31
32 // for jet multiplicity
33 declareProperty("ReqNAllJets", m_reqNAllJets=false);
34 declareProperty("NumberOfAllJets", m_nAllJets=4);
35 declareProperty("AllJetPtCut", m_allJetPtCut=40.*CLHEP::GeV);
36 declareProperty("AllJetEtaCut", m_allJetEtaCut=4.9);
37
38 declareProperty("ReqNCentralJets", m_reqNCentralJets=false);
39 declareProperty("NumberOfCentralJets", m_nCentralJets=2);
40 declareProperty("CentralJetPtCut", m_centralJetPtCut=40.*CLHEP::GeV);
41 declareProperty("CentralJetEtaCut", m_centralJetEtaCut=2.6);
42
43 // for trigger requirement
44 declareProperty("ReqTrigger", m_reqTrigger=false);
45 declareProperty("Triggers", m_triggers=std::vector<std::string>());
46
47 // for Mjj requirement ()
48 declareProperty("ReqVBFMjj", m_reqVbfMjj=false); // logical "ORed" with ReqPhoton
49 declareProperty("MjjCut", m_vbfMjjCut=0.*CLHEP::GeV);
50
51 declareProperty("DoDebug", m_debug=false);
52
53 // photon requirement (p. rose)
54 declareProperty("PhotonContainerKey", m_phSGKey="Photons");
55 declareProperty("ReqPhoton", m_reqPh=false); // logical "ORed" with ReqVBFMjj
56 declareProperty("PhotonPtCut", m_phPtCut=0.*CLHEP::GeV);
57 declareProperty("CentralPhotonEtaCut", m_centralPhEtaCut=2.6);
58
59}
60
61// Destructor
64
65// Athena initialize and finalize
67{
68 ATH_MSG_VERBOSE("initialize() ...");
69
70 // trigger decision tool
71 if(m_trigDecisionTool.retrieve(DisableTool{!m_reqTrigger}).isFailure()) {
72 ATH_MSG_FATAL("Failed to retrieve tool: " << m_trigDecisionTool);
73 return StatusCode::FAILURE;
74 }
75 ATH_MSG_INFO("Retrieved tool: " << m_trigDecisionTool);
76
77 return StatusCode::SUCCESS;
78}
79
81{
82 ATH_MSG_VERBOSE("finalize() ...");
83 ATH_MSG_INFO("Processed " << m_ntot << " events, " << m_npass << " events passed filter ");
84 return StatusCode::SUCCESS;
85}
86
87// The filter itself
89{
90 m_ntot++;
91 bool acceptEvent(true);
92
93 StatusCode sc(StatusCode::SUCCESS);
94
95 bool isTriggerFired(m_triggers.empty());
96 for(unsigned int i(0); i<m_triggers.size(); i++) {
97 if(m_trigDecisionTool->isPassed(m_triggers.at(i))) {
98 isTriggerFired = true;
99 break;
100 }
101 }
102
103 // (1) Count Jet Multiplicity
104 std::vector<TLorentzVector> goodAllJets;
105 std::vector<TLorentzVector> goodCentralJets;
106
107 const xAOD::JetContainer *jets(nullptr);
108 ATH_CHECK(evtStore()->retrieve(jets, m_jetSGKey), false);
109 xAOD::JetContainer::const_iterator jet_itr(jets->begin());
110 xAOD::JetContainer::const_iterator jet_end(jets->end());
111 for(; jet_itr != jet_end; ++jet_itr) {
112 TLorentzVector jetP4 = getCalibedJets( (*jet_itr) );
113 if(this->checkAllJetQuality(jetP4)) { goodAllJets.push_back(jetP4); }
114 if(this->checkCentralJetQuality(jetP4)) { goodCentralJets.push_back(jetP4); }
115 }
116
117 // (2) evaluate maximum Mjj in the event
118 double maxM = 0.;
119 for(unsigned int jet_i = 0; jet_i<goodAllJets.size(); jet_i++) {
120 const TLorentzVector& iP4 = goodAllJets.at(jet_i);
121
122 for(unsigned int jet_k=jet_i+1; jet_k<goodAllJets.size(); jet_k++) {
123 const TLorentzVector& kP4 = goodAllJets.at(jet_k);
124
125 const TLorentzVector jjP4 = iP4 + kP4;
126 const double jjM = jjP4.M();
127
128 if (maxM<jjM) {maxM=jjM;}
129 }
130 }
131
132 //get max pt of any photon in the event (p. rose)
133 double maxPhPt=0.;
134 if(m_reqPh){
135 const xAOD::PhotonContainer *phots(nullptr);
136 ATH_CHECK(evtStore()->retrieve(phots, m_phSGKey), false);
137 for (const auto* ph : *phots){
138 if(abs(ph->eta())<m_centralPhEtaCut)
139 if(ph->pt()>maxPhPt) maxPhPt = ph->pt();
140 }//for
141 }//if
142
143
144
145 const bool passNAllJet = (goodAllJets.size()>=m_nAllJets);
146 const bool passNCentralJet = (goodCentralJets.size()>=m_nCentralJets);
147 const bool passMjjCut = (maxM>m_vbfMjjCut);
148 const bool passPhPtCut = (maxPhPt>m_phPtCut);
149
150 if (m_reqNAllJets) { if (not passNAllJet) {acceptEvent=false;} }
151 if (m_reqNCentralJets) { if (not passNCentralJet) {acceptEvent=false;} }
152 //if (m_reqVbfMjj) { if (not passMjjCut) {acceptEvent=false;} }
153 if (m_reqTrigger) { if (not isTriggerFired) {acceptEvent=false;} }
154 //vbf+gamma addition -- logical OR of mjj and phpt cut
155 if(m_reqVbfMjj and not m_reqPh) { if (not passMjjCut) {acceptEvent=false;} }
156 if(m_reqPh and not m_reqVbfMjj) { if (not passPhPtCut) {acceptEvent=false;} }
157 if(m_reqPh and m_reqVbfMjj) { if (not (passMjjCut or passPhPtCut) ) {acceptEvent=false;} }
158
159 if (acceptEvent) {m_npass++;}
160
161 if (m_debug) {
162 printf("dbg> L%3d : event accepted [%s] "
163 "(NJets=%2d [%6s], NCentralJets=%2d [%6s] Mjj=%10.1f [%6s], Trigger [%6s]) \n",
164 __LINE__,
165 acceptEvent? "Y" : "N",
166 (int)goodAllJets.size(), passNAllJet ? "PASSED" : "FAILED",
167 (int)goodCentralJets.size(), passNCentralJet ? "PASSED" : "FAILED",
168 maxM, passMjjCut ? "PASSED" : "FAILED",
169 isTriggerFired ? "PASSED" : "FAILED"
170 );
171 }
172
173 return acceptEvent;
174}
175
177{
178 if(jet.Pt()<m_allJetPtCut) return false;
179 if(fabs(jet.Eta())>m_allJetEtaCut) return false;
180
181 return true;
182}
183
185{
186 if(jet.Pt()<m_centralJetPtCut) return false;
187 if(fabs(jet.Eta())>m_centralJetEtaCut) return false;
188
189 return true;
190}
191
192TLorentzVector
194{
195 TLorentzVector rc;
196
197 if(!jet) return rc;
198
203 const float& pt =ptAcc(*jet);
204 const float& eta=etaAcc(*jet);
205 const float& phi=phiAcc(*jet);
206 const float& m =mAcc(*jet);
207
208 rc.SetPtEtaPhiM(pt, eta, phi, m);
209
210 return rc;
211}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
Helper class to provide constant type-safe access to aux data.
static Double_t sc
static Double_t rc
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
TLorentzVector getCalibedJets(const xAOD::Jet *jet) const
SkimmingToolHIGG5VBF(const std::string &t, const std::string &n, const IInterface *p)
Constructor with parameters.
bool checkCentralJetQuality(const TLorentzVector &jet) const
bool checkAllJetQuality(const TLorentzVector &jet) const
virtual bool eventPassesFilter() const override
Check that the current event passes this filter.
ToolHandle< Trig::TrigDecisionTool > m_trigDecisionTool
Helper class to provide constant type-safe access to aux data.
Jet_v1 Jet
Definition of the current "jet version".
PhotonContainer_v1 PhotonContainer
Definition of the current "photon container version".
JetContainer_v1 JetContainer
Definition of the current "jet container version".