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
TopSingleLeptonDQM Class Reference

define MonitorEnsembple to be used More...

#include <TopSingleLeptonDQM.h>

Inheritance diagram for TopSingleLeptonDQM:
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...
 
 TopSingleLeptonDQM (const edm::ParameterSet &cfg)
 default constructor More...
 
 ~TopSingleLeptonDQM ()
 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::PFCandidate > * 
ElectronStep
 
std::vector< SelectionStep
< reco::Jet > * > 
JetSteps
 
SelectionStep< reco::MET > * METStep
 
SelectionStep
< reco::PFCandidate > * 
MuonStep
 
std::vector< SelectionStep
< reco::PFJet > * > 
PFJetSteps
 
SelectionStep< reco::Vertex > * PvStep
 
std::map< std::string,
std::pair< edm::ParameterSet,
TopSingleLepton::MonitorEnsemble * > > 
selection_
 
std::vector< std::string > selectionOrder_
 
std::vector< std::string > triggerPaths_
 trigger paths More...
 
edm::EDGetTokenT
< edm::TriggerResults
triggerTable__
 trigger table More...
 
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.

"DQM/Physics/plugins/TopSingleLeptonDQM.h"

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 278 of file TopSingleLeptonDQM.h.

Constructor & Destructor Documentation

TopSingleLeptonDQM::TopSingleLeptonDQM ( const edm::ParameterSet cfg)

default constructor

Definition at line 770 of file TopSingleLeptonDQM.cc.

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

771  : vertexSelect_(0),
772  beamspot_(""),
773  beamspotSelect_(0),
774  MuonStep(0),
775  ElectronStep(0),
776  PvStep(0),
777  METStep(0) {
778  JetSteps.clear();
779  CaloJetSteps.clear();
780  PFJetSteps.clear();
781  // configure preselection
782  edm::ParameterSet presel =
783  cfg.getParameter<edm::ParameterSet>("preselection");
784  if (presel.existsAs<edm::ParameterSet>("trigger")) {
785  edm::ParameterSet trigger =
786  presel.getParameter<edm::ParameterSet>("trigger");
787  triggerTable__ = consumes<edm::TriggerResults>(
788  trigger.getParameter<edm::InputTag>("src"));
789  triggerPaths_ = trigger.getParameter<std::vector<std::string> >("select");
790  }
791  if (presel.existsAs<edm::ParameterSet>("beamspot")) {
792  edm::ParameterSet beamspot =
793  presel.getParameter<edm::ParameterSet>("beamspot");
794  beamspot_ = beamspot.getParameter<edm::InputTag>("src");
795  beamspot__ =
796  consumes<reco::BeamSpot>(beamspot.getParameter<edm::InputTag>("src"));
798  beamspot.getParameter<std::string>("select"));
799  }
800 
801  // conifgure the selection
802  std::vector<edm::ParameterSet> sel =
803  cfg.getParameter<std::vector<edm::ParameterSet> >("selection");
804  for (unsigned int i = 0; i < sel.size(); ++i) {
805  selectionOrder_.push_back(sel.at(i).getParameter<std::string>("label"));
806  selection_[selectionStep(selectionOrder_.back())] = std::make_pair(
807  sel.at(i),
809  selectionStep(selectionOrder_.back()).c_str(),
811  }
812  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin();
813  selIt != selectionOrder_.end(); ++selIt) {
814  std::string key = selectionStep(*selIt), type = objectType(*selIt);
815  if (selection_.find(key) != selection_.end()) {
816  if (type == "muons") {
819  }
820  if (type == "elecs") {
822  selection_[key].first, consumesCollector());
823  }
824  if (type == "pvs") {
827  }
828  if (type == "jets") {
830  consumesCollector()));
831  }
832  if (type == "jets/pf") {
834  selection_[key].first, consumesCollector()));
835  }
836  if (type == "jets/calo") {
838  selection_[key].first, consumesCollector()));
839  }
840  if (type == "met") {
843  }
844  }
845  }
846 }
type
Definition: HCALResponse.h:21
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
SelectionStep< reco::MET > * METStep
SelectionStep< reco::PFCandidate > * ElectronStep
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
edm::EDGetTokenT< reco::BeamSpot > beamspot__
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * MuonStep
std::map< std::string, std::pair< edm::ParameterSet, TopSingleLepton::MonitorEnsemble * > > selection_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::InputTag beamspot_
beamspot
bool first
Definition: L1TdeRCT.cc:75
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
SelectionStep< reco::Vertex > * PvStep
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
std::string objectType(const std::string &label)
list key
Definition: combine.py:13
std::vector< SelectionStep< reco::Jet > * > JetSteps
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector
TopSingleLeptonDQM::~TopSingleLeptonDQM ( )
inline

