CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
SelectionStepHLT< Object > Class Template Reference

#include <TopHLTDQMHelper.h>

Public Member Functions

bool select (const edm::Event &event)
 apply selection More...
 
bool select (const edm::Event &event, const edm::EventSetup &setup)
 apply selection override for jets More...
 
 SelectionStepHLT (const edm::ParameterSet &cfg)
 default constructor More...
 
bool selectVertex (const edm::Event &event)
 
 ~SelectionStepHLT ()
 default destructor More...
 

Private Attributes

edm::InputTag btagLabel_
 choice for b-tag as extra selection type More...
 
double btagWorkingPoint_
 choice of b-tag working point as extra selection type More...
 
int eidPattern_
 
edm::InputTag electronId_
 electronId label as extra selection type More...
 
std::string jetCorrector_
 jet corrector as extra selection type More...
 
edm::InputTag jetIDLabel_
 jetID as an extra selection type More...
 
StringCutObjectSelector
< reco::JetID > * 
jetIDSelect_
 selection string on the jetID More...
 
int max_
 
int min_
 min/max for object multiplicity More...
 
edm::InputTag pvs_
 
StringCutObjectSelector< Object > select_
 string cut selector More...
 
edm::InputTag src_
 input collection More...
 

Detailed Description

template<typename Object>
class SelectionStepHLT< Object >

Definition at line 173 of file TopHLTDQMHelper.h.

Constructor & Destructor Documentation

template<typename Object >
SelectionStepHLT< Object >::SelectionStepHLT ( const edm::ParameterSet cfg)

default constructor

Definition at line 222 of file TopHLTDQMHelper.h.

References SelectionStepHLT< Object >::btagLabel_, SelectionStepHLT< Object >::btagWorkingPoint_, SelectionStepHLT< Object >::eidPattern_, SelectionStepHLT< Object >::electronId_, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), SelectionStepHLT< Object >::jetCorrector_, bTagSequences_cff::jetID, SelectionStepHLT< Object >::jetIDLabel_, SelectionStepHLT< Object >::jetIDSelect_, SelectionStepHLT< Object >::max_, SelectionStepHLT< Object >::min_, and AlCaHLTBitMon_QueryRunRegistry::string.

222  :
223  src_( cfg.getParameter<edm::InputTag>( "src" )),
224  select_( cfg.getParameter<std::string>("select")),
225  jetIDSelect_(0)
226 {
227  // construct min/max if the corresponding params
228  // exist otherwise they are initialized with -1
229  cfg.exists("min") ? min_= cfg.getParameter<int>("min") : min_= -1;
230  cfg.exists("max") ? max_= cfg.getParameter<int>("max") : max_= -1;
231  // read electron extras if they exist
232  if(cfg.existsAs<edm::ParameterSet>("electronId")){
233  edm::ParameterSet elecId=cfg.getParameter<edm::ParameterSet>("electronId");
234  electronId_= elecId.getParameter<edm::InputTag>("src");
235  eidPattern_= elecId.getParameter<int>("pattern");
236  }
237  // read jet corrector label if it exists
238  if(cfg.exists("jetCorrector")){ jetCorrector_= cfg.getParameter<std::string>("jetCorrector"); }
239  // read btag information if it exists
240  if(cfg.existsAs<edm::ParameterSet>("jetBTagger")){
241  edm::ParameterSet jetBTagger=cfg.getParameter<edm::ParameterSet>("jetBTagger");
242  btagLabel_=jetBTagger.getParameter<edm::InputTag>("label");
243  btagWorkingPoint_=jetBTagger.getParameter<double>("workingPoint");
244  }
245  // read jetID information if it exists
246  if(cfg.existsAs<edm::ParameterSet>("jetID")){
248  jetIDLabel_ =jetID.getParameter<edm::InputTag>("label");
250  }
251 }
edm::InputTag btagLabel_
choice for b-tag as extra selection type
T getParameter(std::string const &) const
edm::InputTag src_
input collection
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
edm::InputTag jetIDLabel_
jetID as an extra selection type
edm::InputTag electronId_
electronId label as extra selection type
int min_
min/max for object multiplicity
bool exists(std::string const &parameterName) const
checks if a parameter exists
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
StringCutObjectSelector< Object > select_
string cut selector
double btagWorkingPoint_
choice of b-tag working point as extra selection type
std::string jetCorrector_
jet corrector as extra selection type
template<typename Object>
SelectionStepHLT< Object >::~SelectionStepHLT ( )
inline

