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 Member Functions | Private Attributes
SingleTopTChannelLeptonDQM Class Reference

define MonitorEnsembple to be used More...

#include "DQM/Physics/plugins/SingleTopTChannelLeptonDQM.h"

Inheritance diagram for SingleTopTChannelLeptonDQM:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &setup)
 do this during the event loop More...
 
 SingleTopTChannelLeptonDQM (const edm::ParameterSet &cfg)
 default constructor More...
 
 ~SingleTopTChannelLeptonDQM ()
 default destructor More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

std::string objectType (const std::string &label)
 
std::string selectionStep (const std::string &label)
 

Private Attributes

edm::InputTag beamspot_
 beamspot More...
 
edm::EDGetTokenT< reco::BeamSpotbeamspot__
 
StringCutObjectSelector
< reco::BeamSpot > * 
beamspotSelect_
 string cut selector More...
 
std::vector< SelectionStep
< reco::CaloJet > * > 
CaloJetSteps
 
SelectionStep
< reco::GsfElectron > * 
ElectronStep
 
std::vector< SelectionStep
< reco::Jet > * > 
JetSteps
 
SelectionStep< reco::MET > * METStep
 
SelectionStep< reco::Muon > * MuonStep
 
SelectionStep
< reco::PFCandidate > * 
PFElectronStep
 
std::vector< SelectionStep
< reco::PFJet > * > 
PFJetSteps
 
SelectionStep
< reco::PFCandidate > * 
PFMuonStep
 
SelectionStep< reco::Vertex > * PvStep
 
std::map< std::string,
std::pair< edm::ParameterSet,
SingleTopTChannelLepton::MonitorEnsemble * > > 
selection_
 
std::vector< std::string > selectionOrder_
 
std::vector< std::string > triggerPaths_
 trigger paths More...
 
edm::EDGetTokenT
< edm::TriggerResults
triggerTable__
 trigger table More...
 
edm::InputTag vertex_
 primary vertex More...
 
edm::EDGetTokenT< reco::Vertexvertex__
 
StringCutObjectSelector
< reco::Vertex > * 
vertexSelect_
 string cut selector More...
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

define MonitorEnsembple to be used

Module to apply a monitored selection of top like events in the semi-leptonic channel.

Plugin to apply a monitored selection of top like events with some minimal flexibility in the number and definition of the selection steps. To achieve this flexibility it employes the SelectionStep class. The MonitorEnsemble class is used to provide a well defined set of histograms to be monitored after each selection step. The SelectionStep class provides a flexible and intuitive selection via the StringCutParser. SelectionStep and MonitorEnsemble classes are interleaved. The monitoring starts after a preselection step (which is not monitored in the context of this module) with an instance of the MonitorEnsemble class. The following objects are supported for selection:

These types have to be present as prefix of the selection step paramter label separated from the rest of the label by a ':' (e.g. in the form "jets:step0"). The class expects selection labels of this type. They will be disentangled by the private helper functions objectType and seletionStep as declared below.

Definition at line 235 of file SingleTopTChannelLeptonDQM.h.

Constructor & Destructor Documentation

SingleTopTChannelLeptonDQM::SingleTopTChannelLeptonDQM ( const edm::ParameterSet cfg)

default constructor

Definition at line 797 of file SingleTopTChannelLeptonDQM.cc.

References beamspot_, beamspot__, beamspotSelect_, CaloJetSteps, edm::EDConsumerBase::consumesCollector(), ElectronStep, edm::ParameterSet::existsAs(), first, edm::ParameterSet::getParameter(), i, JetSteps, combine::key, METStep, MuonStep, objectType(), PFElectronStep, PFJetSteps, PFMuonStep, PvStep, EgammaValidation_Wenu_cff::sel, selection_, selectionOrder_, selectionStep(), AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, triggerTable__, vertex_, vertex__, and vertexSelect_.

