CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
TopSingleLeptonDQM Class Reference

define MonitorEnsembple to be used More...

#include <TopSingleLeptonDQM.h>

Inheritance diagram for TopSingleLeptonDQM:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 do this during the event loop More...
 
 TopSingleLeptonDQM (const edm::ParameterSet &cfg)
 default constructor More...
 
 ~TopSingleLeptonDQM ()
 default destructor More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

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__
 
std::unique_ptr< StringCutObjectSelector< reco::BeamSpot > > beamspotSelect_
 string cut selector More...
 
std::vector< std::unique_ptr< SelectionStep< reco::CaloJet > > > CaloJetSteps
 
std::unique_ptr< SelectionStep< reco::PFCandidate > > ElectronStep
 
std::vector< std::unique_ptr< SelectionStep< reco::Jet > > > JetSteps
 
std::unique_ptr< SelectionStep< reco::MET > > METStep
 
std::unique_ptr< SelectionStep< reco::PFCandidate > > MuonStep
 
std::vector< std::unique_ptr< SelectionStep< reco::PFJet > > > PFJetSteps
 
std::unique_ptr< SelectionStep< reco::Vertex > > PvStep
 
std::vector< edm::ParameterSetsel_
 
std::map< std::string, std::pair< edm::ParameterSet, std::unique_ptr< TopSingleLepton::MonitorEnsemble > > > selection_
 
std::vector< std::string > selectionOrder_
 
edm::ParameterSet setup_
 
std::vector< std::string > triggerPaths_
 trigger paths More...
 
edm::EDGetTokenT< edm::TriggerResultstriggerTable__
 trigger table More...
 
std::unique_ptr< StringCutObjectSelector< reco::Vertex > > vertexSelect_
 string cut selector More...
 

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

Constructor & Destructor Documentation

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

default constructor

Definition at line 680 of file TopSingleLeptonDQM.cc.

References beamspot_, beamspot__, beamspotSelect_, CaloJetSteps, ElectronStep, edm::ParameterSet::existsAs(), plotBeamSpotDB::first, edm::ParameterSet::getParameter(), mps_fire::i, JetSteps, crabWrapper::key, METStep, MuonStep, objectType(), PFJetSteps, PvStep, sel_, selection_, selectionOrder_, selectionStep(), setup_, AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, and triggerTable__.

681  : vertexSelect_(nullptr),
682  beamspot_(""),
683  beamspotSelect_(nullptr),
684  MuonStep(nullptr),
685  ElectronStep(nullptr),
686  PvStep(nullptr),
687  METStep(nullptr) {
688  JetSteps.clear();
689  CaloJetSteps.clear();
690  PFJetSteps.clear();
691  // configure preselection
692  edm::ParameterSet presel =
693  cfg.getParameter<edm::ParameterSet>("preselection");
694  if (presel.existsAs<edm::ParameterSet>("trigger")) {
696  presel.getParameter<edm::ParameterSet>("trigger");
697  triggerTable__ = consumes<edm::TriggerResults>(
698  trigger.getParameter<edm::InputTag>("src"));
699  triggerPaths_ = trigger.getParameter<std::vector<std::string> >("select");
700  }
701  if (presel.existsAs<edm::ParameterSet>("beamspot")) {
703  presel.getParameter<edm::ParameterSet>("beamspot");
704  beamspot_ = beamspot.getParameter<edm::InputTag>("src");
705  beamspot__ =
706  consumes<reco::BeamSpot>(beamspot.getParameter<edm::InputTag>("src"));
708  beamspot.getParameter<std::string>("select")));
709  }
710 
711  // configure the selection
712  sel_ = cfg.getParameter<std::vector<edm::ParameterSet> >("selection");
713  setup_ = cfg.getParameter<edm::ParameterSet>("setup");
714  for (unsigned int i = 0; i < sel_.size(); ++i) {
715  selectionOrder_.push_back(sel_.at(i).getParameter<std::string>("label"));
716  selection_[selectionStep(selectionOrder_.back())] = std::make_pair(
717  sel_.at(i),
718  std::unique_ptr<TopSingleLepton::MonitorEnsemble>(
720  selectionStep(selectionOrder_.back()).c_str(),
721  setup_, consumesCollector())));
722  }
723  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin();
724  selIt != selectionOrder_.end(); ++selIt) {
725  std::string key = selectionStep(*selIt), type = objectType(*selIt);
726  if (selection_.find(key) != selection_.end()) {
727  if (type == "muons") {
729  consumesCollector()));
730  }
731  if (type == "elecs") {
733  selection_[key].first, consumesCollector()));
734  }
735  if (type == "pvs") {
736  PvStep.reset(new SelectionStep<reco::Vertex>(selection_[key].first,
737  consumesCollector()));
738  }
739  if (type == "jets") {
740  JetSteps.push_back(std::unique_ptr<SelectionStep<reco::Jet>>(
741  new SelectionStep<reco::Jet>(selection_[key].first, consumesCollector())));
742  }
743  if (type == "jets/pf") {
744  PFJetSteps.push_back(std::unique_ptr<SelectionStep<reco::PFJet>>(
745  new SelectionStep<reco::PFJet>(selection_[key].first, consumesCollector())));
746  }
747  if (type == "jets/calo") {
748  CaloJetSteps.push_back(std::unique_ptr<SelectionStep<reco::CaloJet>>(
749  new SelectionStep<reco::CaloJet>(selection_[key].first, consumesCollector())));
750  }
751  if (type == "met") {
752  METStep.reset(new SelectionStep<reco::MET>(selection_[key].first,
753  consumesCollector()));
754  }
755  }
756  }
757 }
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:161
std::unique_ptr< SelectionStep< reco::PFCandidate > > ElectronStep
std::unique_ptr< StringCutObjectSelector< reco::BeamSpot > > beamspotSelect_
string cut selector
edm::EDGetTokenT< reco::BeamSpot > beamspot__
std::map< std::string, std::pair< edm::ParameterSet, std::unique_ptr< TopSingleLepton::MonitorEnsemble > > > selection_
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
std::vector< std::unique_ptr< SelectionStep< reco::Jet > > > JetSteps
std::unique_ptr< StringCutObjectSelector< reco::Vertex > > vertexSelect_
string cut selector
edm::InputTag beamspot_
beamspot
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
std::unique_ptr< SelectionStep< reco::PFCandidate > > MuonStep
std::unique_ptr< SelectionStep< reco::MET > > METStep
Templated helper class to allow a selection on a certain object collection.
std::string objectType(const std::string &label)
std::vector< edm::ParameterSet > sel_
edm::ParameterSet setup_
std::vector< std::unique_ptr< SelectionStep< reco::PFJet > > > PFJetSteps
std::unique_ptr< SelectionStep< reco::Vertex > > PvStep
std::vector< std::unique_ptr< SelectionStep< reco::CaloJet > > > CaloJetSteps
TopSingleLeptonDQM::~TopSingleLeptonDQM ( )
inline

