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>
15 using namespace std;
25 inline bool accept(const edm::Event& event,
26  const edm::TriggerResults& triggerTable,
27  const std::string& triggerPath) {
28  bool passed = false;
29  const edm::TriggerNames& triggerNames = event.triggerNames(triggerTable);
30  for (unsigned int i = 0; i < triggerNames.triggerNames().size(); ++i) {
31  if (triggerNames.triggerNames()[i] == triggerPath) {
32  if (triggerTable.accept(i)) {
33  passed = true;
34  break;
35  }
36  }
37  }
38  return passed;
39 }
40 
41 inline bool accept(const edm::Event& event,
42  const edm::TriggerResults& triggerTable,
43  const std::vector<std::string>& triggerPaths) {
44  bool passed = false;
45  for (unsigned int j = 0; j < triggerPaths.size(); ++j) {
46  if (accept(event, triggerTable, triggerPaths[j])) {
47  passed = true;
48  break;
49  }
50  }
51  return passed;
52 }
53 
56 
68 class Calculate {
69  public:
71  Calculate(int maxNJets, double wMass);
73  ~Calculate() {};
74 
76  double massWBoson(const std::vector<reco::Jet>& jets);
78  double massTopQuark(const std::vector<reco::Jet>& jets);
80  // double massBTopQuark(const std::vector<reco::Jet>& jets, std::vector<bool>
81  // bjet);
82  double massBTopQuark(const std::vector<reco::Jet>& jets,
83  std::vector<double> VbtagWP, double btagWP_);
84 
86  double tmassWBoson(reco::RecoCandidate* lep, const reco::MET& met,
87  const reco::Jet& b);
89  double tmassTopQuark(reco::RecoCandidate* lep, const reco::MET& met,
90  const reco::Jet& b);
91 
92  private:
96  void operator()(const std::vector<reco::Jet>& jets);
98  void operator2(const std::vector<reco::Jet>&, std::vector<double>, double);
100  void operator()(const reco::Jet& bJet, reco::RecoCandidate* lepton,
101  const reco::MET& met);
102 
103  private:
105  bool failed_;
109  double wMass_;
111  double massWBoson_;
117  double tmassWBoson_;
120 };
121 
132 
181 template <typename Object>
183  public:
188 
190  bool select(const edm::Event& event);
191  bool select(const edm::Event& event, const std::string& type);
193  bool select(const edm::Event& event, const edm::EventSetup& setup);
194  bool selectVertex(const edm::Event& event);
195 
196  private:
200  int min_, max_;
214  double eidCutValue_;
223 
230 };
231 
233 template <typename Object>
236  : select_(cfg.getParameter<std::string>("select")), jetIDSelect_(0) {
237 
238  src_ =
239  iC.consumes<edm::View<Object> >(cfg.getParameter<edm::InputTag>("src"));
240  // exist otherwise they are initialized with -1
241  cfg.exists("min") ? min_ = cfg.getParameter<int>("min") : min_ = -1;
242  cfg.exists("max") ? max_ = cfg.getParameter<int>("max") : max_ = -1;
243  std::string mygSF = "gedGsfElectrons";
244  gsfEs_ = iC.consumes<edm::View<reco::GsfElectron> >(
245  cfg.getUntrackedParameter<edm::InputTag>("myGSF", mygSF));
246  if (cfg.existsAs<edm::ParameterSet>("electronId")) {
247  edm::ParameterSet elecId =
248  cfg.getParameter<edm::ParameterSet>("electronId");
249  electronId_ = iC.consumes<edm::ValueMap<float> >(
250  elecId.getParameter<edm::InputTag>("src"));
251  eidCutValue_ = elecId.getParameter<double>("cutValue");
252  }
253  if (cfg.exists("jetCorrector")) {
254  jetCorrector_ = cfg.getParameter<std::string>("jetCorrector");
255  }
256  if (cfg.existsAs<edm::ParameterSet>("jetBTagger")) {
257  edm::ParameterSet jetBTagger =
258  cfg.getParameter<edm::ParameterSet>("jetBTagger");
259  btagLabel_ = iC.consumes<reco::JetTagCollection>(
260  jetBTagger.getParameter<edm::InputTag>("label"));
261  btagWorkingPoint_ = jetBTagger.getParameter<double>("workingPoint");
262  }
263  if (cfg.existsAs<edm::ParameterSet>("jetID")) {
265  jetIDLabel_ = iC.consumes<reco::JetIDValueMap>(
266  jetID.getParameter<edm::InputTag>("label"));
268  jetID.getParameter<std::string>("select"));
269  }
270 }
271 
273 template <typename Object>
275 
276  // fetch input collection
278  if (!event.getByToken(src_, src)) return false;
279 
280  // load electronId value map if configured such
282  if (!electronId_.isUninitialized()) {
283  if (!event.getByToken(electronId_, electronId)) return false;
284  }
285 
286  // determine multiplicity of selected objects
287  int n = 0;
288  for (typename edm::View<Object>::const_iterator obj = src->begin();
289  obj != src->end(); ++obj) {
290  // special treatment for electrons
291  if (dynamic_cast<const reco::GsfElectron*>(&*obj)) {
292  unsigned int idx = obj - src->begin();
293  if (electronId_.isUninitialized()
294  ? true
295  : ((double)(*electronId)[src->refAt(idx)] >= eidCutValue_)) {
296  if (select_(*obj)) ++n;
297  }
298  }
299  // normal treatment
300  else {
301  if (select_(*obj)) ++n;
302  }
303  }
304  bool accept =
305  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
306  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
307 }
308 
310 template <typename Object>
312  const std::string& type) {
313  // fetch input collection
315  if (!event.getByToken(src_, src)) return false;
316 
317  // special for gsfElectron
319 
320  // load electronId value map if configured such
322  if (!electronId_.isUninitialized()) {
323  if (!event.getByToken(electronId_, electronId)) return false;
324  }
325 
326  // determine multiplicity of selected objects
327  int n = 0;
328  for (typename edm::View<Object>::const_iterator obj = src->begin();
329  obj != src->end(); ++obj) {
330 
331  // special treatment for PF candidates
332  if (dynamic_cast<const reco::PFCandidate*>(&*obj)) {
333  reco::PFCandidate objtmp = dynamic_cast<const reco::PFCandidate&>(*obj);
334 
335  if (objtmp.muonRef().isNonnull() && type == "muon") {
336  if (select_(*obj)) {
337  ++n;
338  }
339  } else if (objtmp.gsfElectronRef().isNonnull() && type == "electron") {
340  if (select_(*obj)) {
341  if (electronId_.isUninitialized()) {
342  ++n;
343  } else if (((double)(*electronId)[obj->gsfElectronRef()] >=
344  eidCutValue_)) {
345  ++n;
346  }
347  }
348  // idx_gsf++;
349  }
350  }
351 
352  // special treatment for electrons
353  else if (dynamic_cast<const reco::GsfElectron*>(&*obj)) {
354  unsigned int idx = obj - src->begin();
355  if (electronId_.isUninitialized()
356  ? true
357  : ((double)(*electronId)[src->refAt(idx)] >= eidCutValue_)) {
358  if (select_(*obj)) ++n;
359  }
360  }
361 
362  // normal treatment
363  else {
364  if (select_(*obj)) ++n;
365  }
366  }
367  bool accept =
368  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
369  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
370 }
371 
372 template <typename Object>
374  // fetch input collection
376  if (!event.getByToken(src_, src)) return false;
377 
378  // load electronId value map if configured such
380  if (!electronId_.isUninitialized()) {
381  if (!event.getByToken(electronId_, electronId)) return false;
382  }
383 
384  // determine multiplicity of selected objects
385  int n = 0;
386  for (typename edm::View<Object>::const_iterator obj = src->begin();
387  obj != src->end(); ++obj) {
388 
389  if (select_(*obj)) ++n;
390  }
391  bool accept =
392  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
393  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
394 }
395 
397 template <typename Object>
399  const edm::EventSetup& setup) {
400  // fetch input collection
402  if (!event.getByToken(src_, src)) return false;
403 
404  // load btag collection if configured such
405  // NOTE that the JetTagCollection needs an
406  // edm::View to reco::Jets; we have to add
407  // another Handle bjets for this purpose
411  if (!btagLabel_.isUninitialized()) {
412  if (!event.getByToken(src_, bjets)) return false;
413  if (!event.getByToken(btagLabel_, btagger)) return false;
414  if (!event.getByToken(pvs_, pvertex)) return false;
415  }
416 
417  // load jetID value map if configured such
419  if (jetIDSelect_) {
420  if (!event.getByToken(jetIDLabel_, jetID)) return false;
421  }
422 
423  // load jet corrector if configured such
424  const JetCorrector* corrector = 0;
425  if (!jetCorrector_.empty()) {
426  // check whether a jet correcto is in the event setup or not
428  JetCorrectionsRecord>())) {
429  corrector = JetCorrector::getJetCorrector(jetCorrector_, setup);
430  } else {
431  edm::LogVerbatim("TopDQMHelpers")
432  << "\n"
433  << "---------------------------------------------------------------\n"
434  << " No JetCorrectionsRecord available from EventSetup:\n"
435  << " - Jets will not be corrected.\n"
436  << " - If you want to change this add the following lines to your "
437  "cfg file:\n"
438  << "\n"
439  << " ## load jet corrections\n"
440  << " process.load(\""
441  "JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff"
442  "\")\n"
443  << " process.prefer(\"ak5CaloL2L3\")\n"
444  << "---------------------------------------------------------------"
445  "\n";
446  }
447  }
448  // determine multiplicity of selected objects
449  int n = 0;
450  for (typename edm::View<Object>::const_iterator obj = src->begin();
451  obj != src->end(); ++obj) {
452  // check for chosen btag discriminator to be above the
453  // corresponding working point if configured such
454  unsigned int idx = obj - src->begin();
455  if (btagLabel_.isUninitialized() ? true : (*btagger)[bjets->refAt(idx)] >
456  btagWorkingPoint_) {
457  bool passedJetID = true;
458  // check jetID for calo jets
459  if (jetIDSelect_ &&
460  dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())) {
461  passedJetID = (*jetIDSelect_)((*jetID)[src->refAt(idx)]);
462  }
463  if (passedJetID) {
464  // scale jet energy if configured such
465  Object jet = *obj;
466  jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.);
467  if (select_(jet)) ++n;
468  }
469  }
470  }
471  bool accept =
472  (min_ >= 0 ? n >= min_ : true) && (max_ >= 0 ? n <= max_ : true);
473  return (min_ < 0 && max_ < 0) ? (n > 0) : accept;
474 }
475 
476 #endif
477 
478 /* Local Variables: */
479 /* show-trailing-whitespace: t */
480 /* truncate-lines: t */
481 /* 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:250
tuple cfg
Definition: looper.py:237
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
bool select(const edm::Event &event)
apply selection
~Calculate()
default destructor
Definition: TopDQMHelpers.h:73
SelectionStep(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
default constructor
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
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
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:25
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
~SelectionStep()
default destructor
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:90
Helper class for the calculation of a top and a W boson mass estime.
Definition: TopDQMHelpers.h:68
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:50
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
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
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.
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:81
static EventSetupRecordKey makeKey()
edm::EDGetTokenT< edm::View< Object > > src_
input collection
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
static std::string const triggerPaths
Definition: EdmProvDump.cc:42
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