default destructor

Definition at line 178 of file TopHLTDQMHelper.h.

178 {};

Member Function Documentation

template<typename Object >
bool SelectionStepHLT< Object >::select ( const edm::Event event)

apply selection

Definition at line 255 of file TopHLTDQMHelper.h.

References accept(), edm::Event::getByLabel(), customizeTrackingMonitorSeedNumber::idx, n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.

Referenced by TopHLTDiLeptonOfflineDQM::analyze().

256 {
257  // fetch input collection
259  if( !event.getByLabel(src_, src) ) return false;
260 
261  // load electronId value map if configured such
263  if(!electronId_.label().empty()) {
264  if( !event.getByLabel(electronId_, electronId) ) return false;
265  }
266 
267  // determine multiplicity of selected objects
268  int n=0;
269  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
270  // special treatment for electrons
271  if(dynamic_cast<const reco::GsfElectron*>(&*obj)){
272  unsigned int idx = obj-src->begin();
273  if( electronId_.label().empty() ? true : ((int)(*electronId)[src->refAt(idx)] & eidPattern_) ){
274  if(select_(*obj))++n;
275  }
276  }
277  // normal treatment
278  else{
279  if(select_(*obj))++n;
280  }
281  }
282  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
283  return (min_<0 && max_<0) ? (n>0):accept;
284 }
edm::InputTag src_
input collection
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::InputTag electronId_
electronId label as extra selection type
int min_
min/max for object multiplicity
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
StringCutObjectSelector< Object > select_
string cut selector
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::string const & label() const
Definition: InputTag.h:42
template<typename Object >
bool SelectionStepHLT< Object >::select ( const edm::Event event,
const edm::EventSetup setup 
)

apply selection override for jets

apply selection (w/o using the template class Object), override for jets

Definition at line 310 of file TopHLTDQMHelper.h.

References accept(), JetCorrector::correction(), edm::EventSetup::find(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), customizeTrackingMonitorSeedNumber::idx, metsig::jet, bTagSequences_cff::jetID, n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.

311 {
312  // fetch input collection
314  if( !event.getByLabel(src_, src) ) return false;
315 
316  // load btag collection if configured such
317  // NOTE that the JetTagCollection needs an
318  // edm::View to reco::Jets; we have to add
319  // another Handle bjets for this purpose
323  if(!btagLabel_.label().empty()){
324  if( !event.getByLabel(src_, bjets) ) return false;
325  if( !event.getByLabel(btagLabel_, btagger) ) return false;
326  if( !event.getByLabel(pvs_, pvertex) ) return false;
327  }
328 
329  // load jetID value map if configured such
331  if(jetIDSelect_){
332  if( !event.getByLabel(jetIDLabel_, jetID) ) return false;
333 
334  }
335 
336  // load jet corrector if configured such
337  const JetCorrector* corrector=0;
338  if(!jetCorrector_.empty()){
339  // check whether a jet correcto is in the event setup or not
340  if(setup.find( edm::eventsetup::EventSetupRecordKey::makeKey<JetCorrectionsRecord>() )){
341  corrector = JetCorrector::getJetCorrector(jetCorrector_, setup);
342  }
343  else{
344  edm::LogVerbatim( "TopDQMHelpers" )
345  << "\n"
346  << "------------------------------------------------------------------------------------- \n"
347  << " No JetCorrectionsRecord available from EventSetup: \n"
348  << " - Jets will not be corrected. \n"
349  << " - If you want to change this add the following lines to your cfg file \n"
350  << " \n"
351  << " ## load jet corrections \n"
352  << " process.load(\"JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff\") \n"
353  << " process.prefer(\"ak5CaloL2L3\") \n"
354  << " \n"
355  << "------------------------------------------------------------------------------------- \n";
356  }
357  }
358  // determine multiplicity of selected objects
359  int n=0;
360  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
361  // check for chosen btag discriminator to be above the
362  // corresponding working point if configured such
363  unsigned int idx = obj-src->begin();
364  if( btagLabel_.label().empty() ? true : (*btagger)[bjets->refAt(idx)]>btagWorkingPoint_ ){
365  bool passedJetID=true;
366  // check jetID for calo jets
367  if( jetIDSelect_ && dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())){
368  passedJetID=(*jetIDSelect_)((*jetID)[src->refAt(idx)]);
369  }
370  if(passedJetID){
371  // scale jet energy if configured such
372  Object jet=*obj; jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.);
373  if(select_(jet))++n;
374  }
375  }
376  }
377  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
378  return (min_<0 && max_<0) ? (n>0):accept;
379 }
edm::InputTag pvs_
edm::InputTag btagLabel_
choice for b-tag as extra selection type
edm::InputTag src_
input collection
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::InputTag jetIDLabel_
jetID as an extra selection type
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
int min_
min/max for object multiplicity
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:90
StringCutObjectSelector< reco::JetID > * jetIDSelect_
selection string on the jetID
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
StringCutObjectSelector< Object > select_
string cut selector
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:51
std::string const & label() const
Definition: InputTag.h:42
double btagWorkingPoint_
choice of b-tag working point as extra selection type
std::string jetCorrector_
jet corrector as extra selection type
template<typename Object >
bool SelectionStepHLT< Object >::selectVertex ( const edm::Event event)

