ATLAS Offline Software
Loading...
Searching...
No Matches
EtaJESCalibStep.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2026 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 = 14000;
63 Emax = getEmaxJES(varEta);
64 }
65
66 // Extract JES from the text handling tool
67 double jesCorrection = getJES(varE, varEta, Emax);
68 int ieta=getEtaBin(varEta);
69
70
71 xAOD::JetFourMom_t calibP4 = jetStartP4 * jesCorrection;
72
73
74 const float etaCorr = calibP4.eta() + getEtaCorr(jesCorrection*varE, varEta) ; //m_textTool_Eta->getValue(*jet, jc);
75 ATH_MSG_DEBUG("eta = "<<etaCorr);
76
77 // Apply the eta correction, use TLV from ROOT to do some math for us
78 TLorentzVector TLVjet;
79 TLVjet.SetPtEtaPhiM( calibP4.P()/cosh(etaCorr), etaCorr, calibP4.phi(), calibP4.M() );
80 calibP4.SetPxPyPzE( TLVjet.Px(), TLVjet.Py(), TLVjet.Pz(), TLVjet.E() );
81 ATH_MSG_DEBUG("JES = "<<jesCorrection << " e="<<varE << " eta="<<jetStartP4.Eta()<< " ieta="<< ieta << " post_pt= ="<< calibP4.Pt() << " etaCorr="<< etaCorr);
82
83 // Set the decorations of this scale
84 jesScaleMomAcc.setAttribute(*jet, calibP4);
85 jet->setJetP4(calibP4);
86 }
87
88 return StatusCode::SUCCESS;
89}
90
91
93{
94 // Open the input file
95 std::string local_path=static_cast<std::string> (m_constantFileName);
96 std::string fileName = PathResolverFindCalibFile(m_constantFileName);
97 if(fileName=="") return false;
98
99 TEnv config(fileName.c_str());
100
101 std::string jetAlgo=static_cast<std::string> (m_jetAlgo);
102
103 std::vector<double> etaBins = static_cast<std::vector<double>> (m_etaBins);
104 if (etaBins.size()==0){ // default binning
105 for (int i=0;i<=90; i++)
106 etaBins.push_back(0.1*i-4.5);
107 }
108
109 ATH_MSG_DEBUG("Number eta bins: " << etaBins.size());
110
111 m_etaBinAxis = new TAxis(etaBins.size()-1,&etaBins[0]);// this is to search for eta bin
112
113 for (uint ieta=0; ieta<etaBins.size()-1; ++ieta)
114 {
115 TString key=Form("JES.%s_Bin%d",jetAlgo.c_str(),ieta);
116 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
117 std::vector<double> params = VectorizeD(config.GetValue(key,"")," ");
118 m_nPar = params.size();
119 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
120 for (uint ipar=0;ipar<m_nPar;++ipar) m_JESFactors[ieta][ipar] = params[ipar];
121 if(m_lowPtExtrap > 0) {
122 //Calculate the slope of the response curve at the minPt for each eta bin
123 //Used in the GetLowPtJES method when Pt < minPt
124 const double *factors = m_JESFactors[ieta];
125 double Ecutoff = m_minPt_JES*cosh(etaBins[ieta]);
126 if(m_useSecondaryminPt_JES && std::abs(etaBins[ieta]) >= m_etaSecondaryminPt_JES) Ecutoff = m_secondaryminPt_JES*cosh(etaBins[ieta]);
127 const double Rcutoff = getLogPolN(factors,Ecutoff);
128 const double Slope = getLogPolNSlope(factors,Ecutoff);
129 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] );
130
131 m_JES_MinPt_E[ieta] = Ecutoff;
132 m_JES_MinPt_R[ieta] = Rcutoff;
133 m_JES_MinPt_Slopes[ieta] = Slope;
134
135 //Calculate the parameters for a 2nd order polynomial extension to the calibration curve below minimum ET
136 //Used in the GetLowPtJES method when Pt < minPt
137 if(m_lowPtExtrap == 2) {
138 ATH_MSG_ERROR("LowPtJESExtrapolationMethod==2 not supported yet");
139 return false;
140 }
141 }
142
143 key=Form("EtaCorr.%s_Bin%d",jetAlgo.c_str(),ieta);
144 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
145 params = VectorizeD(config.GetValue(key,"")," ");
146 m_nPar = params.size();
147 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
148 for (uint ipar=0;ipar<m_nPar;++ipar) m_etaCorrFactors[ieta][ipar] = params[ipar];
149
151 key=Form("EmaxJES.%s_Bin%d",jetAlgo.c_str(),ieta);
152 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
153 params = VectorizeD(config.GetValue(key,"")," ");
154 m_energyFreezeJES[ieta] = params[0];
155 }
156 }
157 return true;
158}
159
160
161
163{
164 // Open the input file
165 std::string local_path=static_cast<std::string> (m_constantFileName);
166 std::string fileName = PathResolverFindCalibFile(m_constantFileName);
167
168 TEnv config(fileName.c_str());
169 if(fileName=="") return false;
170
171 std::string jetAlgo=static_cast<std::string> (m_jetAlgo);
172
173 std::vector<double> etaBins = static_cast<std::vector<double>> (m_etaBins);
174 if (etaBins.size()==0){ // default binning
175 for (int i=0;i<=90; i++)
176 etaBins.push_back(0.1*i-4.5);
177 }
178
179 ATH_MSG_DEBUG("Number eta bins: " << etaBins.size());
180
181 m_etaBinAxis = new TAxis(etaBins.size()-1,&etaBins[0]);// this is to search for eta bin
182
183 std::string calibHistFile = PathResolverFindCalibFile(m_histoFileName);
184 loadSplineHists(calibHistFile, "etaJes");
185
186 //Protections for high order extrapolation methods at low Et (Et < _minPt_JES)
187 if(m_lowPtExtrap != 1) {
188 ATH_MSG_ERROR("Only linear extrapolations are supported for p-splines currently. Please change the config file to reflect this");
189 return false;
190 }
191
192 for (uint ieta=0;ieta<etaBins.size()-1;++ieta) {
193 //Calculate the slope of the response curve at the minPt for each eta bin
194 //Used in the GetLowPtJES method when Pt < minPt
195 double Ecutoff= m_minPt_JES*cosh(etaBins[ieta]);
196 if(m_useSecondaryminPt_JES && std::abs(etaBins[ieta]) >= m_etaSecondaryminPt_JES) Ecutoff = m_secondaryminPt_JES*cosh(etaBins[ieta]);
197 const double Rcutoff = getSplineCorr(ieta, Ecutoff);
198 const double Slope = getSplineSlope(ieta, Ecutoff);
199 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] );
200
201 m_JES_MinPt_E[ieta] = Ecutoff;
202 m_JES_MinPt_R[ieta] = Rcutoff;
203 m_JES_MinPt_Slopes[ieta] = Slope;
204
205 TString key=Form("EtaCorr.%s_Bin%d",jetAlgo.c_str(),ieta);
206 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
207 std::vector<double> params = VectorizeD(config.GetValue(key,"")," ");
208 m_nPar = params.size();
209 ATH_MSG_VERBOSE("Number of parameters: " << m_nPar);
210 for (uint ipar=0;ipar<m_nPar;++ipar) m_etaCorrFactors[ieta][ipar] = params[ipar];
211
213 key=Form("EmaxJES.%s_Bin%d",jetAlgo.c_str(),ieta);
214 ATH_MSG_VERBOSE("reading: " << key << " = "<< config.GetValue(key,""));
215 params = VectorizeD(config.GetValue(key,"")," ");
216 m_energyFreezeJES[ieta] = params[0];
217 }
218
219 }
220 return true;
221}
222
223
224double EtaJESCalibStep::getJES(const double X, const double Y, const double Emax) const
225{
226
228 if ( X/cosh(Y) < m_minPt_JES ) { // WARNING !! Won't work if X is actually pT
229 double R = getLowPtJES(X,Y);
230 return 1.0/R;
231 }
232 }
233 else{
234 if(std::abs(Y) < m_etaSecondaryminPt_JES && X/cosh(Y) < m_minPt_JES){
235 double R = getLowPtJES(X,Y);
236 return 1.0/R;
237 }
238 if(std::abs(Y) >= m_etaSecondaryminPt_JES && X/cosh(Y) < m_secondaryminPt_JES){
239 double R = getLowPtJES(X,Y);
240 return 1.0/R;
241 }
242 }
243
244 double JES_R;
245 int binEta = getEtaBin(Y);
246 const double *factors = m_JESFactors[binEta];
247
248
249 double E = X;
250 if( m_freezeJESatHighE && (E>Emax) && (Emax!=-1)) {
251 E = Emax;
252 }
253
254 double R = 1.;
255
256 if(m_useSpline){
257 R = getSplineCorr(binEta, E);
258 return 1.0/R;
259 } else {
260 R = getLogPolN(factors,E);
261 }
262
263 JES_R = 1/R;
264
265 return JES_R;
266}
267
268
269double EtaJESCalibStep::getLowPtJES(double E_uncorr, double eta_det) const {
270 int ieta = getEtaBin(eta_det);
271 double R=1;
272 // This correspond to m_lowPtExtrap == 0 in the old EtaJESCorrection tool. Not supporting other cases yet.
273 if (m_lowPtExtrap == 0) {
274 const double *factors = m_JESFactors[ieta];
275 double E = m_minPt_JES*cosh(eta_det);
276 R= getLogPolN(factors,E);
277 } else if (m_lowPtExtrap == 1) {
278 double Ecutoff = m_JES_MinPt_E[ieta];
279 double Rcutoff = m_JES_MinPt_R[ieta];
280 double slope = m_JES_MinPt_Slopes[ieta];
281 R = slope*(E_uncorr-Ecutoff)+Rcutoff;
282 }
283 else ATH_MSG_WARNING("Incorrect specification of low Pt JES extrapolation, please check the value of the LowPtJESExtrapolationMethod config flag.");
284
285 return R;
286}
287
288
289double EtaJESCalibStep::getEtaCorr( double X, double Y) const
290{
291 int binEta = getEtaBin(Y);
292 const double *factors = m_etaCorrFactors[binEta];
293
294 if ( X < m_minPt_EtaCorr*cosh(Y) )
295 X = m_minPt_EtaCorr*cosh(Y);
297
298 double eta_corr = getLogPolN(factors,X);
299
300 return -eta_corr;
301}
302
303double EtaJESCalibStep::getEmaxJES(const double Y) const
304{
305 int binEta = getEtaBin(Y);
306 double emaxJES = m_energyFreezeJES[binEta];
307
308 return emaxJES;
309}
310
311double EtaJESCalibStep::getLogPolN(const double *factors, double x) const
312{
313 double y=0;
314 for ( uint i=0; i<m_nPar; ++i )
315 y += factors[i]*TMath::Power(log(x),Int_t(i));
316 return y;
317}
318
319
320
321int EtaJESCalibStep::getEtaBin(double eta_det) const
322{
323 int bin = std::as_const(m_etaBinAxis)->FindBin(eta_det);
324 if (bin<=0) return 0;
325 if (bin>m_etaBinAxis->GetNbins()) return bin-2; // overflow
326 return bin-1;
327}
328
329double EtaJESCalibStep::getLogPolNSlope(const double *factors, double x) const {
330 double y=0;
331 const double inv_x = 1. / x;
332 for ( uint i=0; i<m_nPar; ++i )
333 y += i*factors[i]*TMath::Power(log(x),Int_t(i-1))*inv_x;
334 return y;
335}
336
337
338void EtaJESCalibStep::loadSplineHists(const std::string & fileName, const std::string &etajes_name)
339{
340 std::unique_ptr<TFile> tmpF(TFile::Open( fileName.c_str() ));
341 TList *etajes_l = static_cast<TList*>( tmpF->Get(etajes_name.c_str()));
342
343 m_etajesFactors.resize( etajes_l->GetSize() );
344 if(etajes_l->GetSize() != m_etaBinAxis->GetNbins()+1){
345 ATH_MSG_WARNING("Do not have the correct number of eta bins for " << fileName << "\t" << etajes_name << "\t" << etajes_l->GetSize() );
346 }
347
348 for(int i=0 ; i<m_etaBinAxis->GetNbins(); i++){
349 auto *pTH1 = dynamic_cast<TH1*>(etajes_l->At(i));
350 if (not pTH1) continue;
351 m_etajesFactors[i].reset(pTH1);
352 m_etajesFactors[i]->SetDirectory(nullptr);
353 }
354 tmpF->Close();
355}
356
357
358double EtaJESCalibStep::getSplineSlope(const int ieta, const double minE) const {
359 // Don't want to use interpolation here, so instead just use the values at the bin centers near the cutoff
360 int minBin = m_etajesFactors[ieta]->FindBin(minE);
361
362 double rFirst = m_etajesFactors[ ieta ]->GetBinContent(minBin);
363 double rSecond = m_etajesFactors[ ieta ]->GetBinContent(minBin+1);
364 double binWidth = m_etajesFactors[ ieta ]->GetBinCenter(minBin+1) - m_etajesFactors[ ieta ]->GetBinCenter(minBin);
365 double slope = (rSecond - rFirst) / binWidth;
366
367 return slope;
368}
369
370
371double EtaJESCalibStep::getSplineCorr(const int etaBin, double E) const {
372 double R = m_etajesFactors[ etaBin ]->Interpolate(E);
373 return R;
374}
#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< double > m_etaSecondaryminPt_JES
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< 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< float > m_secondaryminPt_JES
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
Gaudi::Property< bool > m_useSecondaryminPt_JES
double m_energyFreezeJES[s_nEtaBins]
Gaudi::Property< int > m_lowPtExtrap
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