CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopDQMHelpers.h
Go to the documentation of this file.
1 #ifndef TOPDQMHELPERS
2 #define TOPDQMHELPERS
3 
4 #include <string>
5 #include <vector>
6 #include <iostream>
24 inline bool accept(const edm::Event& event,
25  const edm::TriggerResults& triggerTable,
26  const std::string& triggerPath) {
27  bool passed = false;
28  const edm::TriggerNames& triggerNames = event.triggerNames(triggerTable);
29  for (unsigned int i = 0; i < triggerNames.triggerNames().size(); ++i) {
30  if (triggerNames.triggerNames()[i] == triggerPath) {
31  if (triggerTable.accept(i)) {
32  passed = true;
33  break;
34  }
35  }
36  }
37  return passed;
38 }
39 
40 inline bool accept(const edm::Event& event,
41  const edm::TriggerResults& triggerTable,
42  const std::vector<std::string>& triggerPaths) {
43  bool passed = false;
44  for (unsigned int j = 0; j < triggerPaths.size(); ++j) {
45  if (accept(event, triggerTable, triggerPaths[j])) {
46  passed = true;
47  break;
48  }
49  }
50  return passed;
51 }
52 
55 
67 class Calculate {
68  public:
70  Calculate(int maxNJets, double wMass);
72  ~Calculate() {};
73 
75  double massWBoson(const std::vector<reco::Jet>& jets);
77  double massTopQuark(const std::vector<reco::Jet>& jets);
79  // double massBTopQuark(const std::vector<reco::Jet>& jets, std::vector<bool>
80  // bjet);
81  double massBTopQuark(const std::vector<reco::Jet>& jets,
82  std::vector<double> VbtagWP, double btagWP_);
83 
85  double tmassWBoson(reco::RecoCandidate* lep, const reco::MET& met,
86  const reco::Jet& b);
88  double tmassTopQuark(reco::RecoCandidate* lep, const reco::MET& met,
89  const reco::Jet& b);
90 
91  private:
95  void operator()(const std::vector<reco::Jet>& jets);
97  void operator2(const std::vector<reco::Jet>&, std::vector<double>, double);
99  void operator()(const reco::Jet& bJet, reco::RecoCandidate* lepton,
100  const reco::MET& met);
101 
102  private:
104  bool failed_;
108  double wMass_;
110  double massWBoson_;
116  double tmassWBoson_;
119 };
120 
131 
180 template <typename Object>
182  public:
187 
189  bool select(const edm::Event& event);
190  bool select(const edm::Event& event, const std::string& type);
192  bool select(const edm::Event& event, const edm::EventSetup& setup);
193  bool selectVertex(const edm::Event& event);
194 
195  private:
199  int min_, max_;
213  double eidCutValue_;
222 
229 };
230 
232 template <typename Object>
235  : select_(cfg.getParameter<std::string>("select")), jetIDSelect_(0) {
236 
237  src_ =
238  iC.consumes<edm::View<Object> >(cfg.getParameter<edm::InputTag>("src"));
239  // exist otherwise they are initialized with -1
240  cfg.exists("min") ? min_ = cfg.getParameter<int>("min") : min_ = -1;
241  cfg.exists("max") ? max_ = cfg.getParameter<int>("max") : max_ = -1;
242  std::string mygSF = "gedGsfElectrons";
243  gsfEs_ = iC.consumes<edm::View<reco::GsfElectron> >(
244  cfg.getUntrackedParameter<edm::InputTag>("myGSF", mygSF));
245  if (cfg.existsAs<edm::ParameterSet>("electronId")) {
246  edm::ParameterSet elecId =
247  cfg.getParameter<edm::ParameterSet>("electronId");
248  electronId_ = iC.consumes<edm::ValueMap<float> >(
249  elecId.getParameter<edm::InputTag>("src"));
250  eidCutValue_ = elecId.getParameter<double>("cutValue");
251  }
252  if (cfg.exists("jetCorrector")) {
253  jetCorrector_ = cfg.getParameter<std::string>("jetCorrector");
254  }
255  if (cfg.existsAs<edm::ParameterSet>("jetBTagger")) {
256  edm::ParameterSet jetBTagger =
257  cfg.getParameter<edm::ParameterSet>("jetBTagger");
258  btagLabel_ = iC.consumes<reco::JetTagCollection>(
259  jetBTagger.getParameter<edm::InputTag>("label"));
260  btagWorkingPoint_ = jetBTagger.getParameter<double>("workingPoint");
261  }
262  if (cfg.existsAs<edm::ParameterSet>("jetID")) {
263  edm::ParameterSet jetID = cfg.getParameter<edm::ParameterSet>("jetID");
264  jetIDLabel_ = iC.consumes<reco::JetIDValueMap>(
265  jetID.getParameter<edm::InputTag>("label"));
267  jetID.getParameter<std::string>("select"));
268  }
269 }
270 
272 template <typename Object>
274 
275  // fetch input collection
277  if (!event.getByToken(src_, src)) return false;
278 
279  // load electronId value map if configured such
281  if (!electronId_.isUninitialized()) {
282  if (!event.getByToken(electronId_, electronId)) return false;
283  }
284 
285  // determine multiplicity of selected objects
286  int n = 0;
287  for (typename edm::View<Object>::const_iterator obj = src->begin();
288  obj != src->end(); ++obj) {
289  // special treatment for electrons
290  if (dynamic_cast<const reco::GsfElectron*>(&*obj)) {
291  unsigned int idx = obj - src->begin();
292  if (electronId_.isUninitialized()
293  ? true
294  : ((double)(*electronId)[src->refAt(idx)] >= eidCutValue_)) {
295  if (select_(*obj)) ++n;
296  }
297  }
298  // normal treatment
299  else {
300  if (select_(*obj)) ++n;
301  }
302  }
303  bool accept =
304  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
305  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
306 }
307 
309 template <typename Object>
311  const std::string& type) {
312  // fetch input collection
314  if (!event.getByToken(src_, src)) return false;
315 
316  // special for gsfElectron
318 
319  // load electronId value map if configured such
321  if (!electronId_.isUninitialized()) {
322  if (!event.getByToken(electronId_, electronId)) return false;
323  }
324 
325  // determine multiplicity of selected objects
326  int n = 0;
327  for (typename edm::View<Object>::const_iterator obj = src->begin();
328  obj != src->end(); ++obj) {
329 
330  // special treatment for PF candidates
331  if (dynamic_cast<const reco::PFCandidate*>(&*obj)) {
332  reco::PFCandidate objtmp = dynamic_cast<const reco::PFCandidate&>(*obj);
333 
334  if (objtmp.muonRef().isNonnull() && type == "muon") {
335  if (select_(*obj)) {
336  ++n;
337  }
338  } else if (objtmp.gsfElectronRef().isNonnull() && type == "electron") {
339  if (select_(*obj)) {
340  if (electronId_.isUninitialized()) {
341  ++n;
342  } else if (((double)(*electronId)[obj->gsfElectronRef()] >=
343  eidCutValue_)) {
344  ++n;
345  }
346  }
347  // idx_gsf++;
348  }
349  }
350 
351  // special treatment for electrons
352  else if (dynamic_cast<const reco::GsfElectron*>(&*obj)) {
353  unsigned int idx = obj - src->begin();
354  if (electronId_.isUninitialized()
355  ? true
356  : ((double)(*electronId)[src->refAt(idx)] >= eidCutValue_)) {
357  if (select_(*obj)) ++n;
358  }
359  }
360 
361  // normal treatment
362  else {
363  if (select_(*obj)) ++n;
364  }
365  }
366  bool accept =
367  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
368  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
369 }
370 
371 template <typename Object>
373  // fetch input collection
375  if (!event.getByToken(src_, src)) return false;
376 
377  // load electronId value map if configured such
379  if (!electronId_.isUninitialized()) {
380  if (!event.getByToken(electronId_, electronId)) return false;
381  }
382 
383  // determine multiplicity of selected objects
384  int n = 0;
385  for (typename edm::View<Object>::const_iterator obj = src->begin();
386  obj != src->end(); ++obj) {
387 
388  if (select_(*obj)) ++n;
389  }
390  bool accept =
391  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
392  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
393 }
394 
396 template <typename Object>
398  const edm::EventSetup& setup) {
399  // fetch input collection
401  if (!event.getByToken(src_, src)) return false;
402 
403  // load btag collection if configured such
404  // NOTE that the JetTagCollection needs an
405  // edm::View to reco::Jets; we have to add
406  // another Handle bjets for this purpose
410  if (!btagLabel_.isUninitialized()) {
411  if (!event.getByToken(src_, bjets)) return false;
412  if (!event.getByToken(btagLabel_, btagger)) return false;
413  if (!event.getByToken(pvs_, pvertex)) return false;
414  }
415 
416  // load jetID value map if configured such
418  if (jetIDSelect_) {
419  if (!event.getByToken(jetIDLabel_, jetID)) return false;
420  }
421 
422  // load jet corrector if configured such
423  const JetCorrector* corrector = 0;
424  if (!jetCorrector_.empty()) {
425  // check whether a jet correcto is in the event setup or not
427  JetCorrectionsRecord>())) {
428  corrector = JetCorrector::getJetCorrector(jetCorrector_, setup);
429  } else {
430  edm::LogVerbatim("TopDQMHelpers")
431  << "\n"
432  << "---------------------------------------------------------------\n"
433  << " No JetCorrectionsRecord available from EventSetup:\n"
434  << " - Jets will not be corrected.\n"
435  << " - If you want to change this add the following lines to your "
436  "cfg file:\n"
437  << "\n"
438  << " ## load jet corrections\n"
439  << " process.load(\""
440  "JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff"
441  "\")\n"
442  << " process.prefer(\"ak5CaloL2L3\")\n"
443  << "---------------------------------------------------------------"
444  "\n";
445  }
446  }
447  // determine multiplicity of selected objects
448  int n = 0;
449  for (typename edm::View<Object>::const_iterator obj = src->begin();
450  obj != src->end(); ++obj) {
451  // check for chosen btag discriminator to be above the
452  // corresponding working point if configured such
453  unsigned int idx = obj - src->begin();
454  if (btagLabel_.isUninitialized() ? true : (*btagger)[bjets->refAt(idx)] >
455  btagWorkingPoint_) {
456  bool passedJetID = true;
457  // check jetID for calo jets
458  if (jetIDSelect_ &&
459  dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())) {
460  passedJetID = (*jetIDSelect_)((*jetID)[src->refAt(idx)]);
461  }
462  if (passedJetID) {
463  // scale jet energy if configured such
464  Object jet = *obj;
465  jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.);
466  if (select_(jet)) ++n;
467  }
468  }
469  }
470  bool accept =
471  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
472  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
473 }
474 
475 #endif
476 
477 /* Local Variables: */
478 /* show-trailing-whitespace: t */
479 /* truncate-lines: t */
480 /* End: */
double massWBoson_
cache of w boson mass estimate
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
tuple cfg
Definition: looper.py:259
double massBTopQuark(const std::vector< reco::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
calculate b-tagged t-quark mass estimate
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
bool select(const edm::Event &event)
apply selection
~Calculate()
default destructor
Definition: TopDQMHelpers.h:72
SelectionStep(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
default constructor
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
double massTopQuark_
cache of top quark mass estimate
double tmassTopQuark_
cache of top quark transverse mass estimate
Base class for all types of Jets.
Definition: Jet.h:20
bool selectVertex(const edm::Event &event)
bool accept() const
Has at least one path accepted the event?
edm::EDGetTokenT< edm::ValueMap< float > > electronId_
electronId label as extra selection type
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate t-quark mass estimate
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< edm::View< reco::GsfElectron > > gsfEs_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
~SelectionStep()
default destructor
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:91
void operator2(const std::vector< reco::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet
Helper class for the calculation of a top and a W boson mass estime.
Definition: TopDQMHelpers.h:67
int min_
min/max for object multiplicity
int maxNJets_
max. number of jets to be considered
edm::EDGetTokenT< reco::JetTagCollection > btagLabel_
choice for b-tag as extra selection type
tuple corrector
Definition: mvaPFMET_cff.py:89
Definition: MET.h:42
vector< PseudoJet > jets
reco::GsfElectronRef gsfElectronRef() const
return a reference to the corresponding GsfElectron if any
Definition: PFCandidate.cc:574
int j
Definition: DBlmapReader.cc:9
double eidCutValue_
double wMass_
paramater of the w boson mass
double tmassTopQuark(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate top quark transverse mass estimate
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
reco::MuonRef muonRef() const
Definition: PFCandidate.cc:450
double tmassWBoson(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate W boson transverse mass estimate
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
Definition: JetCorrector.cc:50
double b
Definition: hdecay.h:120
double btagWorkingPoint_
choice of b-tag working point as extra selection type
Templated helper class to allow a selection on a certain object collection.
void operator()(const std::vector< reco::Jet > &jets)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
bool failed_
indicate failed associations
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
static EventSetupRecordKey makeKey()
edm::EDGetTokenT< edm::View< Object > > src_
input collection
Calculate(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:3
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
static std::string const triggerPaths
Definition: EdmProvDump.cc:41
std::string jetCorrector_
jet corrector as extra selection type
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
StringCutObjectSelector< Object > select_
string cut selector
double tmassWBoson_
cache of W boson transverse mass estimate
double massBTopQuark_
cache of b-tagged top quark mass estimate
edm::EDGetTokenT< reco::JetIDValueMap > jetIDLabel_
jetID as an extra selection type