ATLAS Offline Software
Loading...
Searching...
No Matches
FJvtEfficiencyTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
6
9
10namespace CP {
15
18 ATH_CHECK(m_evtInfoKey.initialize());
22 ATH_MSG_ERROR("failed to set up fJvt systematics");
23 return StatusCode::FAILURE;
24 }
25 return StatusCode::SUCCESS;
26 }
27
29 if (sys.find(fJvtEfficiencyUp) != sys.end())
31 else if (sys.find(fJvtEfficiencyDown) != sys.end())
33 else
35 return StatusCode::SUCCESS;
36 }
37
40 if (!isInRange(jet)) {
41 sf = -1;
43 }
45 if (!m_accIsHS->isAvailable(jet)) {
46 ATH_MSG_ERROR("Truth tagging required but not available");
48 }
49 if (!(*m_accIsHS)(jet)) {
50 sf = 1;
51 return CorrectionCode::Ok;
52 }
53 }
54 auto evtInfo = SG::makeHandle(m_evtInfoKey);
55 if (!evtInfo.isValid()) {
56 ATH_MSG_ERROR("Failed to retrieve " << m_evtInfoKey.key());
58 }
59 return getEffImpl(jet.pt(), evtInfo->actualInteractionsPerCrossing(), sf);
60 }
61
64 if (!isInRange(jet)) {
65 sf = -1;
67 }
69 if (!m_accIsHS->isAvailable(jet)) {
70 ATH_MSG_ERROR("Truth tagging required but not available");
72 }
73 if (!(*m_accIsHS)(jet)) {
74 sf = 1;
75 return CorrectionCode::Ok;
76 }
77 }
78 auto evtInfo = SG::makeHandle(m_evtInfoKey);
79 if (!evtInfo.isValid()) {
80 ATH_MSG_ERROR("Failed to retrieve " << m_evtInfoKey.key());
82 }
83 return getIneffImpl(jet.pt(), evtInfo->actualInteractionsPerCrossing(), sf);
84 }
85} // namespace CP
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
Handle class for reading a decoration on an object.
Handle class for reading from StoreGate.
Return value from object correction CP tools.
@ Error
Some error happened during the object correction.
@ OutOfValidityRange
Input object is out of validity range.
@ Ok
The correction was done successfully.
Gaudi::Property< std::string > m_file
virtual CorrectionCode getInefficiencyScaleFactor(const xAOD::Jet &jet, float &sf) const override
Calculate the inefficiency scale factor for the provided jet.
virtual CorrectionCode getEfficiencyScaleFactor(const xAOD::Jet &jet, float &sf) const override
Calculate the efficiency scale factor for the provided jet.
SG::ReadHandleKey< xAOD::EventInfo > m_evtInfoKey
Gaudi::Property< std::string > m_wp
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
FJvtEfficiencyTool(const std::string &name)
virtual StatusCode sysApplySystematicVariation(const CP::SystematicSet &sys) override
effects: configure this tool for the given list of systematic variations.
Gaudi::Property< float > m_minEta
StatusCode initHists(const std::string &file, const std::string &wp)
Read the input histograms. Passing an empty 'file' string uses dummy SFs.
Gaudi::Property< float > m_maxEta
CorrectionCode getIneffImpl(float x, float y, float &sf) const
Gaudi::Property< bool > m_doTruthRequirement
CorrectionCode getEffImpl(float x, float y, float &sf) const
std::optional< SG::AuxElement::ConstAccessor< char > > m_accIsHS
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
bool isInRange(const xAOD::Jet &jet) const
Class to wrap a set of SystematicVariations.
StatusCode addAffectingSystematic(const SystematicVariation &systematic, bool recommended)
effects: add a systematic to the list of registered systematics.
Select isolated Photons, Electrons and Muons.
static const SystematicVariation fJvtEfficiencyDown("JET_fJvtEfficiency", -1)
static const SystematicVariation fJvtEfficiencyUp("JET_fJvtEfficiency", 1)
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
Jet_v1 Jet
Definition of the current "jet version".