CMS 3D CMS Logo

JetResolutionDemo.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <iomanip>
3 
6 
13 
16 
18 
19 #include <TH2.h>
20 
21 //
22 // class declaration
23 //
24 
25 class JetResolutionDemo : public edm::one::EDAnalyzer<edm::one::SharedResources> {
26 public:
27  explicit JetResolutionDemo(const edm::ParameterSet&);
28  ~JetResolutionDemo() override;
29 
30 private:
31  void analyze(const edm::Event&, const edm::EventSetup&) override;
32  void endJob() override;
33 
36 
37  bool m_debug = false;
38  bool m_use_conddb = false;
39 
46 
47  TH2* m_res_vs_eta = nullptr;
48 };
49 //
50 //----------- Class Implementation ------------------------------------------
51 //
52 //---------------------------------------------------------------------------
54  m_jets_token = consumes<std::vector<pat::Jet>>(iConfig.getParameter<edm::InputTag>("jets"));
55  m_rho_token = consumes<double>(iConfig.getParameter<edm::InputTag>("rho"));
56  m_debug = iConfig.getUntrackedParameter<bool>("debug", false);
57  m_use_conddb = iConfig.getUntrackedParameter<bool>("useCondDB", false);
58 
59  if (m_use_conddb) {
60  m_payload = iConfig.getParameter<std::string>("payload");
63  } else {
64  m_resolutions_file = iConfig.getParameter<edm::FileInPath>("resolutionsFile").fullPath();
65  m_scale_factors_file = iConfig.getParameter<edm::FileInPath>("scaleFactorsFile").fullPath();
66  }
67  usesResource(TFileService::kSharedResource);
68 }
69 //---------------------------------------------------------------------------
71 
72 //---------------------------------------------------------------------------
75  iEvent.getByToken(m_jets_token, jets);
76 
78  iEvent.getByToken(m_rho_token, rho);
79 
80  // Access jet resolution and scale factor from the condition database
81  // or from text files
84 
85  // Two differents way to create a class instance
86  if (m_use_conddb) {
87  // First way, using the get() static method
90  } else {
91  // Second way, using the constructor
94  }
95 
96  if (m_debug) {
97  // Dump resolution
98  resolution.dump();
99  }
100 
101  if (!m_res_vs_eta) {
102  // Advanced usage. Create the histogram by retriving the eta binning directly from the JetResolutionObject. This suppose you need that the parametrization only depends on eta.
103 
104  // Get the list of bins of this object
105  const std::vector<JME::Binning>& bins = resolution.getResolutionObject()->getDefinition().getBins();
106 
107  // Check that the first bin is eta
108  if ((!bins.empty()) && (bins[0] == JME::Binning::JetEta)) {
109  const std::vector<JME::JetResolutionObject::Record> records = resolution.getResolutionObject()->getRecords();
110  // Get all records from the object. Each record correspond to a different binning and different parameters
111 
112  std::vector<float> etas;
113  for (const auto& record : records) {
114  if (etas.empty()) {
115  etas.push_back(record.getBinsRange()[0].min);
116  etas.push_back(record.getBinsRange()[0].max);
117  } else {
118  etas.push_back(record.getBinsRange()[0].max);
119  }
120  }
121 
122  std::vector<float> res;
123  for (std::size_t i = 0; i < 40; i++) {
124  res.push_back(i * 0.005);
125  }
126 
127  m_res_vs_eta = fs->make<TH2F>("res_vs_eta", "res_vs_eta", etas.size() - 1, &etas[0], res.size() - 1, &res[0]);
128  }
129  }
130 
131  for (const auto& jet : *jets) {
132  if (m_debug) {
133  std::cout << "New jet; pt=" << jet.pt() << " eta=" << jet.eta() << " phi=" << jet.phi()
134  << " e=" << jet.energy() << " rho=" << *rho << std::endl;
135  }
136 
137  // Get resolution for this jet
138  // The method to use is getResolution(const JME::JetParameters& parameters) from a JetResolution object
139 
140  // You *must* now in advance which variables are needed for getting the resolution. For the moment, only pt and eta are needed, but this will
141  // probably change in the futur when PU dependency is added. Please keep an eye on the twiki page
142  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution
143  // to stay up-to-date. All currently supported parameters (ie, 'set' functions) are available here:
144  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution#List_of_supported_parameters
145 
146  // Three way to create a JetParameters object
147 
148  // First, using 'set' functions
149  JME::JetParameters parameters_1;
150  parameters_1.setJetPt(jet.pt());
151  parameters_1.setJetEta(jet.eta());
152  parameters_1.setRho(*rho);
153 
154  // You can also chain calls
155 
156  JME::JetParameters parameters_2;
157  parameters_2.setJetPt(jet.pt()).setJetEta(jet.eta()).setRho(*rho);
158 
159  // Second, using the set() function
160  JME::JetParameters parameters_3;
161  parameters_3.set(JME::Binning::JetPt, jet.pt());
162  parameters_3.set({JME::Binning::JetEta, jet.eta()});
163  parameters_3.set({JME::Binning::Rho, *rho});
164 
165  // Or
166 
167  JME::JetParameters parameters_4;
168  parameters_4.set(JME::Binning::JetPt, jet.pt()).set(JME::Binning::JetEta, jet.eta()).set(JME::Binning::Rho, *rho);
169 
170  // Third, using a initializer_list
171  JME::JetParameters parameters_5 = {
173 
174  // Now, get the resolution
175 
176  float r = resolution.getResolution(parameters_1);
177  if (m_debug) {
178  std::cout << "Resolution with parameters_1: " << r << std::endl;
179  std::cout << "Resolution with parameters_2: " << resolution.getResolution(parameters_2) << std::endl;
180  std::cout << "Resolution with parameters_3: " << resolution.getResolution(parameters_3) << std::endl;
181  std::cout << "Resolution with parameters_4: " << resolution.getResolution(parameters_4) << std::endl;
182  std::cout << "Resolution with parameters_5: " << resolution.getResolution(parameters_5) << std::endl;
183 
184  // You can also use a shortcut to get the resolution
185  float r2 = resolution.getResolution(
187  std::cout << "Resolution using shortcut : " << r2 << std::endl;
188  }
189 
190  m_res_vs_eta->Fill(jet.eta(), r);
191 
192  // We do the same thing to access the scale factors
193  float sf = res_sf.getScaleFactor({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}});
194 
195  // Access up and down variation of the scale factor
196  float sf_up =
198  float sf_down =
200 
201  if (m_debug) {
202  std::cout << "Scale factors (Nominal / Up / Down) : " << sf << " / " << sf_up << " / " << sf_down << std::endl;
203  }
204  }
205 }
206 //---------------------------------------------------------------------------
208 //---------------------------------------------------------------------------
209 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
std::string m_scale_factors_file
static const JetResolutionScaleFactor get(const edm::EventSetup &, const Token &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void endJob() override
JetParameters & setJetEta(float eta)
JetParameters & setRho(float rho)
JetParameters & set(const Binning &bin, float value)
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: Electron.h:6
T getUntrackedParameter(std::string const &, T const &) const
JME::JetResolution::Token m_resolutions_token
~JetResolutionDemo() override
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< double > m_rho_token
float getScaleFactor(const JetParameters &parameters, Variation variation=Variation::NOMINAL, std::string uncertaintySource="") const
JME::JetResolutionScaleFactor::Token m_scale_factors_token
std::string m_resolutions_file
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< pat::Jet > > m_jets_token
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
JetParameters & setJetPt(float pt)
edm::Service< TFileService > fs
JetResolutionDemo(const edm::ParameterSet &)
static const JetResolution get(const edm::EventSetup &, const Token &)