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 
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 
43 
45 
46  TH2* m_res_vs_eta = nullptr;
47 };
48 //
49 //----------- Class Implementation ------------------------------------------
50 //
51 //---------------------------------------------------------------------------
53 {
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");
61  else {
62  m_resolutions_file = iConfig.getParameter<edm::FileInPath>("resolutionsFile").fullPath();
63  m_scale_factors_file = iConfig.getParameter<edm::FileInPath>("scaleFactorsFile").fullPath();
64  }
65 }
66 //---------------------------------------------------------------------------
68 {
69 
70 }
71 
72 //---------------------------------------------------------------------------
74 {
75 
77  iEvent.getByToken(m_jets_token, jets);
78 
80  iEvent.getByToken(m_rho_token, rho);
81 
82  // Access jet resolution and scale factor from the condition database
83  // or from text files
86 
87  // Two differents way to create a class instance
88  if (m_use_conddb) {
89  // First way, using the get() static method
90  resolution = JME::JetResolution::get(iSetup, m_payload);
92  } else {
93  // Second way, using the constructor
96  }
97 
98  if (m_debug) {
99  // Dump resolution
100  resolution.dump();
101  }
102 
103  if (! m_res_vs_eta) {
104 
105  // 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.
106 
107  // Get the list of bins of this object
108  const std::vector<JME::Binning>& bins = resolution.getResolutionObject()->getDefinition().getBins();
109 
110  // Check that the first bin is eta
111  if ((bins.size()) && (bins[0] == JME::Binning::JetEta)) {
112 
113  const std::vector<JME::JetResolutionObject::Record> records = resolution.getResolutionObject()->getRecords();
114  // Get all records from the object. Each record correspond to a different binning and different parameters
115 
116  std::vector<float> etas;
117  for (const auto& record: records) {
118  if (etas.empty()) {
119  etas.push_back(record.getBinsRange()[0].min);
120  etas.push_back(record.getBinsRange()[0].max);
121  } else {
122  etas.push_back(record.getBinsRange()[0].max);
123  }
124  }
125 
126  std::vector<float> res;
127  for (std::size_t i = 0; i < 40; i++) {
128  res.push_back(i * 0.005);
129  }
130 
131  m_res_vs_eta = fs->make<TH2F>("res_vs_eta", "res_vs_eta", etas.size() - 1, &etas[0], res.size() - 1, &res[0]);
132  }
133  }
134 
135  for (const auto& jet: *jets) {
136  if (m_debug) {
137  std::cout << "New jet; pt=" << jet.pt() << " eta=" << jet.eta() << " phi=" << jet.phi() << " e=" << jet.energy() << " rho=" << *rho << std::endl;
138  }
139 
140 
141  // Get resolution for this jet
142  // The method to use is getResolution(const JME::JetParameters& parameters) from a JetResolution object
143 
144  // 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
145  // probably change in the futur when PU dependency is added. Please keep an eye on the twiki page
146  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution
147  // to stay up-to-date. All currently supported parameters (ie, 'set' functions) are available here:
148  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution#List_of_supported_parameters
149 
150 
151  // Three way to create a JetParameters object
152 
153  // First, using 'set' functions
154  JME::JetParameters parameters_1;
155  parameters_1.setJetPt(jet.pt());
156  parameters_1.setJetEta(jet.eta());
157  parameters_1.setRho(*rho);
158 
159  // You can also chain calls
160 
161  JME::JetParameters parameters_2;
162  parameters_2.setJetPt(jet.pt()).setJetEta(jet.eta()).setRho(*rho);
163 
164  // Second, using the set() function
165  JME::JetParameters parameters_3;
166  parameters_3.set(JME::Binning::JetPt, jet.pt());
167  parameters_3.set({JME::Binning::JetEta, jet.eta()});
168  parameters_3.set({JME::Binning::Rho, *rho});
169 
170  // Or
171 
172  JME::JetParameters parameters_4;
173  parameters_4.set(JME::Binning::JetPt, jet.pt()).set(JME::Binning::JetEta, jet.eta()).set(JME::Binning::Rho, *rho);
174 
175  // Third, using a initializer_list
176  JME::JetParameters parameters_5 = {{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}, {JME::Binning::Rho, *rho}};
177 
178  // Now, get the resolution
179 
180  float r = resolution.getResolution(parameters_1);
181  if (m_debug) {
182  std::cout << "Resolution with parameters_1: " << r << std::endl;
183  std::cout << "Resolution with parameters_2: " << resolution.getResolution(parameters_2) << std::endl;
184  std::cout << "Resolution with parameters_3: " << resolution.getResolution(parameters_3) << std::endl;
185  std::cout << "Resolution with parameters_4: " << resolution.getResolution(parameters_4) << std::endl;
186  std::cout << "Resolution with parameters_5: " << resolution.getResolution(parameters_5) << std::endl;
187 
188  // You can also use a shortcut to get the resolution
189  float r2 = resolution.getResolution({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}, {JME::Binning::Rho, *rho}});
190  std::cout << "Resolution using shortcut : " << r2 << std::endl;
191  }
192 
193  m_res_vs_eta->Fill(jet.eta(), r);
194 
195  // We do the same thing to access the scale factors
196  float sf = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}});
197 
198  // Access up and down variation of the scale factor
199  float sf_up = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}}, Variation::UP);
200  float sf_down = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}}, Variation::DOWN);
201 
202  if (m_debug) {
203  std::cout << "Scale factors (Nominal / Up / Down) : " << sf << " / " << sf_up << " / " << sf_down << std::endl;
204  }
205  }
206 }
207 //---------------------------------------------------------------------------
209 {
210 
211 }
212 //---------------------------------------------------------------------------
213 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::string m_scale_factors_file
float getResolution(const JetParameters &parameters) const
void endJob() override
JetParameters & setJetEta(float eta)
JetParameters & setRho(float rho)
const std::vector< Binning > & getBins() const
static const JetResolution get(const edm::EventSetup &, const std::string &)
JetCorrectorParameters::Record record
Definition: classes.h:7
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
JetParameters & set(const Binning &bin, float value)
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: Electron.h:6
~JetResolutionDemo() override
int iEvent
Definition: GenABIO.cc:230
void dump() const
Definition: JetResolution.h:32
const std::vector< Record > & getRecords() const
edm::EDGetTokenT< double > m_rho_token
const Definition & getDefinition() const
vector< PseudoJet > jets
std::string m_resolutions_file
edm::EDGetTokenT< std::vector< pat::Jet > > m_jets_token
const JetResolutionObject * getResolutionObject() const
Definition: JetResolution.h:37
float getScaleFactor(const JetParameters &parameters, Variation variation=Variation::NOMINAL) const
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
JetParameters & setJetPt(float pt)
edm::Service< TFileService > fs
JetResolutionDemo(const edm::ParameterSet &)