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 
146 template <typename Object>
148 public:
150  SelectionStep(const edm::ParameterSet& cfg);
153 
155  bool select(const edm::Event& event);
157  bool select(const edm::Event& event, const edm::EventSetup& setup);
158 
159 private:
163  int min_, max_;
167  std::string jetCorrector_;
174 
179 };
180 
182 template <typename Object>
184  src_( cfg.getParameter<edm::InputTag>( "src" )),
185  select_( cfg.getParameter<std::string>("select")),
186  jetIDSelect_(0)
187 {
188  // construct min/max if the corresponding params
189  // exist otherwise they are initialized with -1
190  cfg.exists("min") ? min_= cfg.getParameter<int>("min") : min_= -1;
191  cfg.exists("max") ? max_= cfg.getParameter<int>("max") : max_= -1;
192  // read electron extras if they exist
193  if(cfg.exists("electronId" )){ electronId_= cfg.getParameter<edm::InputTag>("electronId" ); }
194  // read jet corrector label if it exists
195  if(cfg.exists("jetCorrector")){ jetCorrector_= cfg.getParameter<std::string>("jetCorrector"); }
196  // read btag information if it exists
197  if(cfg.existsAs<edm::ParameterSet>("jetBTagger")){
198  edm::ParameterSet jetBTagger=cfg.getParameter<edm::ParameterSet>("jetBTagger");
199  btagLabel_=jetBTagger.getParameter<edm::InputTag>("label");
200  btagWorkingPoint_=jetBTagger.getParameter<double>("workingPoint");
201  }
202  // read jetID information if it exists
203  if(cfg.existsAs<edm::ParameterSet>("jetID")){
205  jetIDLabel_ =jetID.getParameter<edm::InputTag>("label");
206  jetIDSelect_= new StringCutObjectSelector<reco::JetID>(jetID.getParameter<std::string>("select"));
207  }
208 }
209 
211 template <typename Object>
213 {
214  // fetch input collection
216  if( !event.getByLabel(src_, src) ) return false;
217 
218  // load electronId value map if configured such
220  if(!electronId_.label().empty()) {
221  if( !event.getByLabel(electronId_, electronId) ) return false;
222  }
223 
224  // determine multiplicity of selected objects
225  int n=0;
226  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
227  // special treatment for electrons
228  if(dynamic_cast<const reco::GsfElectron*>(&*obj)){
229  unsigned int idx = obj-src->begin();
230  if( electronId_.label().empty() ? true : (*electronId)[src->refAt(idx)]>0.99 ){
231  if(select_(*obj))++n;
232  }
233  }
234  // normal treatment
235  else{
236  if(select_(*obj))++n;
237  }
238  }
239  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
240  return (min_<0 && max_<0) ? (n>0):accept;
241 }
242 
244 template <typename Object>
246 {
247  // fetch input collection
249  if( !event.getByLabel(src_, src) ) return false;
250 
251  // load btag collection if configured such
252  // NOTE that the JetTagCollection needs an
253  // edm::View to reco::Jets; we have to add
254  // another Handle bjets for this purpose
257  if(!btagLabel_.label().empty()){
258  if( !event.getByLabel(src_, bjets) ) return false;
259  if( !event.getByLabel(btagLabel_, btagger) ) return false;
260  }
261 
262  // load jetID value map if configured such
264  if(jetIDSelect_){
265  if( !event.getByLabel(jetIDLabel_, jetID) ) return false;
266  }
267 
268  // load jet corrector if configured such
269  const JetCorrector* corrector=0;
270  if(!jetCorrector_.empty()){
271  // check whether a jet correcto is in the event setup or not
272  if(setup.find( edm::eventsetup::EventSetupRecordKey::makeKey<JetCorrectionsRecord>() )){
273  corrector = JetCorrector::getJetCorrector(jetCorrector_, setup);
274  }
275  else{
276  edm::LogVerbatim( "TopDQMHelpers" )
277  << "\n"
278  << "------------------------------------------------------------------------------------- \n"
279  << " No JetCorrectionsRecord available from EventSetup: \n"
280  << " - Jets will not be corrected. \n"
281  << " - If you want to change this add the following lines to your cfg file \n"
282  << " \n"
283  << " ## load jet corrections \n"
284  << " process.load(\"JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff\") \n"
285  << " process.prefer(\"ak5CaloL2L3\") \n"
286  << " \n"
287  << "------------------------------------------------------------------------------------- \n";
288  }
289  }
290  // determine multiplicity of selected objects
291  int n=0;
292  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
293  // check for chosen btag discriminator to be above the
294  // corresponding working point if configured such
295  unsigned int idx = obj-src->begin();
296  if( btagLabel_.label().empty() ? true : (*btagger)[bjets->refAt(idx)]>btagWorkingPoint_ ){
297  bool passedJetID=true;
298  // check jetID for calo jets
299  if( jetIDSelect_ && dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())){
300  passedJetID=(*jetIDSelect_)((*jetID)[src->refAt(idx)]);
301  }
302  if(passedJetID){
303  // scale jet energy if configured such
304  Object jet=*obj; jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.);
305  if(select_(jet))++n;
306  }
307  }
308  }
309  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
310  return (min_<0 && max_<0) ? (n>0):accept;
311 }
312 
313 #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:180
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:82
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:359
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 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
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