CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
B2GHadronicHLTValidation Class Reference

#include <HLTriggerOffline/B2G/src/B2GHadronicHLTValidation.cc>

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

Public Member Functions

 B2GHadronicHLTValidation (const edm::ParameterSet &)
 
 ~B2GHadronicHLTValidation () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
std::string monitorPath (const std::string &label) const
 
void triggerBinLabels (const std::vector< std::string > &labels)
 set configurable labels for trigger monitoring histograms More...
 

Private Attributes

double etaJets_
 
MonitorElementhDenJetEta
 
MonitorElementhDenJetPt
 
MonitorElementhDenTriggerMon
 
MonitorElementhNumJetEta
 
MonitorElementhNumJetPt
 
MonitorElementhNumTriggerMon
 
double htMin_
 
bool isAll_ = false
 
bool isSel_ = false
 
edm::Ptr< reco::Jetjet_
 
unsigned int minJets_
 
double ptJets0_
 
double ptJets1_
 
double ptJets_
 
std::string sDir_
 
std::string sJets_
 
std::string sTrigger_
 
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
 
edm::EDGetTokenT< edm::TriggerResultstokTrigger_
 
std::vector< std::string > vsPaths_
 

Detailed Description

Description: compute efficiencies of trigger paths on offline reco selection with respect to pt and eta

Implementation: harvesting

Description:

Description: compute efficiencies of trigger paths on offline reco selection with respect to pt and eta

Implementation: harvesting

Definition at line 48 of file B2GHadronicHLTValidation.h.

Constructor & Destructor Documentation

B2GHadronicHLTValidation::B2GHadronicHLTValidation ( const edm::ParameterSet iConfig)
explicit

Definition at line 113 of file B2GHadronicHLTValidation.h.

References sJets_, sTrigger_, tokJets_, and tokTrigger_.

