ATLAS Offline Software
Loading...
Searching...
No Matches
EtaJESCalibStep.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7// EtaJESCalibStep.cxx
8// Implementation file for class EtaJESCalibStep
9// Author: Max Swiatlowski <mswiatlo@cern.ch>
11
14
15#include "TFile.h"
16
18
19EtaJESCalibStep::EtaJESCalibStep(const std::string& name)
20 : asg::AsgTool( name ){ }
21
22
24 ATH_MSG_DEBUG ("Initializing " << name() << " Use spline="<<m_useSpline);
25
26 ATH_MSG_DEBUG("Reading from " << m_jetInScale << " and writing to " << m_jetOutScale);
27
28 if(! m_useSpline){
29 if(!readMCJESFromText()) {
30 ATH_MSG_ERROR("Problem when reading constant file : "<< m_constantFileName);
31 return StatusCode::FAILURE;
32 }
33 } else { // use spline
34 if(!readMCJESFromHists()) {
35 ATH_MSG_ERROR("Problem when reading constant file : "<< m_histoFileName);
36 return StatusCode::FAILURE;
37 }
38 }
39
40 ATH_CHECK( m_vartoolE.retrieve() );
41 ATH_CHECK( m_vartoolEta.retrieve() );
42
43 return StatusCode::SUCCESS;
44}
45
47 ATH_MSG_DEBUG("Calibrating jet collection.");
48
49
51
52 for(xAOD::Jet* jet: jets){
53
54 const xAOD::JetFourMom_t jetStartP4 = jet->getAttribute<xAOD::JetFourMom_t>(m_jetInScale);
55 jet->setJetP4(jetStartP4);
56
57 // Extract the maximum energy, and store in the context
59 double varE {m_vartoolE->getValue(*jet,jc)};
60 double varEta {m_vartoolEta->getValue(*jet,jc)};
61 double Emax = getEmaxJES(varEta);
62
63 // Extract JES from the text handling tool
64 double jesCorrection = getJES(varE, varEta, Emax);
65 int ieta=getEtaBin(varEta);
66
67
68 xAOD::JetFourMom_t calibP4 = jetStartP4 * jesCorrection;
69
70
71 const float etaCorr = calibP4.eta() + getEtaCorr(jesCorrection*varE, varEta) ; //m_textTool_Eta->getValue(*jet, jc);
72 ATH_MSG_DEBUG("eta = "<<etaCorr);
73
74 // Apply the eta correction, use TLV from ROOT to do some math for us
75 TLorentzVector TLVjet;
76 TLVjet.SetPtEtaPhiM( calibP4.P()/cosh(etaCorr), etaCorr, calibP4.phi(), calibP4.M() );
77 calibP4.SetPxPyPzE( TLVjet.Px(), TLVjet.Py(), TLVjet.Pz(), TLVjet.E() );
78 ATH_MSG_DEBUG("JES = "<<jesCorrection << " e="<<varE << " eta="<<jetStartP4.Eta()<< " ieta="<< ieta << " post_pt= ="<< calibP4.Pt() << " etaCorr="<< etaCorr);
79
80 // Set the decorations of this scale
81 jesScaleMomAcc.setAttribute(*jet, calibP4);
82 jet->setJetP4(calibP4);
83 }
84
85 return StatusCode::SUCCESS;
86}
87
88
90{
91 // Open the input file
92 std::string local_path=static_cast<std::string> (m_constantFileName);
93 std::string fileName = PathResolverFindCalibFile(m_constantFileName);
94 if(fileName=="") return false;
95
96 TEnv config(fileName.c_str());
97
98 std::string jetAlgo=static_cast<std::string> (m_jetAlgo);
99
100 std::vector<double> etaBins = VectorizeD(config.GetValue("JES.EtaBins","")," ");
101 if (etaBins.size()==0){ // default binning
102 for (int i=0;i<=90; i++)
103 etaBins.push_back(0.1*i-4.5);
104 }
105
106 ATH_MSG_DEBUG("Number eta bins: " << etaBins.size());
107
108 m_etaBinAxis = new TAxis(etaBins.size()-1,&etaBins[0]);// this is to search for eta bin
109
110 for (uint ieta=0; ieta<etaBins.size()-1; ++ieta)
111 {
112 TString key=Form("JES.%s_Bin%d",jetAlgo.c_str(),ieta);
113 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
114 std::vector<double> params = VectorizeD(config.GetValue(key,"")," ");
115 m_nPar = params.size();
116 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
117 for (uint ipar=0;ipar<m_nPar;++ipar) m_JESFactors[ieta][ipar] = params[ipar];
118 if(m_lowPtExtrap > 0) {
119 //Calculate the slope of the response curve at the minPt for each eta bin
120 //Used in the GetLowPtJES method when Pt < minPt
121 const double *factors = m_JESFactors[ieta];
122 double Ecutoff = m_minPt_JES*cosh(etaBins[ieta]);
123 const double Rcutoff = getLogPolN(factors,Ecutoff);
124 const double Slope = getLogPolNSlope(factors,Ecutoff);
125 if(Slope > Rcutoff/Ecutoff) ATH_MSG_FATAL("Slope of calibration curve at minimum ET is too steep for the JES factors of etabin " << ieta << ", eta = " << etaBins[ieta] );
126
127 m_JES_MinPt_E[ieta] = Ecutoff;
128 m_JES_MinPt_R[ieta] = Rcutoff;
129 m_JES_MinPt_Slopes[ieta] = Slope;
130
131 //Calculate the parameters for a 2nd order polynomial extension to the calibration curve below minimum ET
132 //Used in the GetLowPtJES method when Pt < minPt
133 if(m_lowPtExtrap == 2) {
134 ATH_MSG_ERROR("LowPtJESExtrapolationMethod==2 not supported yet");
135 return false;
136 }
137 }
138
139 key=Form("EtaCorr.%s_Bin%d",jetAlgo.c_str(),ieta);
140 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
141 params = VectorizeD(config.GetValue(key,"")," ");
142 m_nPar = params.size();
143 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
144 for (uint ipar=0;ipar<m_nPar;++ipar) m_etaCorrFactors[ieta][ipar] = params[ipar];
145
146 key=Form("EmaxJES.%s_Bin%d",jetAlgo.c_str(),ieta);
147 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
148 params = VectorizeD(config.GetValue(key,"")," ");
149 m_energyFreezeJES[ieta] = params[0];
150 }
151 return true;
152}
153
154
155
157{
158 // Open the input file
159 std::string local_path=static_cast<std::string> (m_constantFileName);
160 std::string fileName = PathResolverFindCalibFile(m_constantFileName);
161
162 TEnv config(fileName.c_str());
163 if(fileName=="") return false;
164
165 std::string jetAlgo=static_cast<std::string> (m_jetAlgo);
166
167 std::vector<double> etaBins = static_cast<std::vector<double>> (m_etaBins);
168 if (etaBins.size()==0){ // default binning
169 for (int i=0;i<=90; i++)
170 etaBins.push_back(0.1*i-4.5);
171 }
172
173 ATH_MSG_DEBUG("Number eta bins: " << etaBins.size());
174
175 m_etaBinAxis = new TAxis(etaBins.size()-1,&etaBins[0]);// this is to search for eta bin
176
177 std::string calibHistFile = PathResolverFindCalibFile(m_histoFileName);
178 loadSplineHists(calibHistFile, "etaJes");
179
180 //Protections for high order extrapolation methods at low Et (Et < _minPt_JES)
181 if(m_lowPtExtrap != 1) {
182 ATH_MSG_ERROR("Only linear extrapolations are supported for p-splines currently. Please change the config file to reflect this");
183 return false;
184 }
185
186 for (uint ieta=0;ieta<etaBins.size()-1;++ieta) {
187 //Calculate the slope of the response curve at the minPt for each eta bin
188 //Used in the GetLowPtJES method when Pt < minPt
189 double Ecutoff= m_minPt_JES*cosh(etaBins[ieta]);
190 const double Rcutoff = getSplineCorr(ieta, Ecutoff);
191 const double Slope = getSplineSlope(ieta, Ecutoff);
192 if(Slope > Rcutoff/Ecutoff) ATH_MSG_WARNING("Slope of calibration curve at minimum ET is too steep for the JES factors of etabin " << ieta << ", eta = " << etaBins[ieta] );
193
194 m_JES_MinPt_E[ieta] = Ecutoff;
195 m_JES_MinPt_R[ieta] = Rcutoff;
196 m_JES_MinPt_Slopes[ieta] = Slope;
197
198 TString key=Form("EmaxJES.%s_Bin%d",jetAlgo.c_str(),ieta);
199 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
200 std::vector<double> params = VectorizeD(config.GetValue(key,"")," ");
201 m_energyFreezeJES[ieta] = params[0];
202
203 key=Form("EtaCorr.%s_Bin%d",jetAlgo.c_str(),ieta);
204 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
205 params = VectorizeD(config.GetValue(key,"")," ");
206 m_nPar = params.size();
207 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
208 for (uint ipar=0;ipar<m_nPar;++ipar) m_etaCorrFactors[ieta][ipar] = params[ipar];
209
210
211 }
212 return true;
213}
214
215
216double EtaJESCalibStep::getJES(const double X, const double Y, const double Emax) const
217{
218
219 if ( X/cosh(Y) < m_minPt_JES ) { // WARNING !! Won't work if X is actually pT
220 double R = getLowPtJES(X,Y);
221 return 1.0/R;
222 }
223
224 double JES_R;
225 int binEta = getEtaBin(Y);
226 const double *factors = m_JESFactors[binEta];
227
228
229 double E = X;
230 if( m_freezeJESatHighE && (E>Emax) && (Emax!=-1)) {
231 E = Emax;
232 }
233
234 double R = 1.;
235
236 if(m_useSpline){
237 R = getSplineCorr(binEta, E);
238 return 1.0/R;
239 } else {
240 R = getLogPolN(factors,E);
241 }
242
243 JES_R = 1/R;
244
245 return JES_R;
246}
247
248
249double EtaJESCalibStep::getLowPtJES(double E_uncorr, double eta_det) const {
250 int ieta = getEtaBin(eta_det);
251 double R=1;
252 // This correspond to m_lowPtExtrap == 0 in the old EtaJESCorrection tool. Not supporting other cases yet.
253 if (m_lowPtExtrap == 0) {
254 const double *factors = m_JESFactors[ieta];
255 double E = m_minPt_JES*cosh(eta_det);
256 R= getLogPolN(factors,E);
257 } if (m_lowPtExtrap == 1) {
258 double Ecutoff = m_JES_MinPt_E[ieta];
259 double Rcutoff = m_JES_MinPt_R[ieta];
260 double slope = m_JES_MinPt_Slopes[ieta];
261 R = slope*(E_uncorr-Ecutoff)+Rcutoff;
262 }
263 else ATH_MSG_WARNING("Incorrect specification of low Pt JES extrapolation, please check the value of the LowPtJESExtrapolationMethod config flag.");
264
265 return R;
266}
267
268
269double EtaJESCalibStep::getEtaCorr( double X, double Y) const
270{
271 int binEta = getEtaBin(Y);
272 const double *factors = m_etaCorrFactors[binEta];
273
274 if ( X < m_minPt_EtaCorr*cosh(Y) )
275 X = m_minPt_EtaCorr*cosh(Y);
277
278 double eta_corr = getLogPolN(factors,X);
279
280 return -eta_corr;
281}
282
283double EtaJESCalibStep::getEmaxJES(const double Y) const
284{
285 int binEta = getEtaBin(Y);
286 double emaxJES = m_energyFreezeJES[binEta];
287
288 return emaxJES;
289}
290
291double EtaJESCalibStep::getLogPolN(const double *factors, double x) const
292{
293 double y=0;
294 for ( uint i=0; i<m_nPar; ++i )
295 y += factors[i]*TMath::Power(log(x),Int_t(i));
296 return y;
297}
298
299
300
301int EtaJESCalibStep::getEtaBin(double eta_det) const
302{
303 int bin = std::as_const(m_etaBinAxis)->FindBin(eta_det);
304 if (bin<=0) return 0;
305 if (bin>m_etaBinAxis->GetNbins()) return bin-2; // overflow
306 return bin-1;
307}
308
309double EtaJESCalibStep::getLogPolNSlope(const double *factors, double x) const {
310 double y=0;
311 const double inv_x = 1. / x;
312 for ( uint i=0; i<m_nPar; ++i )
313 y += i*factors[i]*TMath::Power(log(x),Int_t(i-1))*inv_x;
314 return y;
315}
316
317
318void EtaJESCalibStep::loadSplineHists(const std::string & fileName, const std::string &etajes_name)
319{
320 std::unique_ptr<TFile> tmpF(TFile::Open( fileName.c_str() ));
321 TList *etajes_l = static_cast<TList*>( tmpF->Get(etajes_name.c_str()));
322
323 m_etajesFactors.resize( etajes_l->GetSize() );
324 if(etajes_l->GetSize() != m_etaBinAxis->GetNbins()+1){
325 ATH_MSG_WARNING("Do not have the correct number of eta bins for " << fileName << "\t" << etajes_name << "\t" << etajes_l->GetSize() );
326 }
327
328 for(int i=0 ; i<m_etaBinAxis->GetNbins(); i++){
329 auto *pTH1 = dynamic_cast<TH1*>(etajes_l->At(i));
330 if (not pTH1) continue;
331 m_etajesFactors[i].reset(pTH1);
332 m_etajesFactors[i]->SetDirectory(nullptr);
333 }
334 tmpF->Close();
335}
336
337
338double EtaJESCalibStep::getSplineSlope(const int ieta, const double minE) const {
339 // Don't want to use interpolation here, so instead just use the values at the bin centers near the cutoff
340 int minBin = m_etajesFactors[ieta]->FindBin(minE);
341
342 double rFirst = m_etajesFactors[ ieta ]->GetBinContent(minBin);
343 double rSecond = m_etajesFactors[ ieta ]->GetBinContent(minBin+1);
344 double binWidth = m_etajesFactors[ ieta ]->GetBinCenter(minBin+1) - m_etajesFactors[ ieta ]->GetBinCenter(minBin);
345 double slope = (rSecond - rFirst) / binWidth;
346
347 return slope;
348}
349
350
351double EtaJESCalibStep::getSplineCorr(const int etaBin, double E) const {
352 double R = m_etajesFactors[ etaBin ]->Interpolate(E);
353 return R;
354}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
VecD VectorizeD(const TString &str, const TString &sep=" ")
unsigned int uint
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
#define y
#define x
virtual StatusCode calibrate(xAOD::JetContainer &) const override
Apply calibration to a jet container.
double m_etaCorrFactors[s_nEtaBins][s_nParMax]
double getSplineCorr(const int etaBin, double E) const
double m_JES_MinPt_Slopes[s_nEtaBins]
Gaudi::Property< float > m_maxE_EtaCorr
Gaudi::Property< std::vector< double > > m_etaBins
Eta bins if other than 0.1 steps from -4.5 to 4.5.
double getLogPolN(const double *factors, double x) const
double getEmaxJES(const double Y) const
return Emax
double getLogPolNSlope(const double *factors, double x) const
Gaudi::Property< std::string > m_jetAlgo
jet collection to be calibrated
Gaudi::Property< std::string > m_histoFileName
ToolHandle< JetHelper::IVarTool > m_vartoolEta
interface for xAOD::jet variable to be defined by user, this must correspond to jet Eta in currect ve...
Gaudi::Property< std::string > m_jetInScale
double getEtaCorr(double X, double Y=0) const
return Eta correction
std::vector< std::unique_ptr< TH1 > > m_etajesFactors
double m_JES_MinPt_R[s_nEtaBins]
double getSplineSlope(const int ieta, const double minE) const
Gaudi::Property< float > m_minPt_JES
double m_JESFactors[s_nEtaBins][s_nParMax]
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< float > m_lowPtExtrap
Gaudi::Property< bool > m_freezeJESatHighE
EtaJESCalibStep(const std::string &name="EtaJESCalibStep")
Gaudi::Property< std::string > m_constantFileName
name of the text file
void loadSplineHists(const std::string &fileName, const std::string &etajes_name)
Gaudi::Property< std::string > m_jetOutScale
double getLowPtJES(double E_uncorr, double eta_det) const
deal with low pt jets
unsigned int m_nPar
double m_JES_MinPt_E[s_nEtaBins]
double getJES(const double X, const double Y=0, const double Emax=-1) const
return MCJES calibration factor
Gaudi::Property< bool > m_useSpline
Gaudi::Property< float > m_minPt_EtaCorr
ToolHandle< JetHelper::IVarTool > m_vartoolE
double m_energyFreezeJES[s_nEtaBins]
int getEtaBin(double eta_det) const
Class JetContext Designed to read AOD information related to the event, N vertices,...
Definition JetContext.h:24
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
void setAttribute(SG::AuxElement &p, const TYPE &v) const
void binWidth(TH1 *h)
Definition listroot.cxx:80
VecD VectorizeD(const TString &str, const TString &sep=" ")
Jet_v1 Jet
Definition of the current "jet version".
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