ATLAS Offline Software
Loading...
Searching...
No Matches
CorrelationMatrix.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8#include "TH2D.h"
9
10namespace jet
11{
12
14// //
15// Constructor/destructor/initialization //
16// //
18
20 : asg::AsgMessaging(name)
21 , m_isInit(false)
22 , m_name(name.c_str())
23 , m_numBins(-1)
24 , m_minVal(1e10)
25 , m_maxVal(-1e10)
26 , m_fixedVal1(0)
27 , m_fixedVal2(0)
28 , m_corrMat(nullptr)
29 , m_jets(nullptr)
30 , m_eInfos(nullptr)
31 , m_jet1(nullptr)
32 , m_jet2(nullptr)
33 , m_eInfo(nullptr)
34{
36}
37
38CorrelationMatrix::CorrelationMatrix(const TString& name, const int numBins, const double minVal, const double maxVal, const double fixedVal1, const double fixedVal2)
39 : asg::AsgMessaging(name.Data())
40 , m_isInit(false)
41 , m_name(name)
42 , m_numBins(numBins)
43 , m_minVal(minVal)
44 , m_maxVal(maxVal)
45 , m_fixedVal1(fixedVal1)
46 , m_fixedVal2(fixedVal2)
47 , m_corrMat(nullptr)
48 , m_jets(nullptr)
49 , m_eInfos(nullptr)
50 , m_jet1(nullptr)
51 , m_jet2(nullptr)
52 , m_eInfo(nullptr)
53{
54 ATH_MSG_DEBUG("Creating CorrelationMatrix named " << m_name.Data());
55}
56
57
62
63
65{
66 // Check for initialization and tool status
67 if (checkInitialization(uncTool).isFailure())
68 return StatusCode::FAILURE;
69 m_isInit = true;
70
71 // Determine the center of mass energy for kinematic limits if possible
72 // If not possible, then ignores kinematic limits
73 // pt * cosh(eta) = sqrtS/2 for kinematic limit
74 const double sqrtS = uncTool.getSqrtS();
75 const double kinLimit1 = sqrtS > 0 ? (sqrtS/2.) / cosh(m_fixedVal1) : 1e99;
76 const double kinLimit2 = sqrtS > 0 ? (sqrtS/2.) / cosh(m_fixedVal2) : 1e99;
77
78 std::cout << "Initializing for Pt with fixedVal1, fixedVal2 " << m_fixedVal1 << " " << m_fixedVal2 << std::endl;
79 std::cout << "Corresponding to kinematic limits " << kinLimit1 << " " << kinLimit2 << std::endl;
80
81
82 // Build the correlation matrix
84
85 // Prepare the objects for use
86 if (createStore().isFailure())
87 return StatusCode::FAILURE;
88
89 // Set default values
90 if (setDefaultProperties(uncTool).isFailure())
91 return StatusCode::FAILURE;
92
93 // Fill the correlation matrix
94 for (int binX = 1; binX <= m_corrMat->GetNbinsX(); ++binX)
95 {
96 // Set the jet values, (pT, eta, phi, m)
97 // note that only pT and eta matter right now
98 const double valX = m_corrMat->GetXaxis()->GetBinCenter(binX);
99 m_jet1->setJetP4(xAOD::JetFourMom_t(valX,m_fixedVal1,0,0));
100
101 for (int binY = 1; binY <= m_corrMat->GetNbinsY(); ++binY)
102 {
103 const double valY = m_corrMat->GetYaxis()->GetBinCenter(binY);
104
105 // Check for the kinematic limit (simple case)
106 if (valX > kinLimit1 || valY > kinLimit2 ) {
107// std::cout << "Setting error code in space 2" << std::endl;
108// std::cout << "Values, kin limits are " << valX << " " << kinLimit1 << " " << valY << " " << kinLimit2 << std::endl;
109 m_corrMat->SetBinContent(binX,binY,JESUNC_ERROR_CODE);
110 } else
111 {
112 // Set the jet values, (pT, eta, phi, m)
113 // note that only pT and eta matter right now
114 m_jet2->setJetP4(xAOD::JetFourMom_t(valY,m_fixedVal2,0,0));
115
116 // Now determine the correlation for the two jets
117 m_corrMat->SetBinContent(binX,binY,getCorrelation(uncTool));
118 }
119 }
120 }
121
122 // Clear up the objects we used
123 if (clearStore().isFailure())
124 return StatusCode::FAILURE;
125
126 // Done
127 return StatusCode::SUCCESS;
128}
129
131{
132 // Check for initialization and tool status
133 if (checkInitialization(uncTool).isFailure())
134 return StatusCode::FAILURE;
135 m_isInit = true;
136
137 // Determine the center of mass energy for kinematic limits if possible
138 // If not possible, then ignores kinematic limits
139 // pt * cosh(eta) = sqrtS/2 for kinematic limit
140 const double sqrtS = uncTool.getSqrtS();
141 const double kinLimit1 = sqrtS > 0 ? acosh((sqrtS/2.) / m_fixedVal1) : 1e99;
142 const double kinLimit2 = sqrtS > 0 ? acosh((sqrtS/2.) / m_fixedVal2) : 1e99;
143
144 // Build the correlation matrix
146
147 // Prepare the objects for use
148 if (createStore().isFailure())
149 return StatusCode::FAILURE;
150
151 // Set default values
152 if (setDefaultProperties(uncTool).isFailure())
153 return StatusCode::FAILURE;
154
155 // Fill the correlation matrix
156 for (int binX = 1; binX <= m_corrMat->GetNbinsX(); ++binX)
157 {
158 // Set the jet values, (pT, eta, phi, m)
159 // note that only pT and eta matter right now
160 const double valX = m_corrMat->GetXaxis()->GetBinCenter(binX);
161 m_jet1->setJetP4(xAOD::JetFourMom_t(m_fixedVal1,valX,0,0));
162
163 for (int binY = 1; binY <= m_corrMat->GetNbinsY(); ++binY)
164 {
165 const double valY = m_corrMat->GetYaxis()->GetBinCenter(binY);
166
167 // Check for the kinematic limit (simple case)
168 if (valX > kinLimit1 || valY > kinLimit2) {
169 m_corrMat->SetBinContent(binX,binY,JESUNC_ERROR_CODE);
170 } else
171 {
172 // Set the jet values, (pT, eta, phi, m)
173 // note that only pT and eta matter right now
174 m_jet2->setJetP4(xAOD::JetFourMom_t(m_fixedVal2,valY,0,0));
175
176 // Now determine the correlation for the two jets
177 m_corrMat->SetBinContent(binX,binY,getCorrelation(uncTool));
178 }
179 }
180 }
181
182 // Clear up the objects we used
183 if (clearStore().isFailure())
184 return StatusCode::FAILURE;
185
186 // Done
187 return StatusCode::SUCCESS;
188}
189
190
192// //
193// Helper methods for building the matrices //
194// //
196
198{
199 // Prevent double-initialization
200 if (m_isInit)
201 {
202 ATH_MSG_ERROR("CorrelationMatrix is already initialized: " << m_name.Data());
203 return StatusCode::FAILURE;
204 }
205
206 // Warn about any non-four-vector terms which will be ignored
207 bool scalesFourVec = false;
208 for (size_t iComp = 0; iComp < uncTool.getNumComponents(); ++iComp)
209 if (!uncTool.getComponentScalesFourVec(iComp))
210 ATH_MSG_WARNING("CorrelationMatrix will ignore component which does not scale the jet four-vector: " << uncTool.getComponentName(iComp));
211 else
212 scalesFourVec = true;
213
214 // If no components scale the four-vector, end now
215 if (!scalesFourVec)
216 {
217 ATH_MSG_ERROR("No components scale the four-vector for the CorrelationMatrix: " << m_name.Data());
218 return StatusCode::FAILURE;
219 }
220
221 return StatusCode::SUCCESS;
222}
223
225{
226 // Build a jet container and a pair of jets for us to manipulate later
228 m_jets->setStore(new xAOD::JetAuxContainer());
229 m_jets->push_back(new xAOD::Jet());
230 m_jets->push_back(new xAOD::Jet());
231 m_jet1 = m_jets->at(0);
232 m_jet2 = m_jets->at(1);
233
234 // Build an EventInfo object for us to manipulate later
236 m_eInfos->setStore(new xAOD::EventInfoAuxContainer());
237 m_eInfos->push_back(new xAOD::EventInfo());
238 m_eInfo = m_eInfos->at(0);
239
240 return StatusCode::SUCCESS;
241}
242
244{
245 m_jet1 = nullptr;
246 m_jet2 = nullptr;
247 m_eInfo = nullptr;
248
249 //for (size_t iJet = 0; iJet < m_jets->size(); ++iJet)
250 // JESUNC_SAFE_DELETE((*m_jets)[iJet]);
251 m_jets->clear();
253
254 //for (size_t iInfo = 0; iInfo < m_eInfos->size(); ++iInfo)
255 // JESUNC_SAFE_DELETE((*m_eInfos)[iInfo]);
256 m_eInfos->clear();
258
259 return StatusCode::SUCCESS;
260}
261
263{
264
265 // TODO make this part of a common file
266 static const SG::AuxElement::Accessor<int> Nsegments("GhostMuonSegmentCount");
267 static const SG::AuxElement::Accessor<char> IsBjet("IsBjet");
268 static const SG::AuxElement::Accessor<float> mu("averageInteractionsPerCrossing");
269 static const SG::AuxElement::Accessor<float> NPV("NPV");
270
271 // 25 segments is about average for jets receiving a correction
272 // This is modulated by the probability of punchthrough
273 // TODO add punch-through Nsegments/etc dependence on probability
274 Nsegments(*m_jet1) = 0; //25;
275 Nsegments(*m_jet2) = 0; //25;
276 IsBjet(*m_jet1) = false;
277 IsBjet(*m_jet2) = false;
278
279 float sigmaMu = 0;
280 float sigmaNPV = 0;
281
282 const TString config = uncTool.getConfigFile();
283 if (config.Contains("_2011/"))
284 {
285 // Dag, night of Febuary 4/5, 2013
286 sigmaMu = 3.0;
287 sigmaNPV = 3.0;
288 }
289 else if (config.Contains("_2012/") || config.Contains("_2015/Prerec"))
290 {
291 // Craig Sawyer, Jan 22 2013
292 sigmaMu = 5.593*1.11;
293 sigmaNPV = 3.672;
294 }
295 else if (config.Contains("_2015"))
296 {
297 // Kate, Jan 31 2016
298 sigmaMu = 1.9;
299 sigmaNPV = 2.9;
300 }
301 else if (config.Contains("_2016"))
302 {
303 // Kate, Nov 2016
304 // via Eric Corrigan's pileup studies
305 // Scaling term taken from
306 // https://indico.cern.ch/event/437993/contributions/1925644/attachments/1138739/1630981/spagan_MuRescaling.pdf
307 sigmaMu = 5.35;
308 sigmaNPV = 3.49;
309 }
310 else
311 {
312 ATH_MSG_ERROR("Unexpected year for setPileupShiftsForYear in config: " << config.Data());
313 return StatusCode::FAILURE;
314 }
315
316 mu(*m_eInfo) = uncTool.getRefMu()+sigmaMu;
317 NPV(*m_eInfo) = uncTool.getRefNPV()+sigmaNPV;
318
319 return StatusCode::SUCCESS;
320}
321
322TH2D* CorrelationMatrix::buildMatrix(const std::vector<double>& bins) const
323{
324 TH2D* matrix = new TH2D(m_name,m_name,bins.size()-1,&bins[0],bins.size()-1,&bins[0]);
325 matrix->GetZaxis()->SetRangeUser(-1,1);
326 return matrix;
327}
328
330{
331 return getCovariance(uncTool,m_jet1,m_jet2)/sqrt(getCovariance(uncTool,m_jet1,m_jet1)*getCovariance(uncTool,m_jet2,m_jet2));
332}
333
334double CorrelationMatrix::getCovariance(const JetUncertaintiesTool& uncTool, const xAOD::Jet* jet1, const xAOD::Jet* jet2) const
335{
336 double covariance = 0;
337 for (size_t iComp = 0; iComp < uncTool.getNumComponents(); ++iComp)
338 covariance += uncTool.getUncertainty(iComp,*jet1,*m_eInfo)*uncTool.getUncertainty(iComp,*jet2,*m_eInfo);
339
340 return covariance;
341}
342
343} // end jet namespace
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
@ Data
Definition BaseObject.h:11
static const std::vector< std::string > bins
virtual bool getComponentScalesFourVec(const size_t index) const
virtual std::string getConfigFile() const
virtual float getRefNPV() const
virtual double getUncertainty(size_t index, const xAOD::Jet &jet) const
virtual float getRefMu() const
virtual float getSqrtS() const
virtual size_t getNumComponents() const
virtual std::string getComponentName(const size_t index) const
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
AsgMessaging(const std::string &name)
Constructor with a name.
double getCovariance(const JetUncertaintiesTool &uncTool, const xAOD::Jet *jet1, const xAOD::Jet *jet2) const
virtual StatusCode initializeForEta(const JetUncertaintiesTool &uncTool)
TH2D * buildMatrix(const std::vector< double > &bins) const
CorrelationMatrix(const TString &name, const int numBins, const double minVal, const double maxVal, const double fixedVal1, const double fixedVal2)
StatusCode setDefaultProperties(const JetUncertaintiesTool &uncTool)
double getCorrelation(const JetUncertaintiesTool &uncTool) const
virtual StatusCode initializeForPt(const JetUncertaintiesTool &uncTool)
xAOD::EventInfoContainer * m_eInfos
StatusCode checkInitialization(const JetUncertaintiesTool &uncTool) const
xAOD::EventInfo * m_eInfo
xAOD::JetContainer * m_jets
std::vector< double > getLogBins(const size_t numBins, const double minVal, const double maxVal)
std::vector< double > getUniformBins(const size_t numBins, const double minVal, const double maxVal)
Jet_v1 Jet
Definition of the current "jet version".
EventInfoContainer_v1 EventInfoContainer
Define the latest version of the container.
EventInfo_v1 EventInfo
Definition of the latest event info version.
EventInfoAuxContainer_v1 EventInfoAuxContainer
Define the latest version of the auxiliary container.
JetAuxContainer_v1 JetAuxContainer
Definition of the current jet auxiliary container.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17