ATLAS Offline Software
Loading...
Searching...
No Matches
EtaJESCalibStep.h
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.h
8// Header file for class EtaJESCalibStep
9// Author: Max Swiatlowski <mswiatlo@cern.ch>
11#ifndef JETCALIBTOOLS_ETAJESCALIBSTEP_H
12#define JETCALIBTOOLS_ETAJESCALIBSTEP_H 1
13
14#include <string.h>
15
16#include <TString.h>
17#include <TEnv.h>
18
19#include "AsgTools/AsgTool.h"
21#include "AsgTools/ToolHandle.h"
23
25
28
30
32 : public asg::AsgTool,
33 virtual public IJetCalibStep {
34
36
37 public:
38 EtaJESCalibStep(const std::string& name = "EtaJESCalibStep");
39
40 virtual StatusCode initialize() override;
41 virtual StatusCode calibrate(xAOD::JetContainer&) const override;
42
43private:
44 // support functions to extract information from text files
45 double getLogPolN(const double *factors, double x) const;
46 double getLogPolNSlope(const double *factors, double x) const ;
47 int getEtaBin(double eta_det) const;
48 bool readMCJESFromText() ;
49 bool readMCJESFromHists() ;
50 void loadSplineHists(const std::string & fileName, const std::string &etajes_name) ;
51
53 double getJES(const double X, const double Y=0, const double Emax=-1) const;
55 double getEtaCorr( double X, double Y=0) const;
57 double getEmaxJES(const double Y) const;
59 double getLowPtJES(double E_uncorr, double eta_det) const ;
60
61 double getSplineCorr(const int etaBin, double E) const;
62 double getSplineSlope(const int ieta, const double minE) const;
63
64 Gaudi::Property<std::string> m_jetInScale {this, "InScale", "JetPileupScaleMomentum", "Starting jet scale"};
65 Gaudi::Property<std::string> m_jetOutScale {this, "OutScale", "JetEtaJESScaleMomentum", "Ending jet scale"};
66
68 Gaudi::Property< std::string > m_constantFileName { this, "CalibConstantFile", "/afs/cern.ch/work/s/stapiaar/JetDev4/athena/JetToolHelpers/data/file_JES.config", "text file containing constants" };
70 Gaudi::Property< std::string > m_jetAlgo { this, "JetAlgo", "AntiKt4EMPFlow", "jet collection" };
71 Gaudi::Property< float > m_minPt_JES = {this, "MinPtForETAJES",15, "min pT"};
72 Gaudi::Property< bool > m_freezeJESatHighE = {this, "FreezeJEScorrectionatHighE",false, " freeze at high e"};
73 Gaudi::Property< float > m_lowPtExtrap = {this, "LowPtJESExtrapolationMethod", 0, " low pt etrap"};
74 Gaudi::Property< float > m_lowPtMinR = {this, "LowPtJESExtrapolationMinimumResponse", 0.25, " low pt etrap min"};
75
76 Gaudi::Property< float > m_minPt_EtaCorr = {this, "MinPtForEtaCorr" ,8. , ""};
77 Gaudi::Property< float > m_maxE_EtaCorr = {this, "MaxEForEtaCorr" ,2500. , ""};
78
79 Gaudi::Property< std::string > m_histoFileName { this, "HistoFile", "none", "root file containing histos for spline calib" };
80 Gaudi::Property< bool > m_useSpline = {this, "UseSpline",false, " use spline"};
81
83 Gaudi::Property<std::vector<double> > m_etaBins {this, "EtaBins", {} ,""};
84
85 ToolHandle<JetHelper::IVarTool> m_vartoolE{this, "VarToolE", "VarTool", "InputVariable instance E (or pT?)" };
87 ToolHandle<JetHelper::IVarTool> m_vartoolEta{this, "VarToolEta", "VarTool", "InputVariable instance eta (or rapididty?)" };
88
89 // protected:
90 unsigned int m_nPar{}; // number of parameters in config file
91 const static unsigned int s_nEtaBins = 90;
92 const static unsigned int s_nParMax = 9;
99
100
101 TAxis * m_etaBinAxis{};
102
103 std::vector<std::unique_ptr<TH1> > m_etajesFactors;
104
105
106};
107
108#endif
109
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
#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)
static const unsigned int s_nEtaBins
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
static const unsigned int s_nParMax
Gaudi::Property< float > m_lowPtMinR
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
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
JetContainer_v1 JetContainer
Definition of the current "jet container version".