default destructor

Definition at line 283 of file TopSingleLeptonDQM.h.

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

283  {
284  if (vertexSelect_) delete vertexSelect_;
285  if (beamspotSelect_) delete beamspotSelect_;
286  if (MuonStep) delete MuonStep;
287  if (ElectronStep) delete ElectronStep;
288  if (PvStep) delete PvStep;
289  if (METStep) delete METStep;
290  for (unsigned int i = 0; i < JetSteps.size(); i++)
291  if (JetSteps[i]) delete JetSteps[i];
292  for (unsigned int i = 0; i < CaloJetSteps.size(); i++)
293  if (CaloJetSteps[i]) delete CaloJetSteps[i];
294  for (unsigned int i = 0; i < PFJetSteps.size(); i++)
295  if (PFJetSteps[i]) delete PFJetSteps[i];
296  };
int i
Definition: DBlmapReader.cc:9
SelectionStep< reco::MET > * METStep
SelectionStep< reco::PFCandidate > * ElectronStep
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * MuonStep
SelectionStep< reco::Vertex > * PvStep
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
std::vector< SelectionStep< reco::Jet > * > JetSteps
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector

Member Function Documentation

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

do this during the event loop

Implements edm::EDAnalyzer.

Definition at line 848 of file TopSingleLeptonDQM.cc.

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

849  {
850 
853  if (!event.getByToken(triggerTable__, triggerTable)) return;
854  if (!accept(event, *triggerTable, triggerPaths_)) return;
855  }
856  if (!beamspot__.isUninitialized()) {
858  if (!event.getByToken(beamspot__, beamspot)) return;
859  if (!(*beamspotSelect_)(*beamspot)) return;
860  }
861  // cout<<" apply selection steps"<<endl;
862  unsigned int passed = 0;
863  unsigned int nJetSteps = -1;
864  unsigned int nPFJetSteps = -1;
865  unsigned int nCaloJetSteps = -1;
866  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin();
867  selIt != selectionOrder_.end(); ++selIt) {
868  std::string key = selectionStep(*selIt), type = objectType(*selIt);
869  if (selection_.find(key) != selection_.end()) {
870  if (type == "empty") {
871  selection_[key].second->fill(event, setup);
872  }
873  if (type == "muons" && MuonStep != 0) {
874  if (MuonStep->select(event)) {
875  ++passed;
876  // cout<<"selected event! "<<selection_[key].second<<endl;
877  selection_[key].second->fill(event, setup);
878  } else
879  break;
880  }
881  // cout<<" apply selection steps 2"<<endl;
882  if (type == "elecs" && ElectronStep != 0) {
883  // cout<<"In electrons ..."<<endl;
884  if (ElectronStep->select(event, "electron")) {
885  ++passed;
886  selection_[key].second->fill(event, setup);
887  } else
888  break;
889  }
890  // cout<<" apply selection steps 3"<<endl;
891  if (type == "pvs" && PvStep != 0) {
892  if (PvStep->selectVertex(event)) {
893  ++passed;
894  selection_[key].second->fill(event, setup);
895  } else
896  break;
897  }
898  // cout<<" apply selection steps 4"<<endl;
899  if (type == "jets") {
900  nJetSteps++;
901  if (JetSteps[nJetSteps] != NULL) {
902  if (JetSteps[nJetSteps]->select(event, setup)) {
903  ++passed;
904  selection_[key].second->fill(event, setup);
905  } else
906  break;
907  }
908  }
909  if (type == "jets/pf") {
910  nPFJetSteps++;
911  if (PFJetSteps[nPFJetSteps] != NULL) {
912  if (PFJetSteps[nPFJetSteps]->select(event, setup)) {
913  ++passed;
914  selection_[key].second->fill(event, setup);
915  } else
916  break;
917  }
918  }
919  if (type == "jets/calo") {
920  nCaloJetSteps++;
921  if (CaloJetSteps[nCaloJetSteps] != NULL) {
922  if (CaloJetSteps[nCaloJetSteps]->select(event, setup)) {
923  ++passed;
924  selection_[key].second->fill(event, setup);
925  } else
926  break;
927  }
928  }
929  if (type == "met" && METStep != 0) {
930  if (METStep->select(event)) {
931  ++passed;
932  selection_[key].second->fill(event, setup);
933  } else
934  break;
935  }
936  }
937  }
938 }
type
Definition: HCALResponse.h:21
bool select(const edm::Event &event)
apply selection
SelectionStep< reco::MET > * METStep
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
bool selectVertex(const edm::Event &event)
SelectionStep< reco::PFCandidate > * ElectronStep
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
#define NULL
Definition: scimark2.h:8
edm::EDGetTokenT< reco::BeamSpot > beamspot__
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:25
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
std::vector< SelectionStep< reco::PFJet > * > PFJetSteps
SelectionStep< reco::PFCandidate > * MuonStep
std::map< std::string, std::pair< edm::ParameterSet, TopSingleLepton::MonitorEnsemble * > > selection_
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
SelectionStep< reco::Vertex > * PvStep
std::vector< SelectionStep< reco::CaloJet > * > CaloJetSteps
std::string objectType(const std::string &label)
list key
Definition: combine.py:13
bool isUninitialized() const
Definition: EDGetToken.h:71
std::vector< SelectionStep< reco::Jet > * > JetSteps
std::string TopSingleLeptonDQM::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 304 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