default destructor

Definition at line 294 of file TopSingleLeptonDQM.h.

References analyze(), bookHistograms(), and GeneralSetup::setup().

294 {};

Member Function Documentation

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

do this during the event loop

Definition at line 765 of file TopSingleLeptonDQM.cc.

References accept(), egammaObjectModificationsInMiniAOD_cff::beamspot, beamspot__, beamspotSelect_, CaloJetSteps, ElectronStep, edm::Event::getByToken(), edm::EDGetTokenT< T >::isUninitialized(), JetSteps, crabWrapper::key, METStep, MuonStep, objectType(), TriggerAnalyzer::passed, PFJetSteps, PvStep, TauGenJetsDecayModeSelectorAllHadrons_cfi::select, selection_, selectionOrder_, selectionStep(), AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, and triggerTable__.

766  {
767 
770  if (!event.getByToken(triggerTable__, triggerTable)) return;
771  if (!accept(event, *triggerTable, triggerPaths_)) return;
772  }
773  if (!beamspot__.isUninitialized()) {
775  if (!event.getByToken(beamspot__, beamspot)) return;
776  if (!(*beamspotSelect_)(*beamspot)) return;
777  }
778  // cout<<" apply selection steps"<<endl;
779  unsigned int passed = 0;
780  unsigned int nJetSteps = -1;
781  unsigned int nPFJetSteps = -1;
782  unsigned int nCaloJetSteps = -1;
783  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin();
784  selIt != selectionOrder_.end(); ++selIt) {
785  std::string key = selectionStep(*selIt), type = objectType(*selIt);
786  if (selection_.find(key) != selection_.end()) {
787  if (type == "empty") {
788  selection_[key].second->fill(event, setup);
789  }
790  if (type == "muons" && MuonStep != nullptr) {
791  if (MuonStep->select(event)) {
792  ++passed;
793  // cout<<"selected event! "<<selection_[key].second<<endl;
794  selection_[key].second->fill(event, setup);
795  } else
796  break;
797  }
798  // cout<<" apply selection steps 2"<<endl;
799  if (type == "elecs" && ElectronStep != nullptr) {
800  // cout<<"In electrons ..."<<endl;
801  if (ElectronStep->select(event, "electron")) {
802  ++passed;
803  selection_[key].second->fill(event, setup);
804  } else
805  break;
806  }
807  // cout<<" apply selection steps 3"<<endl;
808  if (type == "pvs" && PvStep != nullptr) {
809  if (PvStep->selectVertex(event)) {
810  ++passed;
811  selection_[key].second->fill(event, setup);
812  } else
813  break;
814  }
815  // cout<<" apply selection steps 4"<<endl;
816  if (type == "jets") {
817  nJetSteps++;
818  if (JetSteps[nJetSteps] != nullptr) {
819  if (JetSteps[nJetSteps]->select(event, setup)) {
820  ++passed;
821  selection_[key].second->fill(event, setup);
822  } else
823  break;
824  }
825  }
826  if (type == "jets/pf") {
827  nPFJetSteps++;
828  if (PFJetSteps[nPFJetSteps] != nullptr) {
829  if (PFJetSteps[nPFJetSteps]->select(event, setup)) {
830  ++passed;
831  selection_[key].second->fill(event, setup);
832  } else
833  break;
834  }
835  }
836  if (type == "jets/calo") {
837  nCaloJetSteps++;
838  if (CaloJetSteps[nCaloJetSteps] != nullptr) {
839  if (CaloJetSteps[nCaloJetSteps]->select(event, setup)) {
840  ++passed;
841  selection_[key].second->fill(event, setup);
842  } else
843  break;
844  }
845  }
846  if (type == "met" && METStep != nullptr) {
847  if (METStep->select(event)) {
848  ++passed;
849  selection_[key].second->fill(event, setup);
850  } else
851  break;
852  }
853  }
854  }
855 }
type
Definition: HCALResponse.h:21
std::unique_ptr< SelectionStep< reco::PFCandidate > > ElectronStep
std::unique_ptr< StringCutObjectSelector< reco::BeamSpot > > beamspotSelect_
string cut selector
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< reco::BeamSpot > beamspot__
std::map< std::string, std::pair< edm::ParameterSet, std::unique_ptr< TopSingleLepton::MonitorEnsemble > > > selection_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
std::vector< std::unique_ptr< SelectionStep< reco::Jet > > > JetSteps
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
std::unique_ptr< SelectionStep< reco::PFCandidate > > MuonStep
std::unique_ptr< SelectionStep< reco::MET > > METStep
std::string objectType(const std::string &label)
bool isUninitialized() const
Definition: EDGetToken.h:70
std::vector< std::unique_ptr< SelectionStep< reco::PFJet > > > PFJetSteps
std::unique_ptr< SelectionStep< reco::Vertex > > PvStep
std::vector< std::unique_ptr< SelectionStep< reco::CaloJet > > > CaloJetSteps
void TopSingleLeptonDQM::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotected

