CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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&);
29 
30  private:
31  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
32  virtual void endJob() override;
33 
35 
36  bool m_debug = false;
37  bool m_use_conddb = false;
38 
42 
44 
45  TH2* m_res_vs_eta = nullptr;
46 };
47 //
48 //----------- Class Implementation ------------------------------------------
49 //
50 //---------------------------------------------------------------------------
52 {
53  m_jets_token = consumes<std::vector<pat::Jet>>(iConfig.getParameter<edm::InputTag>("jets"));
54  m_debug = iConfig.getUntrackedParameter<bool>("debug", false);
55  m_use_conddb = iConfig.getUntrackedParameter<bool>("useCondDB", false);
56 
57  if (m_use_conddb)
58  m_payload = iConfig.getParameter<std::string>("payload");
59  else {
60  m_resolutions_file = iConfig.getParameter<edm::FileInPath>("resolutionsFile").fullPath();
61  m_scale_factors_file = iConfig.getParameter<edm::FileInPath>("scaleFactorsFile").fullPath();
62  }
63 }
64 //---------------------------------------------------------------------------
66 {
67 
68 }
69 
70 //---------------------------------------------------------------------------
72 {
73 
75  iEvent.getByToken(m_jets_token, jets);
76 
77  // Access jet resolution and scale factor from the condition database
78  // or from text files
81 
82  // Two differents way to create a class instance
83  if (m_use_conddb) {
84  // First way, using the get() static method
85  resolution = JME::JetResolution::get(iSetup, m_payload);
87  } else {
88  // Second way, using the constructor
91  }
92 
93  if (m_debug) {
94  // Dump resolution
95  resolution.dump();
96  }
97 
98  if (! m_res_vs_eta) {
99 
100  // 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.
101 
102  // Get the list of bins of this object
103  const std::vector<JME::Binning>& bins = resolution.getResolutionObject()->getDefinition().getBins();
104 
105  // Check that the first bin is eta
106  if ((bins.size()) && (bins[0] == JME::Binning::JetEta)) {
107 
108  const std::vector<JME::JetResolutionObject::Record> records = resolution.getResolutionObject()->getRecords();
109  // Get all records from the object. Each record correspond to a different binning and different parameters
110 
111  std::vector<float> etas;
112  for (const auto& record: records) {
113  if (etas.empty()) {
114  etas.push_back(record.getBinsRange()[0].min);
115  etas.push_back(record.getBinsRange()[0].max);
116  } else {
117  etas.push_back(record.getBinsRange()[0].max);
118  }
119  }
120 
121  std::vector<float> res;
122  for (std::size_t i = 0; i < 40; i++) {
123  res.push_back(i * 0.005);
124  }
125 
126  m_res_vs_eta = fs->make<TH2F>("res_vs_eta", "res_vs_eta", etas.size() - 1, &etas[0], res.size() - 1, &res[0]);
127  }
128  }
129 
130  for (const auto& jet: *jets) {
131  if (m_debug) {
132  std::cout << "New jet; pt=" << jet.pt() << " eta=" << jet.eta() << " phi=" << jet.phi() << " e=" << jet.energy() << std::endl;
133  }
134 
135 
136  // Get resolution for this jet
137  // The method to use is getResolution(const JME::JetParameters& parameters) from a JetResolution object
138 
139  // 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
140  // probably change in the futur when PU dependency is added. Please keep an eye on the twiki page
141  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution
142  // to stay up-to-date. All currently supported parameters (ie, 'set' functions) are available here:
143  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyResolution#List_of_supported_parameters
144 
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 
153  // You can also chain calls
154 
155  JME::JetParameters parameters_2;
156  parameters_2.setJetPt(jet.pt()).setJetEta(jet.eta());
157 
158  // Second, using the set() function
159  JME::JetParameters parameters_3;
160  parameters_3.set(JME::Binning::JetPt, jet.pt());
161  parameters_3.set({JME::Binning::JetEta, jet.eta()});
162 
163  // Or
164 
165  JME::JetParameters parameters_4;
166  parameters_4.set(JME::Binning::JetPt, jet.pt()).set(JME::Binning::JetEta, jet.eta());
167 
168  // Third, using a initializer_list
169  JME::JetParameters parameters_5 = {{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}};
170 
171  // Now, get the resolution
172 
173  float r = resolution.getResolution(parameters_1);
174  if (m_debug) {
175  std::cout << "Resolution with parameters_1: " << r << std::endl;
176  std::cout << "Resolution with parameters_2: " << resolution.getResolution(parameters_2) << std::endl;
177  std::cout << "Resolution with parameters_3: " << resolution.getResolution(parameters_3) << std::endl;
178  std::cout << "Resolution with parameters_4: " << resolution.getResolution(parameters_4) << std::endl;
179  std::cout << "Resolution with parameters_5: " << resolution.getResolution(parameters_5) << std::endl;
180 
181  // You can also use a shortcut to get the resolution
182  float r2 = resolution.getResolution({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}});
183  std::cout << "Resolution using shortcut : " << r2 << std::endl;
184  }
185 
186  m_res_vs_eta->Fill(jet.eta(), r);
187 
188  // We do the same thing to access the scale factors
189  float sf = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}});
190 
191  // Access up and down variation of the scale factor
192  float sf_up = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}}, Variation::UP);
193  float sf_down = res_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}}, Variation::DOWN);
194 
195  if (m_debug) {
196  std::cout << "Scale factors (Nominal / Up / Down) : " << sf << " / " << sf_up << " / " << sf_down << std::endl;
197  }
198  }
199 }
200 //---------------------------------------------------------------------------
202 {
203 
204 }
205 //---------------------------------------------------------------------------
206 //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
int i
Definition: DBlmapReader.cc:9
float getResolution(const JetParameters &parameters) const
virtual void endJob() override
JetParameters & setJetEta(float eta)
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:464
#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)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:230
void dump() const
Definition: JetResolution.h:32
const std::vector< Record > & getRecords() const
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 &)
tuple cout
Definition: gather_cfg.py:121
JetParameters & setJetPt(float pt)
edm::Service< TFileService > fs
JetResolutionDemo(const edm::ParameterSet &)