CMS 3D CMS Logo

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

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &setup)
 do this during the event loop
 TopSingleLeptonDQM (const edm::ParameterSet &cfg)
 default constructor
 ~TopSingleLeptonDQM ()
 default destructor

Private Member Functions

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

Private Attributes

edm::InputTag beamspot_
 beamspot
StringCutObjectSelector
< reco::BeamSpot > * 
beamspotSelect_
 string cut selector
std::map< std::string,
std::pair< edm::ParameterSet,
TopSingleLepton::MonitorEnsemble * > > 
selection_
std::vector< std::string > selectionOrder_
std::vector< std::string > triggerPaths_
 trigger paths
edm::InputTag triggerTable_
 trigger table
edm::InputTag vertex_
 primary vertex
StringCutObjectSelector
< reco::Vertex > * 
vertexSelect_
 string cut selector

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(), AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, triggerTable_, vertex_, and vertexSelect_.

                                                                : triggerTable_(""), vertexSelect_(0), beamspot_(""), beamspotSelect_(0)
{
  // configure preselection
  edm::ParameterSet presel=cfg.getParameter<edm::ParameterSet>("preselection");
  if( presel.existsAs<edm::ParameterSet>("trigger") ){
    edm::ParameterSet trigger=presel.getParameter<edm::ParameterSet>("trigger");
    triggerTable_=trigger.getParameter<edm::InputTag>("src");
    triggerPaths_=trigger.getParameter<std::vector<std::string> >("select");
  } 
  if( presel.existsAs<edm::ParameterSet>("vertex" ) ){
    edm::ParameterSet vertex=presel.getParameter<edm::ParameterSet>("vertex");
    vertex_= vertex.getParameter<edm::InputTag>("src");
    vertexSelect_= new StringCutObjectSelector<reco::Vertex>(vertex.getParameter<std::string>("select"));
  }
  if( presel.existsAs<edm::ParameterSet>("beamspot" ) ){
    edm::ParameterSet beamspot=presel.getParameter<edm::ParameterSet>("beamspot");
    beamspot_= beamspot.getParameter<edm::InputTag>("src");
    beamspotSelect_= new StringCutObjectSelector<reco::BeamSpot>(beamspot.getParameter<std::string>("select"));
  }

  // conifgure the selection
  std::vector<edm::ParameterSet> sel=cfg.getParameter<std::vector<edm::ParameterSet> >("selection");
  for(unsigned int i=0; i<sel.size(); ++i){
    selectionOrder_.push_back(sel.at(i).getParameter<std::string>("label"));
    selection_[selectionStep(selectionOrder_.back())] = std::make_pair(sel.at(i), new TopSingleLepton::MonitorEnsemble(selectionStep(selectionOrder_.back()).c_str(), cfg.getParameter<edm::ParameterSet>("setup")));
  }
}
TopSingleLeptonDQM::~TopSingleLeptonDQM ( ) [inline]

default destructor

Definition at line 217 of file TopSingleLeptonDQM.h.

References beamspotSelect_, and vertexSelect_.

                       {
    if( vertexSelect_ ) delete vertexSelect_;
    if( beamspotSelect_ ) delete beamspotSelect_;
  };

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(), relval_parameters_module::step, AlCaHLTBitMon_QueryRunRegistry::string, triggerPaths_, and triggerTable_.

{ 
  if(!triggerTable_.label().empty()){
    edm::Handle<edm::TriggerResults> triggerTable;
    if( !event.getByLabel(triggerTable_, triggerTable) ) return;
    if(!accept(event, *triggerTable, triggerPaths_)) return;
  }
  //cout<<"trig passed"<<endl;
  if(!beamspot_.label().empty()){
    edm::Handle<reco::BeamSpot> beamspot;
    if( !event.getByLabel(beamspot_, beamspot) ) return;
    if(!(*beamspotSelect_)(*beamspot)) return;
  }
  // apply selection steps
  unsigned int passed=0;
  for(std::vector<std::string>::const_iterator selIt=selectionOrder_.begin(); selIt!=selectionOrder_.end(); ++selIt){
    std::string key = selectionStep(*selIt), type = objectType(*selIt);
    if(selection_.find(key)!=selection_.end()){
      if(type=="empty"){
        selection_[key].second->fill(event, setup);
      }
      if(type=="muons"){
        SelectionStep<reco::Muon> step(selection_[key].first);
        if(step.select(event)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="elecs"){
        SelectionStep<reco::GsfElectron> step(selection_[key].first);
        if(step.select(event)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="pvs" ){
        SelectionStep<reco::Vertex> step(selection_[key].first);
        if(step.selectVertex(event)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="jets" ){
        SelectionStep<reco::Jet> step(selection_[key].first);
        if(step.select(event, setup)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="jets/pf" ){
        SelectionStep<reco::PFJet> step(selection_[key].first);
        if(step.select(event, setup)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="jets/calo" ){
        SelectionStep<reco::CaloJet> step(selection_[key].first);
        if(step.select(event, setup)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
      if(type=="met" ){
        SelectionStep<reco::MET> step(selection_[key].first);
        if(step.select(event)){ ++passed;
          selection_[key].second->fill(event, setup);
        } else break;
      }
    }
  }
}
std::string TopSingleLeptonDQM::objectType ( const std::string &  label) [inline, private]

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().

{ return label.substr(0, label.find(':')); };  
std::string TopSingleLeptonDQM::selectionStep ( const std::string &  label) [inline, private]

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().

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

Member Data Documentation

beamspot

Definition at line 244 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

string cut selector

Definition at line 246 of file TopSingleLeptonDQM.h.

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

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().

trigger table

Definition at line 231 of file TopSingleLeptonDQM.h.

Referenced by analyze(), and TopSingleLeptonDQM().

primary vertex

Definition at line 239 of file TopSingleLeptonDQM.h.

Referenced by TopSingleLeptonDQM().

string cut selector

Definition at line 241 of file TopSingleLeptonDQM.h.

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