CMS 3D CMS Logo

ObjMonitor.cc
Go to the documentation of this file.
2 
4 
6 
8 
13 
14 // -----------------------------
15 // constructors and destructor
16 // -----------------------------
17 
19  folderName_ ( iConfig.getParameter<std::string>("FolderName") )
20  , metToken_ ( consumes<reco::PFMETCollection> (iConfig.getParameter<edm::InputTag>("met") ) )
21  , jetToken_ ( mayConsume<reco::PFJetCollection> (iConfig.getParameter<edm::InputTag>("jets") ) )
22  , eleToken_ ( mayConsume<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electrons") ) )
23  , muoToken_ ( mayConsume<reco::MuonCollection> (iConfig.getParameter<edm::InputTag>("muons") ) )
24  , do_met_ (iConfig.getParameter<bool>("doMETHistos") )
25  , do_jet_ (iConfig.getParameter<bool>("doJetHistos") )
26  , do_ht_ (iConfig.getParameter<bool>("doHTHistos") )
27  , num_genTriggerEventFlag_(std::make_unique<GenericTriggerEventFlag>(iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"),consumesCollector(), *this))
28  , den_genTriggerEventFlag_(std::make_unique<GenericTriggerEventFlag>(iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"),consumesCollector(), *this))
29  , metSelection_ ( iConfig.getParameter<std::string>("metSelection") )
30  , jetSelection_ ( iConfig.getParameter<std::string>("jetSelection") )
31  , jetId_ ( iConfig.getParameter<std::string>("jetId") )
32  , htjetSelection_ ( iConfig.getParameter<std::string>("htjetSelection"))
33  , eleSelection_ ( iConfig.getParameter<std::string>("eleSelection") )
34  , muoSelection_ ( iConfig.getParameter<std::string>("muoSelection") )
35  , njets_ ( iConfig.getParameter<int>("njets" ) )
36  , nelectrons_ ( iConfig.getParameter<int>("nelectrons" ) )
37  , nmuons_ ( iConfig.getParameter<int>("nmuons" ) )
38 {
39  if (do_met_){
40  metDQM_.initialise(iConfig);
41  }
42  if (do_jet_){
43  jetDQM_.initialise(iConfig);
44  }
45  if (do_ht_ ){
46  htDQM_.initialise(iConfig);
47  }
48 }
49 
50 ObjMonitor::~ObjMonitor() = default;
51 
53  edm::Run const & iRun,
54  edm::EventSetup const & iSetup)
55 {
56 
57  std::string currentFolder = folderName_ ;
58  ibooker.setCurrentFolder(currentFolder.c_str());
59 
60  if (do_met_) metDQM_.bookHistograms(ibooker);
61  if (do_jet_) jetDQM_.bookHistograms(ibooker);
62  if (do_ht_ ) htDQM_.bookHistograms(ibooker);
63 
64  // Initialize the GenericTriggerEventFlag
65  if ( num_genTriggerEventFlag_ && num_genTriggerEventFlag_->on() ) num_genTriggerEventFlag_->initRun( iRun, iSetup );
66  if ( den_genTriggerEventFlag_ && den_genTriggerEventFlag_->on() ) den_genTriggerEventFlag_->initRun( iRun, iSetup );
67 
68 }
69 
71 
72  // Filter out events if Trigger Filtering is requested
73  if (den_genTriggerEventFlag_->on() && ! den_genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
74 
76  iEvent.getByToken( metToken_, metHandle );
77  reco::PFMET pfmet = metHandle->front();
78  if ( ! metSelection_( pfmet ) ) return;
79 
80  float met = pfmet.pt();
81  float phi = pfmet.phi();
82 
84  iEvent.getByToken( jetToken_, jetHandle );
85  std::vector<reco::PFJet> jets;
86  std::vector<reco::PFJet> htjets;
87  if ( int(jetHandle->size()) < njets_ ) return;
88  for ( auto const & j : *jetHandle ) {
89  if ( jetSelection_( j ) ) {
90  if (jetId_=="loose" || jetId_ =="tight"){
91  double abseta = abs(j.eta());
92  double NHF = j.neutralHadronEnergyFraction();
93  double NEMF = j.neutralEmEnergyFraction();
94  double CHF = j.chargedHadronEnergyFraction();
95  double CEMF = j.chargedEmEnergyFraction();
96  unsigned NumNeutralParticles =j.neutralMultiplicity();
97  unsigned CHM = j.chargedMultiplicity();
98  bool passId = (jetId_=="loose" && looseJetId(abseta,NHF,NEMF,CHF,CEMF,NumNeutralParticles,CHM)) || (jetId_=="tight" && tightJetId(abseta,NHF,NEMF,CHF,CEMF,NumNeutralParticles,CHM));
99  if (passId) jets.push_back(j);
100  }
101  else jets.push_back(j);
102  }
103  if ( htjetSelection_( j ) ) htjets.push_back(j);
104  }
105  if ( int(jets.size()) < njets_ ) return;
106 
108  iEvent.getByToken( eleToken_, eleHandle );
109  std::vector<reco::GsfElectron> electrons;
110  if ( int(eleHandle->size()) < nelectrons_ ) return;
111  for ( auto const & e : *eleHandle ) {
112  if ( eleSelection_( e ) ) electrons.push_back(e);
113  }
114  if ( int(electrons.size()) < nelectrons_ ) return;
115 
117  iEvent.getByToken( muoToken_, muoHandle );
118  if ( int(muoHandle->size()) < nmuons_ ) return;
119  std::vector<reco::Muon> muons;
120  for ( auto const & m : *muoHandle ) {
121  if ( muoSelection_( m ) ) muons.push_back(m);
122  }
123  if ( int(muons.size()) < nmuons_ ) return;
124 
125 
126  bool passNumCond = num_genTriggerEventFlag_->off() || num_genTriggerEventFlag_->accept( iEvent, iSetup);
127  int ls = iEvent.id().luminosityBlock();
128 
129  if (do_met_) metDQM_.fillHistograms(met,phi,ls,passNumCond);
130  if (do_jet_) jetDQM_.fillHistograms(jets,pfmet,ls,passNumCond);
131  if (do_ht_ ) htDQM_.fillHistograms(htjets,met,ls,passNumCond);
132 
133 }
134 
135 bool ObjMonitor::looseJetId(const double & abseta,
136  const double & NHF,
137  const double & NEMF,
138  const double & CHF,
139  const double & CEMF,
140  const unsigned & NumNeutralParticles,
141  const unsigned & CHM)
142 {
143  if (abseta<=2.7){
144  unsigned NumConst = CHM+NumNeutralParticles;
145 
146  return ((NumConst>1 && NHF<0.99 && NEMF<0.99) && ((abseta<=2.4 && CHF>0 && CHM>0 && CEMF<0.99) || abseta>2.4));
147  }
148  else if (abseta<=3){
149  return (NumNeutralParticles>2 && NEMF>0.01 && NHF<0.98);
150  }
151  else {
152  return NumNeutralParticles>10 && NEMF<0.90;
153  }
154 }
155 bool ObjMonitor::tightJetId(const double & abseta,
156  const double & NHF,
157  const double & NEMF,
158  const double & CHF,
159  const double & CEMF,
160  const unsigned & NumNeutralParticles,
161  const unsigned & CHM)
162 {
163  if (abseta<=2.7){
164  unsigned NumConst = CHM+NumNeutralParticles;
165  return (NumConst>1 && NHF<0.90 && NEMF<0.90 ) && ((abseta<=2.4 && CHF>0 && CHM>0 && CEMF<0.99) || abseta>2.4);
166  }
167  else if (abseta<=3){
168  return (NHF<0.98 && NEMF>0.01 && NumNeutralParticles>2);
169  }
170  else {
171  return (NEMF<0.90 && NumNeutralParticles>10);
172  }
173 }
174 
176 {
178  desc.add<std::string> ( "FolderName", "HLT/OBJ" );
179 
180  desc.add<edm::InputTag>( "met", edm::InputTag("pfMet") );
181  desc.add<edm::InputTag>( "jets", edm::InputTag("ak4PFJetsCHS") );
182  desc.add<edm::InputTag>( "electrons",edm::InputTag("gedGsfElectrons") );
183  desc.add<edm::InputTag>( "muons", edm::InputTag("muons") );
184  desc.add<std::string>("metSelection", "pt > 0");
185  desc.add<std::string>("jetSelection", "pt > 0");
186  desc.add<std::string>("jetId", "");
187  desc.add<std::string>("htjetSelection", "pt > 30");
188  desc.add<std::string>("eleSelection", "pt > 0");
189  desc.add<std::string>("muoSelection", "pt > 0");
190  desc.add<int>("njets", 0);
191  desc.add<int>("nelectrons", 0);
192  desc.add<int>("nmuons", 0);
193 
194  edm::ParameterSetDescription genericTriggerEventPSet;
195  genericTriggerEventPSet.add<bool>("andOr");
196  genericTriggerEventPSet.add<edm::InputTag>("dcsInputTag", edm::InputTag("scalersRawToDigi") );
197  genericTriggerEventPSet.add<std::vector<int> >("dcsPartitions",{});
198  genericTriggerEventPSet.add<bool>("andOrDcs", false);
199  genericTriggerEventPSet.add<bool>("errorReplyDcs", true);
200  genericTriggerEventPSet.add<std::string>("dbLabel","");
201  genericTriggerEventPSet.add<bool>("andOrHlt", true);
202  genericTriggerEventPSet.add<edm::InputTag>("hltInputTag", edm::InputTag("TriggerResults::HLT") );
203  genericTriggerEventPSet.add<std::vector<std::string> >("hltPaths",{});
204  genericTriggerEventPSet.add<std::string>("hltDBKey","");
205  genericTriggerEventPSet.add<bool>("errorReplyHlt",false);
206  genericTriggerEventPSet.add<unsigned int>("verbosityLevel",1);
207 
208  desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
209  desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
210 
211  desc.add<bool>("doMETHistos", true);
213  METDQM::fillMetDescription(histoPSet);
214  desc.add<bool>("doJetHistos", true);
215  JetDQM::fillJetDescription(histoPSet);
216  desc.add<bool>("doHTHistos", true);
217  HTDQM::fillHtDescription(histoPSet);
218 
219  desc.add<edm::ParameterSetDescription>("histoPSet",histoPSet);
220 
221  descriptions.add("objMonitoring", desc);
222 }
223 
224 // Define this as a plug-in
void fillHistograms(const double &met, const double &phi, const int &ls, const bool passCond)
Definition: METDQM.cc:39
virtual double pt() const final
transverse momentum
void initialise(const edm::ParameterSet &iConfig)
Definition: JetDQM.cc:8
void fillHistograms(const std::vector< reco::PFJet > &jets, const reco::PFMET &pfmet, const int &ls, const bool passCond)
Definition: JetDQM.cc:90
void bookHistograms(DQMStore::IBooker &)
Definition: METDQM.cc:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: ObjMonitor.cc:175
bool looseJetId(const double &abseta, const double &NHF, const double &NEMF, const double &CHF, const double &CEMF, const unsigned &NumNeutralParticles, const unsigned &CHM)
Definition: ObjMonitor.cc:135
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static void fillJetDescription(edm::ParameterSetDescription &histoPSet)
Definition: JetDQM.cc:180
void bookHistograms(DQMStore::IBooker &)
Definition: HTDQM.cc:16
static void fillMetDescription(edm::ParameterSetDescription &histoPSet)
Definition: METDQM.cc:62
edm::EDGetTokenT< reco::GsfElectronCollection > eleToken_
Definition: ObjMonitor.h:88
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
void initialise(const edm::ParameterSet &iConfig)
Definition: HTDQM.cc:7
bool tightJetId(const double &abseta, const double &NHF, const double &NEMF, const double &CHF, const double &CEMF, const unsigned &NumNeutralParticles, const unsigned &CHM)
Definition: ObjMonitor.cc:155
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
StringCutObjectSelector< reco::MET, true > metSelection_
Definition: ObjMonitor.h:104
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
virtual double phi() const final
momentum azimuthal angle
int iEvent
Definition: GenABIO.cc:230
std::string jetId_
Definition: ObjMonitor.h:106
std::unique_ptr< GenericTriggerEventFlag > num_genTriggerEventFlag_
Definition: ObjMonitor.h:101
edm::EDGetTokenT< reco::PFMETCollection > metToken_
Definition: ObjMonitor.h:86
void initialise(const edm::ParameterSet &iConfig)
Definition: METDQM.cc:7
void bookHistograms(DQMStore::IBooker &)
Definition: JetDQM.cc:23
vector< PseudoJet > jets
~ObjMonitor() override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
StringCutObjectSelector< reco::PFJet, true > jetSelection_
Definition: ObjMonitor.h:105
std::unique_ptr< GenericTriggerEventFlag > den_genTriggerEventFlag_
Definition: ObjMonitor.h:102
std::string folderName_
Definition: ObjMonitor.h:83
bool do_met_
Definition: ObjMonitor.h:93
static void fillHtDescription(edm::ParameterSetDescription &histoPSet)
Definition: HTDQM.cc:61
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool do_jet_
Definition: ObjMonitor.h:95
int nelectrons_
Definition: ObjMonitor.h:111
def ls(path, rec=False)
Definition: eostools.py:348
bool do_ht_
Definition: ObjMonitor.h:97
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
Definition: ObjMonitor.cc:70
met
===> hadronic RAZOR
JetDQM jetDQM_
Definition: ObjMonitor.h:96
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EventID id() const
Definition: EventBase.h:60
std::vector< PFJet > PFJetCollection
collection of PFJet objects
edm::EDGetTokenT< reco::PFJetCollection > jetToken_
Definition: ObjMonitor.h:87
fixed size matrix
HLT enums.
StringCutObjectSelector< reco::Muon, true > muoSelection_
Definition: ObjMonitor.h:109
StringCutObjectSelector< reco::PFJet, true > htjetSelection_
Definition: ObjMonitor.h:107
void fillHistograms(const std::vector< reco::PFJet > &htjets, const double &met, const int &ls, const bool passCond)
Definition: HTDQM.cc:35
edm::EDGetTokenT< reco::MuonCollection > muoToken_
Definition: ObjMonitor.h:89
METDQM metDQM_
Definition: ObjMonitor.h:94
ObjMonitor(const edm::ParameterSet &)
Definition: ObjMonitor.cc:18
StringCutObjectSelector< reco::GsfElectron, true > eleSelection_
Definition: ObjMonitor.h:108
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ObjMonitor.cc:52
Definition: Run.h:42
HTDQM htDQM_
Definition: ObjMonitor.h:98
Collection of PF MET.