ATLAS Offline Software
Loading...
Searching...
No Matches
GlobalNNCalibration.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6#ifndef JetCalibTools_GlobalNNCalibration_H
7#define JetCalibTools_GlobalNNCalibration_H
8
9#include "TH1.h"
10#include "TString.h"
11#include "lwtnn/LightweightGraph.hh"
12
13// Local includes
15#include <string>
16#include <vector>
17#include <map>
18#include <memory>
19
20class JetEventInfo;
21class TEnv;
22
23
30
31
32
33
34
36
37 public:
38 // Constructor/destructor/init
43
48 GlobalNNCalibration(const std::string& name);
49
58 GlobalNNCalibration(const std::string& name, TEnv * config, TString jetAlgo, const TString& calibAreaTag, bool dev);
59
63 virtual ~GlobalNNCalibration() = default;
64
69 virtual StatusCode initialize() override;
70
71
72 protected:
73 // @brief Calibrates the jet, and decorates it with the calibration using the name "JetGNNCScaleMomentum"
74 // @param jet_reco The jet
75 // @param jetEventInfo A set of information about the event and jet
76 virtual StatusCode calibrate(xAOD::Jet& jet, JetEventInfo&) const override;
77
78 private:
79
85 double getJetChargedFraction(const xAOD::Jet& jet_reco, int PVindex) const;
86
91 double getJetDetEta(const xAOD::Jet& jet_reco) const;
92
98 int getJetNtrk1000(const xAOD::Jet& jet_reco, int PVindex) const;
99
105 double getJetWtrk1000(const xAOD::Jet& jet_reco, int PVindex) const;
106
111 double getJESPt(const xAOD::Jet& jet_reco) const;
112
118 std::map<std::string,double> getJetFeatures(const xAOD::Jet& jet_reco, JetEventInfo& jetEventInfo) const;
119
125 double getSplineSlope(const int ieta, const double minPt) const;
126
132 void loadSplineHists(const TString & fileName, const std::string &etajes_name = "etaJes");
133
139 double getSplineCorr(const int etaBin, double E) const;
140
141
147 int getEtaBin(const xAOD::Jet& jet_reco, const std::vector<double>& etaBins) const;
148
149
150 std::vector<std::unique_ptr<lwt::LightweightGraph> > m_lwnns;
151 std::vector<std::unique_ptr<TH1> > m_ptCorrFactors;
152 std::vector<double> m_nnEtaBins;
153 std::vector<double> m_closureEtaBins;
154 std::vector<TString> m_NNInputs;
155
156
157 TEnv * m_config{};
158 TString m_jetAlgo;
159 std::string m_calibAreaTag;
160 bool m_dev{};
161 bool m_doSplineCorr{true};
165
166 std::vector<double> m_JPtS_MinPt_Slopes;
167 std::vector<double> m_JPtS_MinPt_Pt;
168 std::vector<double> m_JPtS_MinPt_R;
169
170
171}; // Class GlobalNNCalibration
172
173
174#endif
175
176
std::vector< std::unique_ptr< TH1 > > m_ptCorrFactors
double getJetDetEta(const xAOD::Jet &jet_reco) const
Returns the detector eta of the jet.
std::vector< double > m_closureEtaBins
std::vector< double > m_JPtS_MinPt_Slopes
std::vector< double > m_JPtS_MinPt_Pt
double getJetWtrk1000(const xAOD::Jet &jet_reco, int PVindex) const
Returns the jet width.
std::vector< TString > m_NNInputs
GlobalNNCalibration()
The constructor.
int getJetNtrk1000(const xAOD::Jet &jet_reco, int PVindex) const
Returns the number of tracks with pT > 1 GeV associated to the jet.
virtual StatusCode initialize() override
Returns the charged fraction of a jet.
void loadSplineHists(const TString &fileName, const std::string &etajes_name="etaJes")
Reads the spline histograms from the file given in the config, and stores them in m_ptCorrFactors.
virtual ~GlobalNNCalibration()=default
The destructor.
double getJetChargedFraction(const xAOD::Jet &jet_reco, int PVindex) const
Returns the charged fraction of a jet.
double getJESPt(const xAOD::Jet &jet_reco) const
Returns the jet pT after the MCJES calibration.
virtual StatusCode calibrate(xAOD::Jet &jet, JetEventInfo &) const override
double getSplineCorr(const int etaBin, double E) const
Returns the correction from spline histogram, which should be applied after the NN correction.
std::map< std::string, double > getJetFeatures(const xAOD::Jet &jet_reco, JetEventInfo &jetEventInfo) const
Returns a map of possible inputs to the NN, and their corresponding values for this jet.
std::vector< double > m_nnEtaBins
int getEtaBin(const xAOD::Jet &jet_reco, const std::vector< double > &etaBins) const
Returns the eta bin, as determined by a list of bin edges.
double getSplineSlope(const int ieta, const double minPt) const
Gets the slope of the spline histogram for a given eta bin, for extrapolation of the calibration.
std::vector< double > m_JPtS_MinPt_R
std::vector< std::unique_ptr< lwt::LightweightGraph > > m_lwnns
JetCalibrationStep(const char *name="JetCalibrationStep")
Jet_v1 Jet
Definition of the current "jet version".