CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopHLTOfflineDQMHelper.h
Go to the documentation of this file.
1 #ifndef TOPHLTOFFLINEDQMHELPER
2 #define TOPHLTOFFLINEDQMHELPER
3 
4 #include <string>
5 #include <vector>
6 //#include <math>
7 #include "TString.h"
14 
15 /*Originally from DQM/Physics package, written by Roger Wolf and Jeremy Andrea*/
25 inline bool
26 acceptHLT(const edm::Event& event, const edm::TriggerResults& triggerTable, const std::string& triggerPath)
27 {
28  bool passed=false;
29  const edm::TriggerNames& triggerNames = event.triggerNames(triggerTable);
30  for(unsigned int i=0; i<triggerNames.triggerNames().size(); ++i){
31  TString name = triggerNames.triggerNames()[i].c_str();
32  //std::cout << name << " " << triggerTable.accept(i) << std::endl;
33  if(name.Contains(TString(triggerPath.c_str()), TString::kIgnoreCase)) {
34  if(triggerTable.accept(i)){
35  passed=true;
36  break;
37  }
38  }
39  }
40  return passed;
41 }
42 
43 inline bool
44 acceptHLT(const edm::Event& event, const edm::TriggerResults& triggerTable, const std::vector<std::string>& triggerPaths)
45 {
46  bool passed=false;
47  for(unsigned int j=0; j<triggerPaths.size(); ++j){
48  if(acceptHLT(event, triggerTable, triggerPaths[j])){
49  passed=true;
50  break;
51  }
52  }
53  return passed;
54 }
55 
56 
59 
72 class CalculateHLT {
73  public:
75  CalculateHLT(int maxNJets, double wMass);
78 
80  double massWBoson(const std::vector<reco::Jet>& jets);
82  double massTopQuark(const std::vector<reco::Jet>& jets);
84 /* double tmassWBoson(const T& mu, const reco::CaloMET& met, const reco::Jet& b);
86  double tmassTopQuark(const T& mu, const reco::CaloMET& met, const reco::Jet& b);
88  double masslb(const T& mu, const reco::CaloMET& met, const reco::Jet& b);*/
89 
91  double tmassWBoson(reco::RecoCandidate* mu, const reco::MET& met, const reco::Jet& b);
93  double tmassTopQuark(reco::RecoCandidate* mu, const reco::MET& met, const reco::Jet& b);
95  double masslb(reco::RecoCandidate* mu, const reco::MET& met, const reco::Jet& b);
96 
97  private:
101  void operator()(const std::vector<reco::Jet>& jets);
102  void operator()(const reco::Jet& bJet, reco::RecoCandidate* lepton, const reco::MET& met);
103  private:
105  bool failed_;
109  double wMass_;
111  double massWBoson_;
115  double tmassWBoson_;
119  double mlb_;
120 
121 
122 };
123 
124 
136 
177  public:
178  virtual bool select(const edm::Event& event) {
179  return false;
180  };
181 
182  virtual bool select(const edm::Event& event, const edm::EventSetup& setup) {
183  assert(false);
184  return false;
185  };
186 };
187 
188 template <typename Object>
190 public:
194  virtual ~SelectionStepHLT(){};
195 
197  virtual bool select(const edm::Event& event);
199  virtual bool select(const edm::Event& event, const edm::EventSetup& setup);
200  bool selectVertex(const edm::Event& event);
201 private:
205  int min_, max_;
227 
229 
234 };
235 
237 template <typename Object>
239  src_( iC.consumes< edm::View<Object> >(cfg.getParameter<edm::InputTag>( "src" ))),
240  select_( cfg.getParameter<std::string>("select")),
241  jetIDSelect_(0)
242 {
243  // construct min/max if the corresponding params
244  // exist otherwise they are initialized with -1
245  cfg.exists("min") ? min_= cfg.getParameter<int>("min") : min_= -1;
246  cfg.exists("max") ? max_= cfg.getParameter<int>("max") : max_= -1;
247  // read electron extras if they exist
248  if(cfg.existsAs<edm::ParameterSet>("electronId")){
249  edm::ParameterSet elecId=cfg.getParameter<edm::ParameterSet>("electronId");
250  electronId_= iC.consumes< edm::ValueMap<float> >(elecId.getParameter<edm::InputTag>("src"));
251  eidPattern_= elecId.getParameter<int>("pattern");
252  }
253  // read jet corrector label if it exists
254  if(cfg.exists("jetCorrector")){ jetCorrector_= cfg.getParameter<std::string>("jetCorrector"); }
255  // read btag information if it exists
256  if(cfg.existsAs<edm::ParameterSet>("jetBTagger")){
257  edm::ParameterSet jetBTagger=cfg.getParameter<edm::ParameterSet>("jetBTagger");
258  btagLabel_= iC.consumes< reco::JetTagCollection >(jetBTagger.getParameter<edm::InputTag>("label"));
259  btagWorkingPoint_=jetBTagger.getParameter<double>("workingPoint");
260  }
261  // read jetID information if it exists
262  if(cfg.existsAs<edm::ParameterSet>("jetID")){
264  jetIDLabel_ = iC.consumes< reco::JetIDValueMap >(jetID.getParameter<edm::InputTag>("label"));
266  }
267 }
268 
270 template <typename Object>
272 {
273  // fetch input collection
275  if( !event.getByToken(src_, src) ) return false;
276 
277  // load electronId value map if configured such
279  if(!electronId_.isUninitialized()) {
280  if( !event.getByToken(electronId_, electronId) ) return false;
281  }
282 
283  // determine multiplicity of selected objects
284  int n=0;
285  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
286  // special treatment for electrons
287  if(dynamic_cast<const reco::GsfElectron*>(&*obj)){
288  unsigned int idx = obj-src->begin();
289  if( electronId_.isUninitialized() ? true : ((int)(*electronId)[src->refAt(idx)] & eidPattern_) ){
290  if(select_(*obj))++n;
291  }
292  }
293  // normal treatment
294  else{
295  if(select_(*obj))++n;
296  }
297  }
298  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
299  return (min_<0 && max_<0) ? (n>0):accept;
300 }
301 template <typename Object>
303 {
304  // fetch input collection
306  if( !event.getByToken(src_, src) ) return false;
307 
308  // load electronId value map if configured such
310  if(!electronId_.isUninitialized()) {
311  if( !event.getByToken(electronId_, electronId) ) return false;
312  }
313 
314  // determine multiplicity of selected objects
315  int n=0;
316  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
317 
318  if(select_(*obj))++n;
319  }
320  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
321  return (min_<0 && max_<0) ? (n>0):accept;
322 }
323 
324 template <typename Object>
326 {
327  throw cms::Exception("SelectionStepHLT") << "you fail" << std::endl;
328  return false;
329 }
330 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool failed_
indicate failed associations
tuple cfg
Definition: looper.py:293
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
double masslb(reco::RecoCandidate *mu, const reco::MET &met, const reco::Jet &b)
calculate mlb estimate
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate top quark mass estimate
void operator()(const std::vector< reco::Jet > &jets)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
double wMass_
paramater of the w boson mass
Base class for all types of Jets.
Definition: Jet.h:20
edm::EDGetTokenT< reco::JetTagCollection > btagLabel_
choice for b-tag as extra selection type
int min_
min/max for object multiplicity
bool accept() const
Has at least one path accepted the event?
assert(m_qm.get())
double tmassTopQuark(reco::RecoCandidate *mu, const reco::MET &met, const reco::Jet &b)
calculate top quark transverse mass estimate
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual bool select(const edm::Event &event, const edm::EventSetup &setup)
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
CalculateHLT(int maxNJets, double wMass)
default constructor
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
virtual ~SelectionStepHLT()
default destructor
SelectionStepHLT(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
default constructor
double massWBoson_
cache of w boson mass estimate
bool selectVertex(const edm::Event &event)
double mlb_
cache of mlb estimate
Definition: MET.h:42
vector< PseudoJet > jets
~CalculateHLT()
default destructor
double massTopQuark_
cache of top quark mass estimate
int j
Definition: DBlmapReader.cc:9
const int mu
Definition: Constants.h:22
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
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
virtual bool select(const edm::Event &event)
apply selection
double tmassWBoson(reco::RecoCandidate *mu, const reco::MET &met, const reco::Jet &b)
calculate W boson transverse mass estimate
edm::EDGetTokenT< edm::ValueMap< float > > electronId_
electronId label as extra selection type
StringCutObjectSelector< Object > select_
string cut selector
bool acceptHLT(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
edm::EDGetTokenT< edm::View< Object > > src_
input collection
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double tmassWBoson_
cache of W boson transverse mass estimate
virtual bool select(const edm::Event &event)
double b
Definition: hdecay.h:120
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
int maxNJets_
max. number of jets to be considered
double btagWorkingPoint_
choice of b-tag working point as extra selection type
std::string jetCorrector_
jet corrector as extra selection type
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::EDGetTokenT< reco::JetIDValueMap > jetIDLabel_
jetID as an extra selection type
double tmassTopQuark_
cache of top quark transverse mass estimate
static std::string const triggerPaths
Definition: EdmProvDump.cc:42