799 
800 {
801  JetSteps.clear();
802  CaloJetSteps.clear();
803  PFJetSteps.clear();
804 
805  // configure preselection
806  edm::ParameterSet presel=cfg.getParameter<edm::ParameterSet>("preselection");
807  if( presel.existsAs<edm::ParameterSet>("trigger") ){
808  edm::ParameterSet trigger=presel.getParameter<edm::ParameterSet>("trigger");
809  triggerTable__=consumes<edm::TriggerResults>(trigger.getParameter<edm::InputTag>("src"));
810  triggerPaths_=trigger.getParameter<std::vector<std::string> >("select");
811  }
812  if( presel.existsAs<edm::ParameterSet>("vertex" ) ){
813  edm::ParameterSet vertex=presel.getParameter<edm::ParameterSet>("vertex");
814  vertex_= vertex.getParameter<edm::InputTag>("src");
815  vertex__= consumes<reco::Vertex>(vertex.getParameter<edm::InputTag>("src"));
816  vertexSelect_= new StringCutObjectSelector<reco::Vertex>(vertex.getParameter<std::string>("select"));
817  }
818  if( presel.existsAs<edm::ParameterSet>("beamspot" ) ){
819  edm::ParameterSet beamspot=presel.getParameter<edm::ParameterSet>("beamspot");
820  beamspot_= beamspot.getParameter<edm::InputTag>("src");
821  beamspot__= consumes<reco::BeamSpot>(beamspot.getParameter<edm::InputTag>("src"));
822  beamspotSelect_= new StringCutObjectSelector<reco::BeamSpot>(beamspot.getParameter<std::string>("select"));
823  }
824  // conifgure the selection
825  std::vector<edm::ParameterSet> sel=cfg.getParameter<std::vector<edm::ParameterSet> >("selection");
826  for(unsigned int i=0; i<sel.size(); ++i){
827  selectionOrder_.push_back(sel.at(i).getParameter<std::string>("label"));
828  selection_[selectionStep(selectionOrder_.back())] = std::make_pair(sel.at(i), new SingleTopTChannelLepton::MonitorEnsemble(selectionStep(selectionOrder_.back()).c_str(), cfg.getParameter<edm::ParameterSet>("setup"), cfg.getParameter<std::vector<edm::ParameterSet> >("selection"), consumesCollector()));
829  }
830  for(std::vector<std::string>::const_iterator selIt=selectionOrder_.begin(); selIt!=selectionOrder_.end(); ++selIt){
831  std::string key = selectionStep(*selIt), type = objectType(*selIt);
832  if(selection_.find(key)!=selection_.end()){
833  if(type=="muons"){
835  }
836  if(type=="muons/pf"){
838  }
839  if(type=="elecs"){
841  }
842  if(type=="elecs/pf"){
844  }
845  if(type=="pvs"){
847  }
848  if(type=="jets" ){
850  }
851  if(type=="jets/pf" ){
853  }
854  if(type=="jets/calo" ){
856  }
857  if(type=="met"){
859  }
860  }
861  }
862 }
std::vector< std::string > selectionOrder_
type
Definition: HCALResponse.h:21
std::string objectType(const std::string &label)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
std::vector< std::string > triggerPaths_
trigger paths
SelectionStep< reco::MET > * METStep
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
std::map< std::string, std::pair< edm::ParameterSet, SingleTopTChannelLepton::MonitorEnsemble * > > selection_
edm::InputTag vertex_
primary vertex
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector
SelectionStep< reco::Vertex > * PvStep
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
SelectionStep< reco::GsfElectron > * ElectronStep
bool first
Definition: L1TdeRCT.cc:75
std::vector< SelectionStep< reco::Jet > * > JetSteps
std::string selectionStep(const std::string &label)
SelectionStep< reco::Muon > * MuonStep
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * PFMuonStep
list key
Definition: combine.py:13
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
edm::EDGetTokenT< reco::Vertex > vertex__
edm::EDGetTokenT< reco::BeamSpot > beamspot__
SelectionStep< reco::PFCandidate > * PFElectronStep
SingleTopTChannelLeptonDQM::~SingleTopTChannelLeptonDQM ( )
inline

default destructor

Definition at line 240 of file SingleTopTChannelLeptonDQM.h.

References beamspotSelect_, CaloJetSteps, ElectronStep, i, JetSteps, METStep, MuonStep, PFElectronStep, PFJetSteps, PFMuonStep, PvStep, and vertexSelect_.

240  {
241  if( vertexSelect_ ) delete vertexSelect_;
242  if( beamspotSelect_ ) delete beamspotSelect_;
243  if( MuonStep) delete MuonStep;
244  if( PFMuonStep) delete PFMuonStep;
245  if( ElectronStep) delete ElectronStep;
246  if( PFElectronStep) delete PFElectronStep;
247  if( PvStep) delete PvStep;
248  for(unsigned int i = 0; i < JetSteps.size(); i++)
249  if( JetSteps[i]) delete JetSteps[i];
250  for(unsigned int i = 0; i < CaloJetSteps.size(); i++)
251  if( CaloJetSteps[i]) delete CaloJetSteps[i];
252  for(unsigned int i = 0; i < PFJetSteps.size(); i++)
253  if( PFJetSteps[i]) delete PFJetSteps[i];
254  if( METStep) delete METStep;
255 
256 
257 
258 
259  };
int i
Definition: DBlmapReader.cc:9
SelectionStep< reco::MET > * METStep
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector
SelectionStep< reco::Vertex > * PvStep
SelectionStep< reco::GsfElectron > * ElectronStep
std::vector< SelectionStep< reco::Jet > * > JetSteps
SelectionStep< reco::Muon > * MuonStep
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * PFMuonStep
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
SelectionStep< reco::PFCandidate > * PFElectronStep

