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 "DQM/Physics/plugins/TopSingleLeptonDQM.h"

Inheritance diagram for TopSingleLeptonDQM:
edm::EDAnalyzer

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
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

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

Private Attributes

edm::InputTag beamspot_
 beamspot More...
 
StringCutObjectSelector
< reco::BeamSpot > * 
beamspotSelect_
 string cut selector More...
 
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::InputTag triggerTable_
 trigger table More...
 
edm::InputTag vertex_
 primary vertex More...
 
StringCutObjectSelector
< reco::Vertex > * 
vertexSelect_
 string cut selector More...
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

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

Constructor & Destructor Documentation

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

default constructor

Definition at line 530 of file TopSingleLeptonDQM.cc.

References beamspot_, beamspotSelect_, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), i, EgammaValidation_Wenu_cff::sel, selection_, selectionOrder_, selectionStep(), triggerPaths_, triggerTable_, vertex_, and vertexSelect_.

531 {
532  // configure preselection
533  edm::ParameterSet presel=cfg.getParameter<edm::ParameterSet>("preselection");
534  if( presel.existsAs<edm::ParameterSet>("trigger") ){
535  edm::ParameterSet trigger=presel.getParameter<edm::ParameterSet>("trigger");
536  triggerTable_=trigger.getParameter<edm::InputTag>("src");
537  triggerPaths_=trigger.getParameter<std::vector<std::string> >("select");
538  }
539  if( presel.existsAs<edm::ParameterSet>("vertex" ) ){
540  edm::ParameterSet vertex=presel.getParameter<edm::ParameterSet>("vertex");
541  vertex_= vertex.getParameter<edm::InputTag>("src");
542  vertexSelect_= new StringCutObjectSelector<reco::Vertex>(vertex.getParameter<std::string>("select"));
543  }
544  if( presel.existsAs<edm::ParameterSet>("beamspot" ) ){
545  edm::ParameterSet beamspot=presel.getParameter<edm::ParameterSet>("beamspot");
546  beamspot_= beamspot.getParameter<edm::InputTag>("src");
547  beamspotSelect_= new StringCutObjectSelector<reco::BeamSpot>(beamspot.getParameter<std::string>("select"));
548  }
549 
550  // conifgure the selection
551  std::vector<edm::ParameterSet> sel=cfg.getParameter<std::vector<edm::ParameterSet> >("selection");
552  for(unsigned int i=0; i<sel.size(); ++i){
553  selectionOrder_.push_back(sel.at(i).getParameter<std::string>("label"));
554  selection_[selectionStep(selectionOrder_.back())] = std::make_pair(sel.at(i), new TopSingleLepton::MonitorEnsemble(selectionStep(selectionOrder_.back()).c_str(), cfg.getParameter<edm::ParameterSet>("setup")));
555  }
556 }
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:187
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
edm::InputTag vertex_
primary vertex
std::map< std::string, std::pair< edm::ParameterSet, TopSingleLepton::MonitorEnsemble * > > selection_
edm::InputTag beamspot_
beamspot
edm::InputTag triggerTable_
trigger table
StringCutObjectSelector< reco::Vertex > * vertexSelect_
string cut selector
TopSingleLeptonDQM::~TopSingleLeptonDQM ( )
inline

default destructor

Definition at line 217 of file TopSingleLeptonDQM.h.

References beamspotSelect_, and vertexSelect_.

217  {
218  if( vertexSelect_ ) delete vertexSelect_;
219  if( beamspotSelect_ ) delete beamspotSelect_;
220  };
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
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 559 of file TopSingleLeptonDQM.cc.

References accept(), beamspot_, beamspotSelect_, first, edm::Event::getByLabel(), combine::key, edm::InputTag::label(), objectType(), SelectionStep< Object >::select(), selection_, selectionOrder_, selectionStep(), SelectionStep< Object >::selectVertex(), launcher::step, triggerPaths_, and triggerTable_.

560 {
561  if(!triggerTable_.label().empty()){
563  if( !event.getByLabel(triggerTable_, triggerTable) ) return;
564  if(!accept(event, *triggerTable, triggerPaths_)) return;
565  }
566  //cout<<"trig passed"<<endl;
567  if(!beamspot_.label().empty()){
569  if( !event.getByLabel(beamspot_, beamspot) ) return;
570  if(!(*beamspotSelect_)(*beamspot)) return;
571  }
572  // apply selection steps
573  unsigned int passed=0;
574  for(std::vector<std::string>::const_iterator selIt=selectionOrder_.begin(); selIt!=selectionOrder_.end(); ++selIt){
575  std::string key = selectionStep(*selIt), type = objectType(*selIt);
576  if(selection_.find(key)!=selection_.end()){
577  if(type=="empty"){
578  selection_[key].second->fill(event, setup);
579  }
580  if(type=="muons"){
582  if(step.select(event)){ ++passed;
583  selection_[key].second->fill(event, setup);
584  } else break;
585  }
586  if(type=="elecs"){
588  if(step.select(event)){ ++passed;
589  selection_[key].second->fill(event, setup);
590  } else break;
591  }
592  if(type=="pvs" ){
594  if(step.selectVertex(event)){ ++passed;
595  selection_[key].second->fill(event, setup);
596  } else break;
597  }
598  if(type=="jets" ){
600  if(step.select(event, setup)){ ++passed;
601  selection_[key].second->fill(event, setup);
602  } else break;
603  }
604  if(type=="jets/pf" ){
606  if(step.select(event, setup)){ ++passed;
607  selection_[key].second->fill(event, setup);
608  } else break;
609  }
610  if(type=="jets/calo" ){
612  if(step.select(event, setup)){ ++passed;
613  selection_[key].second->fill(event, setup);
614  } else break;
615  }
616  if(type=="met" ){
618  if(step.select(event)){ ++passed;
619  selection_[key].second->fill(event, setup);
620  } else break;
621  }
622  }
623  }
624 }
type
Definition: HCALResponse.h:22
list step
Definition: launcher.py:15
StringCutObjectSelector< reco::BeamSpot > * beamspotSelect_
string cut selector
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
std::vector< std::string > triggerPaths_
trigger paths
std::vector< std::string > selectionOrder_
std::string selectionStep(const std::string &label)
std::map< std::string, std::pair< edm::ParameterSet, TopSingleLepton::MonitorEnsemble * > > selection_
edm::InputTag beamspot_
beamspot
edm::InputTag triggerTable_
trigger table
bool first
Definition: L1TdeRCT.cc:94
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
Templated helper class to allow a selection on a certain object collection.
std::string const & label() const
Definition: InputTag.h:25
std::string objectType(const std::string &label)
list key
Definition: combine.py:13
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 228 of file TopSingleLeptonDQM.h.

Referenced by analyze().

228 { return label.substr(0, label.find(':')); };
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 231 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

Member Data Documentation

edm::InputTag TopSingleLeptonDQM::beamspot_
private

beamspot

Definition at line 244 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

string cut selector

Definition at line 246 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 256 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 250 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

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

trigger paths

Definition at line 237 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

edm::InputTag TopSingleLeptonDQM::triggerTable_
private

trigger table

Definition at line 231 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

edm::InputTag TopSingleLeptonDQM::vertex_
private

primary vertex

Definition at line 239 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

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

string cut selector

Definition at line 241 of file TopSingleLeptonDQM.h.

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