Imperial Analysis
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ICJetSrcHelper.hh
Go to the documentation of this file.
1 #ifndef UserCode_ICHiggsTauTau_ICJetSrcHelper_h
2 #define UserCode_ICHiggsTauTau_ICJetSrcHelper_h
3 
4 #include <memory>
5 #include "boost/functional/hash.hpp"
6 #include "FWCore/ParameterSet/interface/ParameterSet.h"
7 #include "FWCore/Framework/interface/Event.h"
8 #include "FWCore/Framework/interface/EventSetup.h"
9 #include "FWCore/Framework/interface/EDProducer.h"
10 #include "FWCore/Utilities/interface/InputTag.h"
11 #include "DataFormats/JetReco/interface/CaloJet.h"
12 #include "DataFormats/JetReco/interface/PFJet.h"
13 #include "DataFormats/TrackReco/interface/TrackFwd.h"
14 #include "PhysicsTools/SelectorUtils/interface/JetIDSelectionFunctor.h"
15 #if CMSSW_MAJOR_VERSION > 7 || (CMSSW_MAJOR_VERSION == 7 && CMSSW_MINOR_VERSION >= 6)
16 #include "JetMETCorrections/Modules/interface/CorrectedJetProducer.h"
17 #else
18 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
19 #endif
20 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
26 
46 template <class U>
47 struct JetSrcHelper {
48  explicit JetSrcHelper(const edm::ParameterSet &config, edm::ConsumesCollector && collector)
50  config.getParameter<edm::InputTag>("inputJetFlavour")),
51  include_jet_flavour(config.getParameter<bool>("includeJetFlavour")),
52  apply_jec_factors(config.getParameter<bool>("applyJECs")),
53  include_jec_factors(config.getParameter<bool>("includeJECs")),
54  cut(config.getParameter<std::string>("cutAfterJECs")),
55  cut_string(config.getParameter<std::string>("cutAfterJECs")),
56  apply_post_jec_cut(config.getParameter<bool>("applyCutAfterJECs")),
57  input_sv_info((config.getParameter<edm::InputTag>("inputSVInfo"))),
58  include_sv_info_ids(config.getParameter<bool>("requestSVInfo")) {
59  collector.consumes<edm::ValueMap<std::vector<int>>>(input_jet_flavour);
60  collector.consumes<reco::SecondaryVertexTagInfoCollection>(input_sv_info);
62  edm::ParameterSet pset = config.getParameter<edm::ParameterSet>("JECs");
63  std::vector<std::string> vec =
64  pset.getParameterNamesForType<std::string>();
65  for (unsigned i = 0; i < vec.size(); ++i) {
66  jecs.push_back(
67  std::make_pair(vec[i], pset.getParameter<std::string>(vec[i])));
68 #if CMSSW_MAJOR_VERSION > 7 || (CMSSW_MAJOR_VERSION ==7 && CMSSW_MINOR_VERSION >= 6 )
69  collector.consumes<reco::JetCorrector>(jecs[i].second);
70 #else
71  collector.consumes<double>(jecs[i].second);
72 #endif
73  }
74  }
75  edm::ParameterSet btag_pset =
76  config.getParameter<edm::ParameterSet>("BTagDiscriminators");
77  std::vector<std::string> btag_vec =
78  btag_pset.getParameterNamesForType<edm::InputTag>();
79  for (unsigned i = 0; i < btag_vec.size(); ++i) {
80  input_btags.push_back(std::make_pair(
81  btag_vec[i], btag_pset.getParameter<edm::InputTag>(btag_vec[i])));
82  collector.consumes<reco::JetFloatAssociation::Container>(input_btags[i].second);
83  }
84  }
85 
86  void DoSetup(edm::EDProducer * prod) {
87  if (include_sv_info_ids) {
88  prod->produces<reco::SecondaryVertexTagInfoRefVector>("requestedSVInfo");
89  }
90  PrintOptional(1, include_jet_flavour, "includeJetFlavour");
91  PrintOptional(1, apply_jec_factors, "applyJECs");
92  PrintOptional(1, include_jec_factors, "includeJECs");
93  PrintOptional(1, apply_post_jec_cut, "applyCutAfterJECs");
94  PrintOptional(1, include_sv_info_ids, "requestSVInfo");
95  }
97 
98  template <class T>
99  void Produce(edm::Handle<edm::View<U> > const &jets_handle,
100  std::vector<T> *jets, std::vector<unsigned> *passed,
101  edm::Event &event, const edm::EventSetup &setup) {
102  edm::Handle<edm::ValueMap<std::vector<int>> > jet_flavour_handle;
104  event.getByLabel(input_jet_flavour, jet_flavour_handle);
105  #if CMSSW_MAJOR_VERSION > 7 || (CMSSW_MAJOR_VERSION == 7 && CMSSW_MINOR_VERSION >= 6 )
106  std::vector<reco::JetCorrector const *> correctors(jecs.size(), NULL);
107  #else
108  std::vector<JetCorrector const*> correctors(jecs.size(), NULL);
109  #endif
110  for (unsigned i = 0; i < correctors.size(); ++i) {
111  #if CMSSW_MAJOR_VERSION > 7 || (CMSSW_MAJOR_VERSION == 7 && CMSSW_MINOR_VERSION >= 6)
112  edm::Handle<reco::JetCorrector> corrector_handle;
113  event.getByLabel(jecs[i].second,corrector_handle);
114  correctors[i] = corrector_handle.product();
115  #else
116  correctors[i] = JetCorrector::getJetCorrector(jecs[i].second, setup);
117  #endif
118  }
119 
120  edm::Handle<reco::SecondaryVertexTagInfoCollection> sv_info_handle;
121  if (include_sv_info_ids) {
122  event.getByLabel(input_sv_info, sv_info_handle);
123  }
124 
125  std::vector<edm::Handle<reco::JetFloatAssociation::Container> >
126  btag_handles(input_btags.size());
127  for (unsigned i = 0; i < btag_handles.size(); ++i) {
128  event.getByLabel(input_btags[i].second, btag_handles[i]);
129  }
130 
131  std::auto_ptr<reco::SecondaryVertexTagInfoRefVector> sv_info_requests(
132  new reco::SecondaryVertexTagInfoRefVector());
133 
134  for (unsigned i = 0; i < jets_handle->size(); ++i) {
135  U const &src = jets_handle->at(i);
136  double full_jec_factor = 1.0;
137  std::vector<float> factors(correctors.size(), 0.);
138  bool keep_jet = true;
139  // Are there jet energy corrections to apply?
140  if (correctors.size() > 0) {
141  // Yes - make a copy of the source Jet
142  U jet_cpy = jets_handle->at(i);
143  // Loop through each correction and apply
144  for (unsigned j = 0; j < correctors.size(); ++j) {
145  #ifndef CMSSW_4_2_8_patch7
146  #if CMSSW_MAJOR_VERSION > 7 || (CMSSW_MAJOR_VERSION ==7 && CMSSW_MINOR_VERSION >= 6)
147  double factor = correctors[j]->correction(jet_cpy);
148  #else
149  double factor = correctors[j]->correction(jet_cpy, event, setup);
150  #endif
151  #else
152  double factor = correctors[j]->correction(jet_cpy,
153  edm::RefToBase<reco::Jet>(jets_handle->refAt(i)), event, setup);
154  #endif
155  if (apply_jec_factors) full_jec_factor *= factor;
156  factors[j] = factor;
157  jet_cpy.scaleEnergy(factor);
158  }
159  if (apply_post_jec_cut && !cut(jet_cpy)) keep_jet = false;
160  }
161 
162  if (keep_jet) {
163  passed->push_back(i);
164  } else {
165  continue;
166  }
167 
168  jets->push_back(T());
169  T &dest = jets->back();
170 
171  dest.set_id(jet_hasher(&src));
172  dest.set_pt(src.pt() * full_jec_factor);
173  dest.set_eta(src.eta());
174  dest.set_phi(src.phi());
175  dest.set_energy(src.energy() * full_jec_factor);
176  dest.set_charge(src.charge());
177  dest.set_jet_area(src.jetArea());
178  if (include_jet_flavour) {
179  dest.set_parton_flavour((*jet_flavour_handle)[jets_handle->refAt(i)].at(0));
180  dest.set_hadron_flavour((*jet_flavour_handle)[jets_handle->refAt(i)].at(1));
181  }
182  // Only write correction into the output jet if the user wants it
183  if (include_jec_factors) {
184  for (unsigned j = 0; j < correctors.size(); ++j) {
185  dest.SetJecFactor(jecs[j].first, factors[j]);
186  }
187  }
188 
189  if (include_sv_info_ids) {
190  for (unsigned iS = 0; iS < sv_info_handle->size(); ++iS) {
191  if (&src == (*sv_info_handle)[iS].jet().get()) {
192  reco::SecondaryVertexTagInfoRef ref(sv_info_handle, iS);
193  sv_info_requests->push_back(ref);
194  std::vector<std::size_t> jet_sv_ids;
195  for (unsigned iV = 0; iV < (*sv_info_handle)[iS].nVertices();
196  ++iV) {
197  reco::Vertex const &sv =
198  (*sv_info_handle)[iS].secondaryVertex(iV);
199  jet_sv_ids.push_back(vertex_hasher(&sv));
200  }
201  dest.set_secondary_vertices(jet_sv_ids);
202  }
203  }
204  }
205 
206  for (unsigned iB = 0; iB < btag_handles.size(); ++iB) {
207  dest.SetBDiscriminator(input_btags[iB].first,
208  (*btag_handles[iB])[edm::RefToBase<reco::Jet>(
209  jets_handle->refAt(i))]);
210  }
211  }
212  if (include_sv_info_ids) {
213  event.put(sv_info_requests, "requestedSVInfo");
214  }
215  }
216 
217  boost::hash<U const*> jet_hasher;
218  boost::hash<reco::Vertex const*> vertex_hasher;
219 
220  edm::InputTag input_jet_flavour;
222 
223  std::vector<std::pair<std::string, std::string> > jecs;
226 
227  StringCutObjectSelector<U> cut;
228  std::string cut_string;
230 
231  edm::InputTag input_sv_info;
233 
234  std::vector<std::pair<std::string, edm::InputTag> > input_btags;
235 };
236 
246 template <>
247 struct JetSrcHelper<pat::Jet> {
248  explicit JetSrcHelper(const edm::ParameterSet &config, edm::ConsumesCollector && collector)
249  : include_jet_flavour(config.getParameter<bool>("includeJetFlavour")),
250  include_jec_factors(config.getParameter<bool>("includeJECs")),
251  input_sv_info((config.getParameter<edm::InputTag>("inputSVInfo"))),
252  include_sv_info_ids(config.getParameter<bool>("requestSVInfo")),
253  is_slimmed(config.getParameter<bool>("isSlimmed")),
254  slimmed_puid_label(config.getParameter<std::string>("slimmedPileupIDLabel")) {
255  collector.consumes<reco::SecondaryVertexTagInfoCollection>(input_sv_info);
256  }
257  void DoSetup(edm::EDProducer * prod) {
258  if (include_sv_info_ids) {
259  prod->produces<reco::SecondaryVertexTagInfoRefVector>("requestedSVInfo");
260  }
261  PrintOptional(1, include_jet_flavour, "includeJetFlavour");
262  PrintOptional(1, include_jec_factors, "includeJECs");
263  PrintOptional(1, include_sv_info_ids, "requestSVInfo");
264  }
266 
267  template <class T>
268  void Produce(edm::Handle<edm::View<pat::Jet> > const &jets_handle,
269  std::vector<T> *jets, std::vector<unsigned> *passed,
270  edm::Event &event, const edm::EventSetup &setup) {
271  edm::Handle<reco::SecondaryVertexTagInfoCollection> sv_info_handle;
272  if (include_sv_info_ids) {
273  event.getByLabel(input_sv_info, sv_info_handle);
274  }
275 
276  std::auto_ptr<reco::SecondaryVertexTagInfoRefVector> sv_info_requests(
277  new reco::SecondaryVertexTagInfoRefVector());
278 
279  for (unsigned i = 0; i < jets_handle->size(); ++i) {
280  pat::Jet const &src = jets_handle->at(i);
281  passed->push_back(i);
282  jets->push_back(T());
283  T &dest = jets->back();
284 
285  dest.set_id(jet_hasher(&src));
286  dest.set_pt(src.pt());
287  dest.set_eta(src.eta());
288  dest.set_phi(src.phi());
289  dest.set_energy(src.energy());
290  dest.set_charge(src.charge());
291  dest.set_jet_area(src.jetArea());
292  if (include_jet_flavour) {
293  dest.set_parton_flavour(src.partonFlavour());
294  dest.set_hadron_flavour(src.hadronFlavour());
295  }
296  // Only write correction into the output jet if the user wants it
297  if (include_jec_factors && src.availableJECLevels().size() >= 1) {
298  std::vector<std::string> jec_levels = src.availableJECLevels();
299  for (unsigned j = 1; j < jec_levels.size(); ++j) {
300  dest.SetJecFactor(
301  jec_levels[j],
302  src.jecFactor(jec_levels[j]) / src.jecFactor(jec_levels[j - 1]));
303  }
304  }
305 
306  if (include_sv_info_ids) {
307  for (unsigned iS = 0; iS < sv_info_handle->size(); ++iS) {
308  if (&(*(src.originalObjectRef())) ==
309  (*sv_info_handle)[iS].jet().get()) {
310  reco::SecondaryVertexTagInfoRef ref(sv_info_handle, iS);
311  sv_info_requests->push_back(ref);
312  std::vector<std::size_t> jet_sv_ids;
313  for (unsigned iV = 0; iV < (*sv_info_handle)[iS].nVertices();
314  ++iV) {
315  reco::Vertex const &sv =
316  (*sv_info_handle)[iS].secondaryVertex(iV);
317  jet_sv_ids.push_back(vertex_hasher(&sv));
318  }
319  dest.set_secondary_vertices(jet_sv_ids);
320  }
321  }
322  }
323 
324  for (unsigned iB = 0; iB < src.getPairDiscri().size(); ++iB) {
325  dest.SetBDiscriminator(src.getPairDiscri().at(iB).first,
326  src.getPairDiscri().at(iB).second);
327  }
328  }
329  if (include_sv_info_ids) {
330  event.put(sv_info_requests, "requestedSVInfo");
331  }
332  }
333 
334  boost::hash<pat::Jet const*> jet_hasher;
335  boost::hash<reco::Vertex const*> vertex_hasher;
336 
339  edm::InputTag input_sv_info;
342  std::string slimmed_puid_label;
343 };
344 
345 #endif
void DoSetup(edm::EDProducer *prod)
std::vector< std::pair< std::string, edm::InputTag > > input_btags
void DoSetup(edm::EDProducer *prod)
JetSrcHelper(const edm::ParameterSet &config, edm::ConsumesCollector &&collector)
boost::hash< U const * > jet_hasher
A struct to help with the configuration of the ICJetProducer, in particular the aspects that depend o...
edm::InputTag input_sv_info
std::string cut_string
void Produce(edm::Handle< edm::View< pat::Jet > > const &jets_handle, std::vector< T > *jets, std::vector< unsigned > *passed, edm::Event &event, const edm::EventSetup &setup)
boost::hash< reco::Vertex const * > vertex_hasher
std::vector< std::pair< std::string, std::string > > jecs
bool include_jet_flavour
JetSrcHelper(const edm::ParameterSet &config, edm::ConsumesCollector &&collector)
bool include_sv_info_ids
void PrintOptional(unsigned depth, bool value, std::string text)
Definition: Consumes.h:19
boost::hash< reco::Vertex const * > vertex_hasher
bool include_jec_factors
boost::hash< pat::Jet const * > jet_hasher
edm::InputTag input_jet_flavour
void Produce(edm::Handle< edm::View< U > > const &jets_handle, std::vector< T > *jets, std::vector< unsigned > *passed, edm::Event &event, const edm::EventSetup &setup)
StringCutObjectSelector< U > cut