Member Function Documentation

void SingleTopTChannelLeptonDQM::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
virtual

do this during the event loop

Implements edm::EDAnalyzer.

Definition at line 866 of file SingleTopTChannelLeptonDQM.cc.

References accept(), beamspot__, beamspotSelect_, CaloJetSteps, ElectronStep, edm::Event::getByToken(), edm::EDGetTokenT< T >::isUninitialized(), JetSteps, combine::key, METStep, MuonStep, NULL, objectType(), PFElectronStep, PFJetSteps, PFMuonStep, benchmark_cfg::select, SelectionStep< Object >::select(), selection_, selectionOrder_, selectionStep(), AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, triggerTable__, vertex__, and vertexSelect_.

867 {
870  if(!event.getByToken(triggerTable__, triggerTable) ) return;
871  if(!accept(event, *triggerTable, triggerPaths_)) return;
872  }
875  if( !event.getByToken(beamspot__, beamspot) ) return;
876  if(!(*beamspotSelect_)(*beamspot)) return;
877  }
878 
879  if(!vertex__.isUninitialized()){
881  if( !event.getByToken(vertex__, vertex) ) return;
882  edm::View<reco::Vertex>::const_iterator pv = vertex->begin();
883  //if ((pv->isFake()) || (pv->ndof() < 4) || (abs(pv->z())>24.) || (pv->position().Rho() > 2.0))
884  if(!(*vertexSelect_)(*pv)) return;
885  }
886 
887 
888  // apply selection steps
889  unsigned int passed=0;
890  unsigned int nJetSteps = -1;
891  unsigned int nPFJetSteps = -1;
892  unsigned int nCaloJetSteps = -1;
893  for(std::vector<std::string>::const_iterator selIt=selectionOrder_.begin(); selIt!=selectionOrder_.end(); ++selIt){
894  std::string key = selectionStep(*selIt), type = objectType(*selIt);
895  if(selection_.find(key)!=selection_.end()){
896  if(type=="empty"){
897  selection_[key].second->fill(event, setup);
898  }
899  if(type=="presel" ){
900  selection_[key].second->fill(event, setup);
901  }
902  if(type=="elecs" && ElectronStep != 0){
903  if(ElectronStep->select(event)){ ++passed;
904  selection_[key].second->fill(event, setup);
905  } else break;
906  }
907  if(type=="elecs/pf" && PFElectronStep != 0){
908 
909  if(PFElectronStep->select(event, "electron")){ ++passed;
910 
911  selection_[key].second->fill(event, setup);
912 
913  } else break;
914  }
915  if(type=="muons" && MuonStep != 0){
916  if(MuonStep->select(event)){ ++passed;
917  selection_[key].second->fill(event, setup);
918  } else break;
919  }
920  if(type=="muons/pf" && PFMuonStep != 0){
921  // cout << "MUON SELECTION" << endl;
922  if(PFMuonStep->select(event, "muon")){ ++passed;
923  selection_[key].second->fill(event, setup);
924  } else break;
925  }
926  if(type=="jets" ){
927  nJetSteps++;
928  if(JetSteps[nJetSteps] != NULL){
929  if(JetSteps[nJetSteps]->select(event, setup)){ ++passed;
930  selection_[key].second->fill(event, setup);
931  } else break;
932  }
933  }
934  if(type=="jets/pf" ){
935  // cout << "JET SELECTION" << endl;
936  nPFJetSteps++;
937  if(PFJetSteps[nPFJetSteps] != NULL){
938  if(PFJetSteps[nPFJetSteps]->select(event, setup)){ ++passed;
939  selection_[key].second->fill(event, setup);
940  }else break;
941  }
942  }
943  if(type=="jets/calo" ){
944  nCaloJetSteps++;
945  if(CaloJetSteps[nCaloJetSteps] != NULL){
946  if(CaloJetSteps[nCaloJetSteps]->select(event, setup)){ ++passed;
947  selection_[key].second->fill(event, setup);
948  } else break;
949  }
950  }
951  if(type=="met" && METStep != 0 ){
952  if(METStep->select(event)){ ++passed;
953  selection_[key].second->fill(event, setup);
954  } else break;
955  }
956  }
957  }
958 }
std::vector< std::string > selectionOrder_
type
Definition: HCALResponse.h:21
std::string objectType(const std::string &label)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
bool select(const edm::Event &event)
apply selection
std::vector< std::string > triggerPaths_
trigger paths
SelectionStep< reco::MET > * METStep
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
#define NULL
Definition: scimark2.h:8
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:26
std::map< std::string, std::pair< edm::ParameterSet, SingleTopTChannelLepton::MonitorEnsemble * > > selection_
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector
SelectionStep< reco::GsfElectron > * ElectronStep
std::vector< SelectionStep< reco::Jet > * > JetSteps
std::string selectionStep(const std::string &label)
SelectionStep< reco::Muon > * MuonStep
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * PFMuonStep
list key
Definition: combine.py:13
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
edm::EDGetTokenT< reco::Vertex > vertex__
edm::EDGetTokenT< reco::BeamSpot > beamspot__
bool isUninitialized() const
Definition: EDGetToken.h:71
SelectionStep< reco::PFCandidate > * PFElectronStep
std::string SingleTopTChannelLeptonDQM::objectType ( const std::string &  label)
inlineprivate

