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