113  :
114  sDir_(iConfig.getUntrackedParameter<std::string>("sDir","HLTValidation/B2G/Efficiencies/")),
115  sJets_(iConfig.getUntrackedParameter<std::string>("sJets","ak5PFJets")),
116  ptJets_(iConfig.getUntrackedParameter<double>("ptJets",0.)),
117  ptJets0_(iConfig.getUntrackedParameter<double>("ptJets0",0.)),
118  ptJets1_(iConfig.getUntrackedParameter<double>("ptJets1",0.)),
119  etaJets_(iConfig.getUntrackedParameter<double>("etaJets",0.)),
120  minJets_(iConfig.getUntrackedParameter<unsigned int>("minJets",0)),
121  htMin_(iConfig.getUntrackedParameter<double>("htMin",0.0)),
122  sTrigger_(iConfig.getUntrackedParameter<std::string>("sTrigger","TriggerResults")),
123  vsPaths_(iConfig.getUntrackedParameter< std::vector<std::string> >("vsPaths"))
124 
125 {
126  // Jets
127  tokJets_ = consumes< edm::View<reco::Jet> >(edm::InputTag(sJets_));
128  // Trigger
129  tokTrigger_ = consumes<edm::TriggerResults>(edm::InputTag(sTrigger_, "", "HLT"));
130 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::TriggerResults > tokTrigger_
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
std::vector< std::string > vsPaths_
B2GHadronicHLTValidation::~B2GHadronicHLTValidation ( )
override

Definition at line 133 of file B2GHadronicHLTValidation.h.

References DEFINE_FWK_MODULE.

134 {
135 
136  // do anything here that needs to be done at desctruction time
137  // (e.g. close files, deallocate resources etc.)
138 
139 }

Member Function Documentation

void B2GHadronicHLTValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 44 of file B2GHadronicHLTValidation.cc.

References reco::LeafCandidate::eta(), etaJets_, MonitorElement::Fill(), edm::Event::getByToken(), hDenJetEta, hDenJetPt, hDenTriggerMon, hNumJetEta, hNumJetPt, hNumTriggerMon, htMin_, mps_fire::i, training_settings::idx, isAll_, edm::Ptr< T >::isNonnull(), isSel_, jet_, fwrapper::jets, minJets_, dataset::name, reco::LeafCandidate::pt(), ptJets0_, ptJets1_, ptJets_, tokJets_, tokTrigger_, edm::TriggerNames::triggerNames(), edm::Event::triggerNames(), and vsPaths_.

45 {
46  using namespace edm;
47 
48  isAll_ = false; isSel_ = false;
49 
50  // Jets
52  if (!iEvent.getByToken(tokJets_,jets))
53  edm::LogWarning("B2GHadronicHLTValidation") << "Jets collection not found \n";
54  unsigned int nGoodJ = 0;
55  double ht = 0.0;
56  // Check to see if we want asymmetric jet pt cuts
57  if ( ptJets0_ > 0.0 || ptJets1_ > 0.0 ) {
58  if ( ptJets0_ > 0.0 ) {
59  if ( !jets->empty() && jets->at(0).pt() > ptJets0_ ) {
60  nGoodJ++;
61  jet_ = jets->ptrAt(0) ;
62  }
63  }
64  if ( ptJets1_ > 0.0 ) {
65  if ( jets->size() > 1 && jets->at(1).pt() > ptJets1_ ) {
66  nGoodJ++;
67  jet_ = jets->ptrAt(1) ;
68  }
69  }
70  } else if (minJets_ > 0 || htMin_ > 0 ) {
71  for(edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j){
72  if (j->pt() < ptJets_) continue;
73  if (fabs(j->eta()) > etaJets_) continue;
74  nGoodJ++;
75  ht += j->pt();
76  if (nGoodJ == minJets_) jet_ = jets->ptrAt( j -jets->begin() );
77  }
78  }
79 
80  if ( nGoodJ >= minJets_ || ht > htMin_ ) isAll_ = true;
81 
82 
83  //Trigger
84  Handle<edm::TriggerResults> triggerTable;
85  if (!iEvent.getByToken(tokTrigger_,triggerTable))
86  edm::LogWarning("B2GHadronicHLTValidation") << "Trigger collection not found \n";
87  const edm::TriggerNames& triggerNames = iEvent.triggerNames(*triggerTable);
88  bool isInteresting = false;
89  for (unsigned int i=0; i<triggerNames.triggerNames().size(); ++i) {
90 
91  TString name = triggerNames.triggerNames()[i].c_str();
92  for (unsigned int j=0; j<vsPaths_.size(); j++) {
93  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
94  isInteresting = true;
95  break;
96  }
97  }
98  if (isInteresting) break;
99  }
100 
101  if (isAll_ && isInteresting) isSel_ = true;
102  else isSel_ = false;
103 
104  //Histos
105  if (isAll_) {
106  if ( jet_.isNonnull() ) {
107  hDenJetPt->Fill(jet_->pt());
108  hDenJetEta->Fill(jet_->eta());
109  }
110  for(unsigned int idx=0; idx<vsPaths_.size(); ++idx){
111  hDenTriggerMon->Fill(idx+0.5);
112  }
113 
114  }
115  if (isSel_) {
116  hNumJetPt->Fill(jet_->pt());
117  hNumJetEta->Fill(jet_->eta());
118  for (unsigned int i=0; i<triggerNames.triggerNames().size(); ++i) {
119  TString name = triggerNames.triggerNames()[i].c_str();
120  for (unsigned int j=0; j<vsPaths_.size(); j++) {
121  if (name.Contains(TString(vsPaths_[j]), TString::kIgnoreCase)) {
122  hNumTriggerMon->Fill(j+0.5 );
123  }
124  }
125  }
126  }
127 }
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
double pt() const final
transverse momentum
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
void Fill(long long x)
edm::EDGetTokenT< edm::TriggerResults > tokTrigger_
vector< PseudoJet > jets
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
edm::EDGetTokenT< edm::View< reco::Jet > > tokJets_
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
std::vector< std::string > vsPaths_
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
void B2GHadronicHLTValidation::bookHistograms ( DQMStore::IBooker dbe,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprivate

Definition at line 132 of file B2GHadronicHLTValidation.cc.

References DQMStore::IBooker::book1D(), hDenJetEta, hDenJetPt, hDenTriggerMon, hNumJetEta, hNumJetPt, hNumTriggerMon, sDir_, DQMStore::IBooker::setCurrentFolder(), triggerBinLabels(), and vsPaths_.

133 {
134  dbe.setCurrentFolder(sDir_);
135  hDenJetPt = dbe.book1D("PtLastJetAll", "PtLastJetAll", 60, 0., 3000.);
136  hDenJetEta = dbe.book1D("EtaLastJetAll", "EtaLastJetAll", 30, -3., 3.);
137  hNumJetPt = dbe.book1D("PtLastJetSel", "PtLastJetSel", 60, 0., 3000.);
138  hNumJetEta = dbe.book1D("EtaLastJetSel", "EtaLastJetSel", 30, -3., 3.);
139  // determine number of bins for trigger monitoring
140  unsigned int nPaths = vsPaths_.size();
141  // monitored trigger occupancy for single lepton triggers
142  hNumTriggerMon = dbe.book1D("TriggerMonSel" , "TriggerMonSel", nPaths, 0., nPaths);
143  hDenTriggerMon = dbe.book1D("TriggerMonAll" , "TriggerMonAll", nPaths, 0., nPaths);
144  // set bin labels for trigger monitoring
146 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
void triggerBinLabels(const std::vector< std::string > &labels)
set configurable labels for trigger monitoring histograms
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
std::vector< std::string > vsPaths_
void B2GHadronicHLTValidation::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 151 of file B2GHadronicHLTValidation.cc.

References edm::ConfigurationDescriptions::addDefault(), and edm::ParameterSetDescription::setUnknown().

151  {
152  //The following says we do not know what parameters are allowed so do no validation
153  // Please change this to state exactly what you do use, even if it is no parameters
155  desc.setUnknown();
156  descriptions.addDefault(desc);
157 }
void addDefault(ParameterSetDescription const &psetDescription)
std::string B2GHadronicHLTValidation::monitorPath ( const std::string &  label) const
inlineprivate

deduce monitorPath from label, the label is expected to be of type 'selectionPath:monitorPath'

Definition at line 61 of file B2GHadronicHLTValidation.h.

References tablePrinter::labels, and triggerBinLabels().

Referenced by triggerBinLabels().

61 { return label.substr(label.find(':')+1); };
void B2GHadronicHLTValidation::triggerBinLabels ( const std::vector< std::string > &  labels)
inlineprivate

set configurable labels for trigger monitoring histograms

Definition at line 94 of file B2GHadronicHLTValidation.h.

References hDenTriggerMon, hNumTriggerMon, training_settings::idx, monitorPath(), and MonitorElement::setBinLabel().

Referenced by bookHistograms(), and monitorPath().

95 {
96  for(unsigned int idx=0; idx<labels.size(); ++idx){
98  hDenTriggerMon->setBinLabel( idx+1, "["+monitorPath(labels[idx])+"]", 1);
99  }
100 }
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
std::string monitorPath(const std::string &label) const

Member Data Documentation

double B2GHadronicHLTValidation::etaJets_
private

Definition at line 81 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

MonitorElement* B2GHadronicHLTValidation::hDenJetEta
private

Definition at line 71 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* B2GHadronicHLTValidation::hDenJetPt
private

Definition at line 69 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* B2GHadronicHLTValidation::hDenTriggerMon
private

Definition at line 73 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), bookHistograms(), and triggerBinLabels().

MonitorElement* B2GHadronicHLTValidation::hNumJetEta
private

Definition at line 70 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* B2GHadronicHLTValidation::hNumJetPt
private

Definition at line 68 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* B2GHadronicHLTValidation::hNumTriggerMon
private

Definition at line 72 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), bookHistograms(), and triggerBinLabels().

