ATLAS Offline Software
Loading...
Searching...
No Matches
UncertaintyComponent.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8#include "TFile.h"
9
10namespace jet
11{
12
13bool operator < (const UncertaintyComponent& componentA, const UncertaintyComponent& componentB)
14{
15 return componentA.getName() < componentB.getName();
16}
17
18bool operator == (const UncertaintyComponent& componentA, const UncertaintyComponent& componentB)
19{
20 return componentA.getName() == componentB.getName();
21}
22
23
25// //
26// Constructor/destructor/initialization //
27// //
29
31 : asg::AsgMessaging(name)
32 , m_isInit(false)
33 , m_uncHistName("")
34 , m_validHistName("")
35 , m_scaleVar(CompScaleVar::UNKNOWN)
36 , m_topology(JetTopology::UNKNOWN)
37 , m_energyScale(1)
38 , m_interpolate(Interpolate::UNKNOWN)
39 , m_splitNumber(0)
41 , m_uncHist(nullptr)
42 , m_validHist(nullptr)
43{
45}
46
47UncertaintyComponent::UncertaintyComponent(const ComponentHelper& component, const size_t numHist)
48 : asg::AsgMessaging(component.name.Data())
49 , m_isInit(false)
50 , m_uncHistName(!component.uncNames.empty()?component.uncNames.at(0):"NONE")
51 , m_validHistName(component.validName)
52 , m_scaleVar(component.scaleVar)
53 , m_topology(component.topology)
54 , m_energyScale(component.energyScale)
55 , m_interpolate(component.interpolate)
56 , m_splitNumber(component.splitNum)
57 , m_uncHist(nullptr)
58 , m_validHist(nullptr)
59{
60 ATH_MSG_DEBUG("Creating UncertaintyComponent named " << m_uncHistName.Data());
61 if (m_uncHistName == "")
62 ATH_MSG_FATAL("Cannot create an UncertaintyComponent with an empty uncHistName");
64 ATH_MSG_FATAL("Cannot create an UncertaintyComponent scaling an UNKNOWN variable");
65 if (numHist != component.uncNames.size())
66 ATH_MSG_FATAL("Expected " << numHist << " uncertainty histograms, but received " << component.uncNames.size() << " : " << m_uncHistName.Data());
67}
68
70 : asg::AsgMessaging(Form("%s_copy",toCopy.m_uncHistName.Data()))
71 , m_isInit(toCopy.m_isInit)
74 , m_scaleVar(toCopy.m_scaleVar)
75 , m_topology(toCopy.m_topology)
79 , m_uncHist(nullptr)
80 , m_validHist(nullptr)
81{
82
83 if (toCopy.m_uncHist)
85 if (toCopy.m_validHist)
87}
88
95
96StatusCode UncertaintyComponent::initialize(TFile* histFile)
97{
98 // Prevent double-initialization
99 if (m_isInit)
100 {
101 ATH_MSG_ERROR("Component is already initialized: " << getName().Data());
102 return StatusCode::FAILURE;
103 }
104
105 // Create the histograms
107 if (!m_uncHist)
108 {
109 ATH_MSG_ERROR("Failed to create uncertainty histogram for component: " << getName().Data());
110 return StatusCode::FAILURE;
111 }
112 if (m_validHistName != "")
113 {
115 if (!m_validHist)
116 {
117 ATH_MSG_ERROR("Failed to create validity histogram of name \"" << getValidName().Data() << "\" for component: " << getName().Data());
118 return StatusCode::FAILURE;
119 }
120 }
121
122 // Initialize the histograms
123 if (m_uncHist->initialize(histFile).isFailure()) return StatusCode::FAILURE;
124 if (m_validHist && m_validHist->initialize(histFile).isFailure()) return StatusCode::FAILURE;
125
126 // Print a summary
127 ATH_MSG_DEBUG("Successfully initialized UncertaintyComponent named " << getName().Data());
128
129 m_isInit = true;
130 return StatusCode::SUCCESS;
131}
132
133
135// //
136// Methods to test for special cases //
137// //
139
141{
142 if (!m_isInit)
143 {
144 ATH_MSG_ERROR("Cannot call method before initialization, component: " << getName().Data());
145 return false;
146 }
147
148 const TH1* histo = m_uncHist->getHisto();
149 return !(fabs(histo->GetMinimum()) > 1.e-8 || fabs(histo->GetMaximum()) > 1.e-8);
150}
151
152
154// //
155// Uncertainty/validity retrieval methods //
156// //
158
160{
161 if (!m_isInit)
162 {
163 ATH_MSG_ERROR("Component hasn't been initialized: " << getName().Data());
164 return false;
165 }
166 return getValidityImpl(jet,eInfo);
167}
168
170{
171 if (!m_isInit)
172 {
173 ATH_MSG_ERROR("Component hasn't been initialized: " << getName().Data());
174 return JESUNC_ERROR_CODE;
175 }
177}
178
179bool UncertaintyComponent::getValidUncertainty(double& unc, const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const
180{
181 if (getValidity(jet,eInfo))
182 {
183 unc = getUncertainty(jet,eInfo);
184 return true;
185 }
186 return false;
187}
188
189
191// //
192// Split factor methods //
193// //
195
197{
198 // Immediate return for the most common case
199 if (!m_splitNumber)
200 return 1;
201
202 // SplitNumber was specified, we have to calculate the factor
203 double splitFactor = 1;
204 const TH1* histo = m_uncHist->getHisto();
205
206 if (m_splitNumber == 1 || m_splitNumber == -1)
207 {
208 // Straight line in log(pT) from 0 to 1
209 // y = mx+b
210 // m = 1/[log(max)-log(min)]
211 // x = log(min) --> b = -m*log(min)
212
213 const double minPt = histo->GetXaxis()->GetBinLowEdge(1);
214 const double maxPt = histo->GetXaxis()->GetBinLowEdge(histo->GetNbinsX()+1);
215 const double valPt = jet.pt()*m_energyScale;
216
217 const double slope = 1./(log(maxPt)-log(minPt));
218 const double intercept = -slope*log(minPt);
219
220 splitFactor = slope*log(valPt <= minPt ? minPt+1.e-3 : valPt >= maxPt ? maxPt-1.e-3 : valPt)+intercept;
221 }
222 else if (m_splitNumber == 2 || m_splitNumber == -2)
223 {
224 // Straight line in |eta| from 0 to 1
225 // y = mx + b
226 // m = 1/(max - min)
227 // x = min --> b = -m*min
228 const double minEta = 0;
229 const double maxEta = 4.5;
230 const double absEta = fabs(jet.eta());
231
232 const double slope = 1./(maxEta - minEta);
233 const double intercept = -slope*minEta;
234
235 splitFactor = slope*(absEta <= minEta ? minEta+1.e-3 : absEta >= maxEta ? maxEta-1.e-3 : absEta)+intercept;
236 }
237 else if (m_splitNumber == 3 || m_splitNumber == -3)
238 {
239 // Increasing with log(pT) and increasing with |eta|
240 // z = mx + ny + c
241 // z(min,min) = 0
242 // z(max,max) = 1
243 // Linear in both dimensions means need factor of 1/2 in a single dimension
244 // m = 0.5/[log(maxPt)-log(minPt)]
245 // n = 0.5/(maxEta - minEta)
246 // c = -minPt*m - minEta*n
247 // 2*z = (logx-logxmin)/(logxmax-logxmin) + (y-ymin)/(ymax-ymin)
248 const double minPt = histo->GetXaxis()->GetBinLowEdge(1);
249 const double maxPt = histo->GetXaxis()->GetBinLowEdge(histo->GetNbinsX()+1);
250 const double valPt = jet.pt()*m_energyScale;
251
252 const double minEta = 0;
253 const double maxEta = 4.5;
254 const double absEta = fabs(jet.eta());
255
256 const double slopePt = 1./(log(maxPt)-log(minPt));
257 const double slopeEta = 1./(maxEta-minEta);
258
259 const double fixedPt = valPt <= minPt ? minPt+1.e-3 : valPt >= maxPt ? maxPt - 1.e-3 : valPt;
260 const double fixedEta = absEta <= minEta ? minEta+1.e-3 : absEta >= maxEta ? maxEta-1.e-3 : absEta;
261
262 splitFactor = (slopePt*(log(fixedPt)-log(minPt)) + slopeEta*(fixedEta-minEta))/2.;
263 }
264 else if (m_splitNumber == 4 || m_splitNumber == -4)
265 {
266 // Increasing with log(pT) and decreasing with |eta|
267 // See description above, follows similarly
268 // 2*z = (logx-logxmin)/(logxmax-logxmin) + (ymax-y)/(ymax-ymin)
269 const double minPt = histo->GetXaxis()->GetBinLowEdge(1);
270 const double maxPt = histo->GetXaxis()->GetBinLowEdge(histo->GetNbinsX()+1);
271 const double valPt = jet.pt()*m_energyScale;
272
273 const double minEta = 0;
274 const double maxEta = 4.5;
275 const double absEta = fabs(jet.eta());
276
277 const double slopePt = 1./(log(maxPt)-log(minPt));
278 const double slopeEta = 1./(maxEta-minEta);
279
280 const double fixedPt = valPt <= minPt ? minPt+1.e-3 : valPt >= maxPt ? maxPt - 1.e-3 : valPt;
281 const double fixedEta = absEta <= minEta ? minEta+1.e-3 : absEta >= maxEta ? maxEta-1.e-3 : absEta;
282
283 splitFactor = (slopePt*(log(fixedPt)-log(minPt)) + slopeEta*(maxEta-fixedEta))/2.;
284 }
285
286
287 // Now check if this is the functional part or the complementary part
288 if (m_splitNumber < 0)
289 splitFactor = sqrt(1-splitFactor*splitFactor);
290
291 return splitFactor;
292}
293
294
296// //
297// Helper methods //
298// //
300
301bool UncertaintyComponent::getValidBool(const double validity) const
302{
303 if (validity < 1.e-5 && validity > -1.e-5) return false;
304 if (validity < 1+1.e-5 && validity > 1-1.e-5) return true;
305 ATH_MSG_ERROR(Form("Validity value not in expected range: %lf for histogram %s",validity,getValidName().Data()));
306 return false;
307}
308
310{
311 // Check for the simple case (where we want the default four-vector itself)
312 if (massDef == CompMassDef::UNKNOWN || massDef == CompMassDef::FourVecMass)
313 return jet.m();
314
315 // Not the simple case, check for the specified four-vector instead and return it if it is available
317 if (scale.isAvailable(jet))
318 return scale(jet).M();
319
320 // Fall-back on the TA moment as a float if applicable (legacy support)
321 SG::AuxElement::ConstAccessor<float> scaleTAMoment("JetTrackAssistedMassCalibrated");
322 if (massDef == CompMassDef::TAMass && scaleTAMoment.isAvailable(jet))
323 return scaleTAMoment(jet);
324
325 // Fall-back on the calo mass as the 4-vec if applicable (legacy support)
326 if (massDef == CompMassDef::CaloMass)
327 return jet.m();
328
329 // Specified scale is not available, error
330 ATH_MSG_ERROR("Failed to retrieve the " << CompMassDef::enumToString(massDef).Data() << " mass from the jet");
331 return JESUNC_ERROR_CODE;
332}
333
335{
336 // Check for the simple case (where we want the default four-vector itself)
337 if (massDef == CompMassDef::UNKNOWN || massDef == CompMassDef::FourVecMass)
338 return jet.m()/jet.pt();
339
340 // Not the simple case, check for the specified four-vector instead and return it if it is available
342 if (scale.isAvailable(jet))
343 return scale(jet).M()/scale(jet).Pt();
344
345 // Fall-back on the TA moment as a float if applicable (legacy support)
346 SG::AuxElement::ConstAccessor<float> scaleTAMoment("JetTrackAssistedMassCalibrated");
347 if (massDef == CompMassDef::TAMass && scaleTAMoment.isAvailable(jet))
348 return scaleTAMoment(jet)/jet.pt();
349
350 // Fall-back on the calo mass as the 4-vec if applicable (legacy support)
351 if (massDef == CompMassDef::CaloMass)
352 return jet.m()/jet.pt();
353
354 // Specified scale is not available, error
355 ATH_MSG_ERROR("Failed to retrieve the " << CompMassDef::enumToString(massDef).Data() << " mass from the jet");
356 return JESUNC_ERROR_CODE;
357
358}
359
361{
362 // Check for the simple case (where we want the default four-vector itself)
363 if (massDef == CompMassDef::UNKNOWN || massDef == CompMassDef::FourVecMass)
364 return jet.m()/jet.e();
365
366 // Not the simple case, check for the specified four-vector instead and return it if it is available
368 if (scale.isAvailable(jet))
369 return scale(jet).M()/scale(jet).E();
370
371 // Fall-back on the TA moment as a float if applicable (legacy support)
372 SG::AuxElement::ConstAccessor<float> scaleTAMoment("JetTrackAssistedMassCalibrated");
373 if (massDef == CompMassDef::TAMass && scaleTAMoment.isAvailable(jet))
374 return scaleTAMoment(jet)/jet.e();
375
376 // Fall-back on the calo mass as the 4-vec if applicable (legacy support)
377 if (massDef == CompMassDef::CaloMass)
378 return jet.m()/jet.e();
379
380 // Specified scale is not available, error
381 ATH_MSG_ERROR("Failed to retrieve the " << CompMassDef::enumToString(massDef).Data() << " mass from the jet");
382 return JESUNC_ERROR_CODE;
383
384}
385
386} // end jet namespace
387
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
@ Data
Definition BaseObject.h:11
if(febId1==febId2)
static const Attributes_t empty
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:569
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
AsgMessaging(const std::string &name)
Constructor with a name.
std::vector< TString > uncNames
JetFourMomAccessor is an extension of JetAttributeAccessor::AccessorWrapper<xAOD::JetFourMom_t> Acces...
const JetTopology::TypeEnum m_topology
virtual double getUncertainty(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
const Interpolate::TypeEnum m_interpolate
virtual double getMassOverPt(const xAOD::Jet &jet, const CompMassDef::TypeEnum massDef) const
UncertaintyComponent(const ComponentHelper &component, const size_t numHist=1)
virtual TString getName() const
virtual double getUncertaintyImpl(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const =0
virtual TString getValidName() const
virtual double getAbsMass(const xAOD::Jet &jet, const CompMassDef::TypeEnum massDef) const
virtual double getSplitFactor(const xAOD::Jet &jet) const
virtual double getMassOverE(const xAOD::Jet &jet, const CompMassDef::TypeEnum massDef) const
virtual bool getValidity(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
virtual bool getValidBool(const double validity) const
virtual bool getValidUncertainty(double &unc, const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
UncertaintyHistogram * m_validHist
const CompScaleVar::TypeEnum m_scaleVar
virtual bool getValidityImpl(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const =0
virtual StatusCode initialize(TFile *histFile)
UncertaintyHistogram * m_uncHist
TString getJetScaleString(const TypeEnum type)
TString enumToString(const TypeEnum type)
bool operator==(const UncertaintyComponent &componentA, const UncertaintyComponent &componentB)
bool operator<(const UncertaintyComponent &componentA, const UncertaintyComponent &componentB)
Jet_v1 Jet
Definition of the current "jet version".
EventInfo_v1 EventInfo
Definition of the latest event info version.