304  {
305  return label.substr(0, label.find(':'));
306  };
std::string TopSingleLeptonDQM::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 309 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

309  {
310  return label.substr(label.find(':') + 1);
311  };

Member Data Documentation

edm::InputTag TopSingleLeptonDQM::beamspot_
private

beamspot

Definition at line 322 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

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

Definition at line 323 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

string cut selector

Definition at line 325 of file TopSingleLeptonDQM.h.

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

std::vector<SelectionStep<reco::CaloJet>*> TopSingleLeptonDQM::CaloJetSteps
private

Definition at line 343 of file TopSingleLeptonDQM.h.

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

SelectionStep<reco::PFCandidate>* TopSingleLeptonDQM::ElectronStep
private

Definition at line 339 of file TopSingleLeptonDQM.h.

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

std::vector<SelectionStep<reco::Jet>*> TopSingleLeptonDQM::JetSteps
private

Definition at line 342 of file TopSingleLeptonDQM.h.

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

SelectionStep<reco::MET>* TopSingleLeptonDQM::METStep
private

Definition at line 341 of file TopSingleLeptonDQM.h.

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

SelectionStep<reco::PFCandidate>* TopSingleLeptonDQM::MuonStep
private

Definition at line 338 of file TopSingleLeptonDQM.h.

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

std::vector<SelectionStep<reco::PFJet>*> TopSingleLeptonDQM::PFJetSteps
private

Definition at line 344 of file TopSingleLeptonDQM.h.

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

SelectionStep<reco::Vertex>* TopSingleLeptonDQM::PvStep
private

Definition at line 340 of file TopSingleLeptonDQM.h.

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

std::map<std::string, std::pair<edm::ParameterSet, TopSingleLepton::MonitorEnsemble*> > TopSingleLeptonDQM::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 337 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::vector<std::string> TopSingleLeptonDQM::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 329 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

trigger paths

Definition at line 317 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

trigger table

Definition at line 311 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

string cut selector

Definition at line 319 of file TopSingleLeptonDQM.h.

Referenced by ~TopSingleLeptonDQM().