double B2GHadronicHLTValidation::htMin_
private

Definition at line 83 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

bool B2GHadronicHLTValidation::isAll_ = false
private

Definition at line 89 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

bool B2GHadronicHLTValidation::isSel_ = false
private

Definition at line 90 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

edm::Ptr<reco::Jet> B2GHadronicHLTValidation::jet_
private

Definition at line 75 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

unsigned int B2GHadronicHLTValidation::minJets_
private

Definition at line 82 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

double B2GHadronicHLTValidation::ptJets0_
private

Definition at line 79 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

double B2GHadronicHLTValidation::ptJets1_
private

Definition at line 80 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

double B2GHadronicHLTValidation::ptJets_
private

Definition at line 78 of file B2GHadronicHLTValidation.h.

Referenced by analyze().

std::string B2GHadronicHLTValidation::sDir_
private

Definition at line 67 of file B2GHadronicHLTValidation.h.

Referenced by bookHistograms().

std::string B2GHadronicHLTValidation::sJets_
private

Definition at line 76 of file B2GHadronicHLTValidation.h.

Referenced by B2GHadronicHLTValidation().

std::string B2GHadronicHLTValidation::sTrigger_
private

Definition at line 85 of file B2GHadronicHLTValidation.h.

Referenced by B2GHadronicHLTValidation().

edm::EDGetTokenT< edm::View<reco::Jet> > B2GHadronicHLTValidation::tokJets_
private

Definition at line 77 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and B2GHadronicHLTValidation().

edm::EDGetTokenT<edm::TriggerResults> B2GHadronicHLTValidation::tokTrigger_
private

Definition at line 86 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and B2GHadronicHLTValidation().

std::vector<std::string> B2GHadronicHLTValidation::vsPaths_
private

Definition at line 87 of file B2GHadronicHLTValidation.h.

Referenced by analyze(), and bookHistograms().