deduce object type from ParameterSet label, the label is expected to be of type 'objectType:selectionStep'

Definition at line 269 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

269 { return label.substr(0, label.find(':')); };
std::string SingleTopTChannelLeptonDQM::selectionStep ( const std::string &  label)
inlineprivate

deduce selection step from ParameterSet label, the label is expected to be of type 'objectType:selectionStep'

Definition at line 272 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

272 { return label.substr(label.find(':')+1); };

Member Data Documentation

edm::InputTag SingleTopTChannelLeptonDQM::beamspot_
private

beamspot

Definition at line 288 of file SingleTopTChannelLeptonDQM.h.

Referenced by SingleTopTChannelLeptonDQM().

edm::EDGetTokenT<reco::BeamSpot> SingleTopTChannelLeptonDQM::beamspot__
private

Definition at line 289 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

StringCutObjectSelector<reco::BeamSpot>* SingleTopTChannelLeptonDQM::beamspotSelect_
private

string cut selector

Definition at line 291 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), SingleTopTChannelLeptonDQM(), and ~SingleTopTChannelLeptonDQM().

std::vector<SelectionStep<reco::CaloJet> * > SingleTopTChannelLeptonDQM::CaloJetSteps
private
SelectionStep<reco::GsfElectron>* SingleTopTChannelLeptonDQM::ElectronStep
private
std::vector<SelectionStep<reco::Jet> * > SingleTopTChannelLeptonDQM::JetSteps
private
SelectionStep<reco::MET>* SingleTopTChannelLeptonDQM::METStep
private
SelectionStep<reco::Muon>* SingleTopTChannelLeptonDQM::MuonStep
private
SelectionStep<reco::PFCandidate>* SingleTopTChannelLeptonDQM::PFElectronStep
private
std::vector<SelectionStep<reco::PFJet> * > SingleTopTChannelLeptonDQM::PFJetSteps
private
SelectionStep<reco::PFCandidate>* SingleTopTChannelLeptonDQM::PFMuonStep
private
SelectionStep<reco::Vertex>* SingleTopTChannelLeptonDQM::PvStep
private
std::map<std::string, std::pair<edm::ParameterSet, SingleTopTChannelLepton::MonitorEnsemble*> > SingleTopTChannelLeptonDQM::selection_
private

this is the heart component of the plugin; std::string keeps a label the selection step for later identification, edm::ParameterSet keeps the configuration of the selection for the SelectionStep class, MonitoringEnsemble keeps an instance of the MonitorEnsemble class to be filled after each selection step

Definition at line 301 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

std::vector<std::string> SingleTopTChannelLeptonDQM::selectionOrder_
private

needed to guarantee the selection order as defined by the order of ParameterSets in the selection vector as defined in the config

Definition at line 295 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

std::vector<std::string> SingleTopTChannelLeptonDQM::triggerPaths_
private

trigger paths

Definition at line 280 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

edm::EDGetTokenT<edm::TriggerResults> SingleTopTChannelLeptonDQM::triggerTable__
private

trigger table

Definition at line 272 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

edm::InputTag SingleTopTChannelLeptonDQM::vertex_
private

primary vertex

Definition at line 282 of file SingleTopTChannelLeptonDQM.h.

Referenced by SingleTopTChannelLeptonDQM().

edm::EDGetTokenT<reco::Vertex> SingleTopTChannelLeptonDQM::vertex__
private

Definition at line 283 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), and SingleTopTChannelLeptonDQM().

StringCutObjectSelector<reco::Vertex>* SingleTopTChannelLeptonDQM::vertexSelect_
private

string cut selector

Definition at line 285 of file SingleTopTChannelLeptonDQM.h.

Referenced by analyze(), SingleTopTChannelLeptonDQM(), and ~SingleTopTChannelLeptonDQM().