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 
10 
20 inline bool
21 accept(const edm::Event& event, const edm::TriggerResults& triggerTable, const std::string& triggerPath)
22 {
23  bool passed=false;
24  const edm::TriggerNames& triggerNames = event.triggerNames(triggerTable);
25  for(unsigned int i=0; i<triggerNames.triggerNames().size(); ++i){
26  if(triggerNames.triggerNames()[i] == triggerPath) {
27  if(triggerTable.accept(i)){
28  passed=true;
29  break;
30  }
31  }
32  }
33  return passed;
34 }
35 
36 inline bool
37 accept(const edm::Event& event, const edm::TriggerResults& triggerTable, const std::vector<std::string>& triggerPaths)
38 {
39  bool passed=false;
40  for(unsigned int j=0; j<triggerPaths.size(); ++j){
41  if(accept(event, triggerTable, triggerPaths[j])){
42  passed=true;
43  break;
44  }
45  }
46  return passed;
47 }
48 
49 
52 
64 class Calculate {
65  public:
67  Calculate(int maxNJets, double wMass);
70 
72  double massWBoson(const std::vector<reco::Jet>& jets);
74  double massTopQuark(const std::vector<reco::Jet>& jets);
75 
76  private:
80  void operator()(const std::vector<reco::Jet>& jets);
81 
82  private:
84  bool failed_;
86  int maxNJets_;
88  double wMass_;
90  double massWBoson_;
92  double massTopQuark_;
93 };
94 
95 
106 
147 template <typename Object>
149 public:
151  SelectionStep(const edm::ParameterSet& cfg);
154 
156  bool select(const edm::Event& event);
158  bool select(const edm::Event& event, const edm::EventSetup& setup);
159 
160 private:
164  int min_, max_;
179  std::string jetCorrector_;
186 
191 };
192 
194 template <typename Object>
196  src_( cfg.getParameter<edm::InputTag>( "src" )),
197  select_( cfg.getParameter<std::string>("select")),
198  jetIDSelect_(0)
199 {
200  // construct min/max if the corresponding params
201  // exist otherwise they are initialized with -1
202  cfg.exists("min") ? min_= cfg.getParameter<int>("min") : min_= -1;
203  cfg.exists("max") ? max_= cfg.getParameter<int>("max") : max_= -1;
204  // read electron extras if they exist
205  if(cfg.existsAs<edm::ParameterSet>("electronId")){
206  edm::ParameterSet elecId=cfg.getParameter<edm::ParameterSet>("electronId");
207  electronId_= elecId.getParameter<edm::InputTag>("src");
208  eidPattern_= elecId.getParameter<int>("pattern");
209  }
210  // read jet corrector label if it exists
211  if(cfg.exists("jetCorrector")){ jetCorrector_= cfg.getParameter<std::string>("jetCorrector"); }
212  // read btag information if it exists
213  if(cfg.existsAs<edm::ParameterSet>("jetBTagger")){
214  edm::ParameterSet jetBTagger=cfg.getParameter<edm::ParameterSet>("jetBTagger");
215  btagLabel_=jetBTagger.getParameter<edm::InputTag>("label");
216  btagWorkingPoint_=jetBTagger.getParameter<double>("workingPoint");
217  }
218  // read jetID information if it exists
219  if(cfg.existsAs<edm::ParameterSet>("jetID")){
221  jetIDLabel_ =jetID.getParameter<edm::InputTag>("label");
222  jetIDSelect_= new StringCutObjectSelector<reco::JetID>(jetID.getParameter<std::string>("select"));
223  }
224 }
225 
227 template <typename Object>
229 {
230  // fetch input collection
232  if( !event.getByLabel(src_, src) ) return false;
233 
234  // load electronId value map if configured such
236  if(!electronId_.label().empty()) {
237  if( !event.getByLabel(electronId_, electronId) ) return false;
238  }
239 
240  // determine multiplicity of selected objects
241  int n=0;
242  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
243  // special treatment for electrons
244  if(dynamic_cast<const reco::GsfElectron*>(&*obj)){
245  unsigned int idx = obj-src->begin();
246  if( electronId_.label().empty() ? true : ((int)(*electronId)[src->refAt(idx)] & eidPattern_) ){
247  if(select_(*obj))++n;
248  }
249  }
250  // normal treatment
251  else{
252  if(select_(*obj))++n;
253  }
254  }
255  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
256  return (min_<0 && max_<0) ? (n>0):accept;
257 }
258 
260 template <typename Object>
262 {
263  // fetch input collection
265  if( !event.getByLabel(src_, src) ) return false;
266 
267  // load btag collection if configured such
268  // NOTE that the JetTagCollection needs an
269  // edm::View to reco::Jets; we have to add
270  // another Handle bjets for this purpose
273  if(!btagLabel_.label().empty()){
274  if( !event.getByLabel(src_, bjets) ) return false;
275  if( !event.getByLabel(btagLabel_, btagger) ) return false;
276  }
277 
278  // load jetID value map if configured such
280  if(jetIDSelect_){
281  if( !event.getByLabel(jetIDLabel_, jetID) ) return false;
282  }
283 
284  // load jet corrector if configured such
285  const JetCorrector* corrector=0;
286  if(!jetCorrector_.empty()){
287  // check whether a jet correcto is in the event setup or not
288  if(setup.find( edm::eventsetup::EventSetupRecordKey::makeKey<JetCorrectionsRecord>() )){
289  corrector = JetCorrector::getJetCorrector(jetCorrector_, setup);
290  }
291  else{
292  edm::LogVerbatim( "TopDQMHelpers" )
293  << "\n"
294  << "------------------------------------------------------------------------------------- \n"
295  << " No JetCorrectionsRecord available from EventSetup: \n"
296  << " - Jets will not be corrected. \n"
297  << " - If you want to change this add the following lines to your cfg file \n"
298  << " \n"
299  << " ## load jet corrections \n"
300  << " process.load(\"JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff\") \n"
301  << " process.prefer(\"ak5CaloL2L3\") \n"
302  << " \n"
303  << "------------------------------------------------------------------------------------- \n";
304  }
305  }
306  // determine multiplicity of selected objects
307  int n=0;
308  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
309  // check for chosen btag discriminator to be above the
310  // corresponding working point if configured such
311  unsigned int idx = obj-src->begin();
312  if( btagLabel_.label().empty() ? true : (*btagger)[bjets->refAt(idx)]>btagWorkingPoint_ ){
313  bool passedJetID=true;
314  // check jetID for calo jets
315  if( jetIDSelect_ && dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())){
316  passedJetID=(*jetIDSelect_)((*jetID)[src->refAt(idx)]);
317  }
318  if(passedJetID){
319  // scale jet energy if configured such
320  Object jet=*obj; jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.);
321  if(select_(jet))++n;
322  }
323  }
324  }
325  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
326  return (min_<0 && max_<0) ? (n>0):accept;
327 }
328 
329 #endif
double massWBoson_
cache of w boson mass estimate
Definition: TopDQMHelpers.h:90
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:181
bool select(const edm::Event &event)
apply selection
edm::InputTag src_
input collection
~Calculate()
default destructor
Definition: TopDQMHelpers.h:69
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
double massTopQuark_
cache of top quark mass estimate
Definition: TopDQMHelpers.h:92
bool accept() const
Has at least one path accepted the event?
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:21
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
Definition: TopDQMHelpers.cc:9
~SelectionStep()
default destructor
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:90
edm::InputTag jetIDLabel_
jetID as an extra selection type
Helper class for the calculation of a top and a W boson mass estime.
Definition: TopDQMHelpers.h:64
int min_
min/max for object multiplicity
tuple obj
Example code starts here #.
Definition: VarParsing.py:655
int maxNJets_
max. number of jets to be considered
Definition: TopDQMHelpers.h:86
int j
Definition: DBlmapReader.cc:9
double wMass_
paramater of the w boson mass
Definition: TopDQMHelpers.h:88
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
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
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:40
edm::InputTag electronId_
electronId label as extra selection type
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)
bool failed_
indicate failed associations
Definition: TopDQMHelpers.h:84
Calculate(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:3
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
std::string jetCorrector_
jet corrector as extra selection type
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
edm::InputTag btagLabel_
choice for b-tag as extra selection type
StringCutObjectSelector< Object > select_
string cut selector
SelectionStep(const edm::ParameterSet &cfg)
default constructor
tuple src
Definition: align_tpl.py:87