Templated helper class to allow a selection on a certain object collection. More...
#include <DQM/Physics/interface/TopDQMHelpers.h>
Public Member Functions | |
bool | select (const edm::Event &event) |
apply selection | |
bool | select (const edm::Event &event, const edm::EventSetup &setup) |
apply selection override for jets | |
SelectionStep (const edm::ParameterSet &cfg) | |
default constructor | |
bool | selectVertex (const edm::Event &event) |
~SelectionStep () | |
default destructor | |
Private Attributes | |
edm::InputTag | btagLabel_ |
choice for b-tag as extra selection type | |
double | btagWorkingPoint_ |
choice of b-tag working point as extra selection type | |
int | eidPattern_ |
edm::InputTag | electronId_ |
electronId label as extra selection type | |
std::string | jetCorrector_ |
jet corrector as extra selection type | |
edm::InputTag | jetIDLabel_ |
jetID as an extra selection type | |
StringCutObjectSelector < reco::JetID > * | jetIDSelect_ |
selection string on the jetID | |
int | max_ |
int | min_ |
min/max for object multiplicity | |
edm::InputTag | pvs_ |
StringCutObjectSelector< Object > | select_ |
string cut selector | |
edm::InputTag | src_ |
input collection |
Templated helper class to allow a selection on a certain object collection.
Templated helper class to allow a selection on a certain object collection, which may be monitored by a separate class afterwards. The class wraps and slightly extends the features of the StringCutParser to allow also to apply event based selections, according to a minimal or maximal number of elements in the collection after the object selection has been applied. It takes an edm::ParameterSet in the constructor, which should contain the following elements:
The parameters _src_ and _select_ are mandatory. The parameters _min_ and _max_ are optional. The parameters _electronId_ and _jetCorrector_ are optional. They are added to keep the possibility to apply selections on id'ed electrons or on corrected jets. They may be omitted in the PSet for simplification reasons if not needed at any time. They are not effiective for other object collections but electrons or jets. If none of the two parameters _min_ or _max_ is found in the event the select function returns true if at least one object fullfilled the requirements.
The class has one template value, which is the object collection to apply the selection on. This has to be parsed to the StringCutParser class. The function select is overrided for jets to circumvent problems with the template specialisation. Note that for MET not type1 or muon corrections are supported on reco candidates.
Definition at line 149 of file TopDQMHelpers.h.
SelectionStep< Object >::SelectionStep | ( | const edm::ParameterSet & | cfg | ) |
default constructor
Definition at line 198 of file TopDQMHelpers.h.
References SelectionStep< Object >::btagLabel_, SelectionStep< Object >::btagWorkingPoint_, SelectionStep< Object >::eidPattern_, SelectionStep< Object >::electronId_, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), SelectionStep< Object >::jetCorrector_, bTagSequences_cff::jetID, SelectionStep< Object >::jetIDLabel_, SelectionStep< Object >::jetIDSelect_, SelectionStep< Object >::max_, SelectionStep< Object >::min_, and AlCaHLTBitMon_QueryRunRegistry::string.
: src_( cfg.getParameter<edm::InputTag>( "src" )), select_( cfg.getParameter<std::string>("select")), jetIDSelect_(0) { // construct min/max if the corresponding params // exist otherwise they are initialized with -1 cfg.exists("min") ? min_= cfg.getParameter<int>("min") : min_= -1; cfg.exists("max") ? max_= cfg.getParameter<int>("max") : max_= -1; // read electron extras if they exist if(cfg.existsAs<edm::ParameterSet>("electronId")){ edm::ParameterSet elecId=cfg.getParameter<edm::ParameterSet>("electronId"); electronId_= elecId.getParameter<edm::InputTag>("src"); eidPattern_= elecId.getParameter<int>("pattern"); } // read jet corrector label if it exists if(cfg.exists("jetCorrector")){ jetCorrector_= cfg.getParameter<std::string>("jetCorrector"); } // read btag information if it exists if(cfg.existsAs<edm::ParameterSet>("jetBTagger")){ edm::ParameterSet jetBTagger=cfg.getParameter<edm::ParameterSet>("jetBTagger"); btagLabel_=jetBTagger.getParameter<edm::InputTag>("label"); btagWorkingPoint_=jetBTagger.getParameter<double>("workingPoint"); } // read jetID information if it exists if(cfg.existsAs<edm::ParameterSet>("jetID")){ edm::ParameterSet jetID=cfg.getParameter<edm::ParameterSet>("jetID"); jetIDLabel_ =jetID.getParameter<edm::InputTag>("label"); jetIDSelect_= new StringCutObjectSelector<reco::JetID>(jetID.getParameter<std::string>("select")); } }
SelectionStep< Object >::~SelectionStep | ( | ) | [inline] |
bool SelectionStep< Object >::select | ( | const edm::Event & | event | ) |
apply selection
Definition at line 231 of file TopDQMHelpers.h.
References accept(), edm::Event::getByLabel(), UserOptions_cff::idx, n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.
Referenced by TopDiLeptonOfflineDQM::analyze(), TopHLTSingleLeptonDQM::analyze(), and TopSingleLeptonDQM::analyze().
{ // fetch input collection edm::Handle<edm::View<Object> > src; if( !event.getByLabel(src_, src) ) return false; // load electronId value map if configured such edm::Handle<edm::ValueMap<float> > electronId; if(!electronId_.label().empty()) { if( !event.getByLabel(electronId_, electronId) ) return false; } // determine multiplicity of selected objects int n=0; for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){ // special treatment for electrons if(dynamic_cast<const reco::GsfElectron*>(&*obj)){ unsigned int idx = obj-src->begin(); if( electronId_.label().empty() ? true : ((int)(*electronId)[src->refAt(idx)] & eidPattern_) ){ if(select_(*obj))++n; } } // normal treatment else{ if(select_(*obj))++n; } } bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true); return (min_<0 && max_<0) ? (n>0):accept; }
bool SelectionStep< 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 286 of file TopDQMHelpers.h.
References accept(), JetCorrector::correction(), edm::EventSetup::find(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), UserOptions_cff::idx, metsig::jet, bTagSequences_cff::jetID, n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.
{ // fetch input collection edm::Handle<edm::View<Object> > src; if( !event.getByLabel(src_, src) ) return false; // load btag collection if configured such // NOTE that the JetTagCollection needs an // edm::View to reco::Jets; we have to add // another Handle bjets for this purpose edm::Handle<edm::View<reco::Jet> > bjets; edm::Handle<reco::JetTagCollection> btagger; edm::Handle<edm::View<reco::Vertex> > pvertex; if(!btagLabel_.label().empty()){ if( !event.getByLabel(src_, bjets) ) return false; if( !event.getByLabel(btagLabel_, btagger) ) return false; if( !event.getByLabel(pvs_, pvertex) ) return false; } // load jetID value map if configured such edm::Handle<reco::JetIDValueMap> jetID; if(jetIDSelect_){ if( !event.getByLabel(jetIDLabel_, jetID) ) return false; } // load jet corrector if configured such const JetCorrector* corrector=0; if(!jetCorrector_.empty()){ // check whether a jet correcto is in the event setup or not if(setup.find( edm::eventsetup::EventSetupRecordKey::makeKey<JetCorrectionsRecord>() )){ corrector = JetCorrector::getJetCorrector(jetCorrector_, setup); } else{ edm::LogVerbatim( "TopDQMHelpers" ) << "\n" << "------------------------------------------------------------------------------------- \n" << " No JetCorrectionsRecord available from EventSetup: \n" << " - Jets will not be corrected. \n" << " - If you want to change this add the following lines to your cfg file \n" << " \n" << " ## load jet corrections \n" << " process.load(\"JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff\") \n" << " process.prefer(\"ak5CaloL2L3\") \n" << " \n" << "------------------------------------------------------------------------------------- \n"; } } // determine multiplicity of selected objects int n=0; for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){ // check for chosen btag discriminator to be above the // corresponding working point if configured such unsigned int idx = obj-src->begin(); if( btagLabel_.label().empty() ? true : (*btagger)[bjets->refAt(idx)]>btagWorkingPoint_ ){ bool passedJetID=true; // check jetID for calo jets if( jetIDSelect_ && dynamic_cast<const reco::CaloJet*>(src->refAt(idx).get())){ passedJetID=(*jetIDSelect_)((*jetID)[src->refAt(idx)]); } if(passedJetID){ // scale jet energy if configured such Object jet=*obj; jet.scaleEnergy(corrector ? corrector->correction(*obj) : 1.); if(select_(jet))++n; } } } bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true); return (min_<0 && max_<0) ? (n>0):accept; }
bool SelectionStep< Object >::selectVertex | ( | const edm::Event & | event | ) |
Definition at line 262 of file TopDQMHelpers.h.
References accept(), edm::Event::getByLabel(), n, getGTfromDQMFile::obj, and alcazmumu_cfi::src.
Referenced by TopSingleLeptonDQM::analyze().
{ // fetch input collection edm::Handle<edm::View<Object> > src; if( !event.getByLabel(src_, src) ) return false; // load electronId value map if configured such edm::Handle<edm::ValueMap<float> > electronId; if(!electronId_.label().empty()) { if( !event.getByLabel(electronId_, electronId) ) return false; } // determine multiplicity of selected objects int n=0; for(typename edm::View<Object>::const_iterator obj=src->begin(); obj!=src->end(); ++obj){ if(select_(*obj))++n; } bool accept=(min_>=0 ? n>=min_:true) && (max_>=0 ? n<=max_:true); return (min_<0 && max_<0) ? (n>0):accept; }
edm::InputTag SelectionStep< Object >::btagLabel_ [private] |
choice for b-tag as extra selection type
Definition at line 182 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
double SelectionStep< Object >::btagWorkingPoint_ [private] |
choice of b-tag working point as extra selection type
Definition at line 184 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
int SelectionStep< 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 178 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
edm::InputTag SelectionStep< Object >::electronId_ [private] |
electronId label as extra selection type
Definition at line 167 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
std::string SelectionStep< Object >::jetCorrector_ [private] |
jet corrector as extra selection type
Definition at line 180 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
edm::InputTag SelectionStep< Object >::jetIDLabel_ [private] |
jetID as an extra selection type
Definition at line 186 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
StringCutObjectSelector<reco::JetID>* SelectionStep< Object >::jetIDSelect_ [private] |
selection string on the jetID
Definition at line 193 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
int SelectionStep< Object >::max_ [private] |
Definition at line 165 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
int SelectionStep< Object >::min_ [private] |
min/max for object multiplicity
Definition at line 165 of file TopDQMHelpers.h.
Referenced by SelectionStep< Object >::SelectionStep().
edm::InputTag SelectionStep< Object >::pvs_ [private] |
Definition at line 188 of file TopDQMHelpers.h.
StringCutObjectSelector<Object> SelectionStep< Object >::select_ [private] |
string cut selector
Definition at line 191 of file TopDQMHelpers.h.
edm::InputTag SelectionStep< Object >::src_ [private] |
input collection
Definition at line 163 of file TopDQMHelpers.h.