Definition at line 758 of file TopSingleLeptonDQM.cc.

References selection_.

759  {
760 
761  for (auto selIt = selection_.begin(); selIt != selection_.end(); ++selIt) {
762  selIt->second.second->book(ibooker);
763  }
764 }
std::map< std::string, std::pair< edm::ParameterSet, std::unique_ptr< TopSingleLepton::MonitorEnsemble > > > selection_
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 307 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

307  {
308  return label.substr(0, label.find(':'));
309  };
char const * label
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 312 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

312  {
313  return label.substr(label.find(':') + 1);
314  };
char const * label

Member Data Documentation

edm::InputTag TopSingleLeptonDQM::beamspot_
private

beamspot

Definition at line 325 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

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

Definition at line 326 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<StringCutObjectSelector<reco::BeamSpot> > TopSingleLeptonDQM::beamspotSelect_
private

string cut selector

Definition at line 328 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

Definition at line 346 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<SelectionStep<reco::PFCandidate> > TopSingleLeptonDQM::ElectronStep
private

Definition at line 342 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

Definition at line 345 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<SelectionStep<reco::MET> > TopSingleLeptonDQM::METStep
private

Definition at line 344 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<SelectionStep<reco::PFCandidate> > TopSingleLeptonDQM::MuonStep
private

Definition at line 341 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

Definition at line 347 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<SelectionStep<reco::Vertex> > TopSingleLeptonDQM::PvStep
private

Definition at line 343 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::vector<edm::ParameterSet> TopSingleLeptonDQM::sel_
private

Definition at line 349 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

std::map<std::string, std::pair<edm::ParameterSet, std::unique_ptr<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 340 of file TopSingleLeptonDQM.h.

Referenced by analyze(), bookHistograms(), 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 332 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

edm::ParameterSet TopSingleLeptonDQM::setup_
private

Definition at line 350 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

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

trigger paths

Definition at line 320 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

trigger table

Definition at line 314 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

std::unique_ptr<StringCutObjectSelector<reco::Vertex> > TopSingleLeptonDQM::vertexSelect_
private

string cut selector

Definition at line 322 of file TopSingleLeptonDQM.h.