Definition at line 286 of file TopHLTDQMHelper.h.

References accept(), edm::Event::getByLabel(), n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.

287 {
288  // fetch input collection
290  if( !event.getByLabel(src_, src) ) return false;
291 
292  // load electronId value map if configured such
294  if(!electronId_.label().empty()) {
295  if( !event.getByLabel(electronId_, electronId) ) return false;
296  }
297 
298  // determine multiplicity of selected objects
299  int n=0;
300  for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){
301 
302  if(select_(*obj))++n;
303  }
304  bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true);
305  return (min_<0 && max_<0) ? (n>0):accept;
306 }
edm::InputTag src_
input collection
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::InputTag electronId_
electronId label as extra selection type
int min_
min/max for object multiplicity
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
StringCutObjectSelector< Object > select_
string cut selector
std::string const & label() const
Definition: InputTag.h:42

Member Data Documentation

template<typename Object>
edm::InputTag SelectionStepHLT< Object >::btagLabel_
private

choice for b-tag as extra selection type

Definition at line 206 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
double SelectionStepHLT< Object >::btagWorkingPoint_
private

choice of b-tag working point as extra selection type

Definition at line 208 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
int SelectionStepHLT< Object >::eidPattern_
private

electronId pattern we expect the following pattern: 0: fails 1: passes electron ID only 2: passes electron Isolation only 3: passes electron ID and Isolation only 4: passes conversion rejection 5: passes conversion rejection and ID 6: passes conversion rejection and Isolation 7: passes the whole selection As described on https://twiki.cern.ch/twiki/bin/view/CMS/SimpleCutBasedEleID

Definition at line 202 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
edm::InputTag SelectionStepHLT< Object >::electronId_
private

electronId label as extra selection type

Definition at line 191 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
std::string SelectionStepHLT< Object >::jetCorrector_
private

jet corrector as extra selection type

Definition at line 204 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
edm::InputTag SelectionStepHLT< Object >::jetIDLabel_
private

jetID as an extra selection type

Definition at line 210 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
StringCutObjectSelector<reco::JetID>* SelectionStepHLT< Object >::jetIDSelect_
private

selection string on the jetID

Definition at line 217 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
int SelectionStepHLT< Object >::max_
private

Definition at line 189 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
int SelectionStepHLT< Object >::min_
private

min/max for object multiplicity

Definition at line 189 of file TopHLTDQMHelper.h.

Referenced by SelectionStepHLT< Object >::SelectionStepHLT().

template<typename Object>
edm::InputTag SelectionStepHLT< Object >::pvs_
private

Definition at line 212 of file TopHLTDQMHelper.h.

template<typename Object>
StringCutObjectSelector<Object> SelectionStepHLT< Object >::select_
private

string cut selector

Definition at line 215 of file TopHLTDQMHelper.h.

template<typename Object>
edm::InputTag SelectionStepHLT< Object >::src_
private

input collection

Definition at line 187 of file TopHLTDQMHelper.h.