CMS 3D CMS Logo

ObjMonitor.cc
Go to the documentation of this file.
2 
4 
6 
8 
14 
15 // -----------------------------
16 // constructors and destructor
17 // -----------------------------
18 
20  : folderName_(iConfig.getParameter<std::string>("FolderName")),
21  metToken_(consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("met"))),
22  jetToken_(mayConsume<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("jets"))),
23  eleToken_(mayConsume<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electrons"))),
24  muoToken_(mayConsume<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("muons"))),
25  phoToken_(mayConsume<reco::PhotonCollection>(iConfig.getParameter<edm::InputTag>("photons"))),
26  trkToken_(mayConsume<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))),
27  do_met_(iConfig.getParameter<bool>("doMETHistos")),
28  do_jet_(iConfig.getParameter<bool>("doJetHistos")),
29  do_ht_(iConfig.getParameter<bool>("doHTHistos")),
30  do_hmg_(iConfig.getParameter<bool>("doHMesonGammaHistos")),
31  num_genTriggerEventFlag_(std::make_unique<GenericTriggerEventFlag>(
32  iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"), consumesCollector(), *this)),
33  den_genTriggerEventFlag_(std::make_unique<GenericTriggerEventFlag>(
34  iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"), consumesCollector(), *this)),
35  metSelection_(iConfig.getParameter<std::string>("metSelection")),
36  jetSelection_(iConfig.getParameter<std::string>("jetSelection")),
37  jetId_(iConfig.getParameter<std::string>("jetId")),
38  htjetSelection_(iConfig.getParameter<std::string>("htjetSelection")),
39  eleSelection_(iConfig.getParameter<std::string>("eleSelection")),
40  muoSelection_(iConfig.getParameter<std::string>("muoSelection")),
41  phoSelection_(iConfig.getParameter<std::string>("phoSelection")),
42  trkSelection_(iConfig.getParameter<std::string>("trkSelection")),
43  njets_(iConfig.getParameter<int>("njets")),
44  nelectrons_(iConfig.getParameter<int>("nelectrons")),
45  nmuons_(iConfig.getParameter<int>("nmuons")),
46  nphotons_(iConfig.getParameter<int>("nphotons")),
47  nmesons_(iConfig.getParameter<int>("nmesons")) {
48  if (do_met_) {
49  metDQM_.initialise(iConfig);
50  }
51  if (do_jet_) {
52  jetDQM_.initialise(iConfig);
53  }
54  if (do_ht_) {
55  htDQM_.initialise(iConfig);
56  }
57  if (do_hmg_) {
58  hmgDQM_.initialise(iConfig);
59  }
60 }
61 
62 ObjMonitor::~ObjMonitor() = default;
63 
64 void ObjMonitor::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) {
65  std::string currentFolder = folderName_;
66  ibooker.setCurrentFolder(currentFolder);
67 
68  if (do_met_)
69  metDQM_.bookHistograms(ibooker);
70  if (do_jet_)
71  jetDQM_.bookHistograms(ibooker);
72  if (do_ht_)
73  htDQM_.bookHistograms(ibooker);
74  if (do_hmg_)
75  hmgDQM_.bookHistograms(ibooker);
76 
77  // Initialize the GenericTriggerEventFlag
79  num_genTriggerEventFlag_->initRun(iRun, iSetup);
81  den_genTriggerEventFlag_->initRun(iRun, iSetup);
82 }
83 
85  // Filter out events if Trigger Filtering is requested
86  if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
87  return;
88 
90  iEvent.getByToken(metToken_, metHandle);
91  reco::PFMET pfmet = metHandle->front();
92  if (!metSelection_(pfmet))
93  return;
94 
95  float met = pfmet.pt();
96  float phi = pfmet.phi();
97 
99  iEvent.getByToken(jetToken_, jetHandle);
100  std::vector<reco::PFJet> jets;
101  std::vector<reco::PFJet> htjets;
102  if (jetHandle->size() < njets_)
103  return;
104  for (auto const& j : *jetHandle) {
105  if (jetSelection_(j)) {
106  if (jetId_ == "loose" || jetId_ == "tight") {
107  double abseta = abs(j.eta());
108  double NHF = j.neutralHadronEnergyFraction();
109  double NEMF = j.neutralEmEnergyFraction();
110  double CHF = j.chargedHadronEnergyFraction();
111  double CEMF = j.chargedEmEnergyFraction();
112  unsigned NumNeutralParticles = j.neutralMultiplicity();
113  unsigned CHM = j.chargedMultiplicity();
114  bool passId = (jetId_ == "loose" && looseJetId(abseta, NHF, NEMF, CHF, CEMF, NumNeutralParticles, CHM)) ||
115  (jetId_ == "tight" && tightJetId(abseta, NHF, NEMF, CHF, CEMF, NumNeutralParticles, CHM));
116  if (passId)
117  jets.push_back(j);
118  } else
119  jets.push_back(j);
120  }
121  if (htjetSelection_(j))
122  htjets.push_back(j);
123  }
124  if (jets.size() < njets_)
125  return;
126 
128  iEvent.getByToken(eleToken_, eleHandle);
129  std::vector<reco::GsfElectron> electrons;
130  if (eleHandle->size() < nelectrons_)
131  return;
132  for (auto const& e : *eleHandle) {
133  if (eleSelection_(e))
134  electrons.push_back(e);
135  }
136  if (electrons.size() < nelectrons_)
137  return;
138 
140  iEvent.getByToken(muoToken_, muoHandle);
141  if (muoHandle->size() < nmuons_)
142  return;
143  std::vector<reco::Muon> muons;
144  for (auto const& m : *muoHandle) {
145  if (muoSelection_(m))
146  muons.push_back(m);
147  }
148  if (muons.size() < nmuons_)
149  return;
150 
152  iEvent.getByToken(phoToken_, phoHandle);
153  if (phoHandle->size() < nphotons_)
154  return;
155  std::vector<reco::Photon> photons;
156  for (auto const& m : *phoHandle) {
157  if (phoSelection_(m))
158  photons.push_back(m);
159  }
160  if (photons.size() < nphotons_)
161  return;
162 
163  std::vector<TLorentzVector> passedMesons;
164  if (do_hmg_) {
166  iEvent.getByToken(trkToken_, trkHandle);
167  // find isolated mesons (phi or rho)
168  TLorentzVector t1, t2;
169  float hadronMassHyp[2] = {0.1396, 0.4937}; // pi or K mass
170  float loMassLim[2] = {0.5, 0.9}; // rho or phi mass
171  float hiMassLim[2] = {1.0, 1.11}; // rho or phi mass
172 
173  for (size_t i = 0; i < trkHandle->size(); ++i) {
174  const reco::Track trk1 = trkHandle->at(i);
175  if (!trkSelection_(trk1))
176  continue;
177  for (size_t j = i + 1; j < trkHandle->size(); ++j) {
178  const reco::Track trk2 = trkHandle->at(j);
179  if (!trkSelection_(trk2))
180  continue;
181  if (trk1.charge() * trk2.charge() != -1)
182  continue;
183 
184  for (unsigned hyp = 0; hyp < 2; ++hyp) {
185  t1.SetPtEtaPhiM(trk1.pt(), trk1.eta(), trk1.phi(), hadronMassHyp[hyp]);
186  t2.SetPtEtaPhiM(trk2.pt(), trk2.eta(), trk2.phi(), hadronMassHyp[hyp]);
187  TLorentzVector mesCand = t1 + t2;
188 
189  // cuts
190  if (mesCand.M() < loMassLim[hyp] || mesCand.M() > hiMassLim[hyp])
191  continue; //mass
192  if (mesCand.Pt() < 35. || fabs(mesCand.Rapidity()) > 2.1)
193  continue; //pT eta
194 
195  // isolation
196  float absIso = 0.;
197  for (size_t k = 0; k < trkHandle->size(); ++k) {
198  if (k == i || k == j)
199  continue;
200  const reco::Track trkN = trkHandle->at(k);
201  if (trkN.charge() == 0 || trkN.pt() < 0.5 || (trkN.dz() > 0.1) ||
202  deltaR(trkN.eta(), trkN.phi(), mesCand.Eta(), mesCand.Phi()) > 0.5)
203  continue;
204  absIso += trkN.pt();
205  }
206  if (absIso / mesCand.Pt() > 0.2)
207  continue;
208  passedMesons.push_back(mesCand);
209  }
210  }
211  }
212  if (passedMesons.size() < nmesons_)
213  return;
214  }
215 
216  bool passNumCond = num_genTriggerEventFlag_->off() || num_genTriggerEventFlag_->accept(iEvent, iSetup);
217  int ls = iEvent.id().luminosityBlock();
218 
219  if (do_met_)
220  metDQM_.fillHistograms(met, phi, ls, passNumCond);
221  if (do_jet_)
222  jetDQM_.fillHistograms(jets, pfmet, ls, passNumCond);
223  if (do_ht_)
224  htDQM_.fillHistograms(htjets, met, ls, passNumCond);
225  if (do_hmg_)
226  hmgDQM_.fillHistograms(photons, passedMesons, ls, passNumCond);
227 }
228 
229 bool ObjMonitor::looseJetId(const double& abseta,
230  const double& NHF,
231  const double& NEMF,
232  const double& CHF,
233  const double& CEMF,
234  const unsigned& NumNeutralParticles,
235  const unsigned& CHM) {
236  if (abseta <= 2.7) {
237  unsigned NumConst = CHM + NumNeutralParticles;
238 
239  return ((NumConst > 1 && NHF < 0.99 && NEMF < 0.99) &&
240  ((abseta <= 2.4 && CHF > 0 && CHM > 0 && CEMF < 0.99) || abseta > 2.4));
241  } else if (abseta <= 3) {
242  return (NumNeutralParticles > 2 && NEMF > 0.01 && NHF < 0.98);
243  } else {
244  return NumNeutralParticles > 10 && NEMF < 0.90;
245  }
246 }
247 bool ObjMonitor::tightJetId(const double& abseta,
248  const double& NHF,
249  const double& NEMF,
250  const double& CHF,
251  const double& CEMF,
252  const unsigned& NumNeutralParticles,
253  const unsigned& CHM) {
254  if (abseta <= 2.7) {
255  unsigned NumConst = CHM + NumNeutralParticles;
256  return (NumConst > 1 && NHF < 0.90 && NEMF < 0.90) &&
257  ((abseta <= 2.4 && CHF > 0 && CHM > 0 && CEMF < 0.99) || abseta > 2.4);
258  } else if (abseta <= 3) {
259  return (NHF < 0.98 && NEMF > 0.01 && NumNeutralParticles > 2);
260  } else {
261  return (NEMF < 0.90 && NumNeutralParticles > 10);
262  }
263 }
264 
267  desc.add<std::string>("FolderName", "HLT/OBJ");
268 
269  desc.add<edm::InputTag>("met", edm::InputTag("pfMet"));
270  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
271  desc.add<edm::InputTag>("electrons", edm::InputTag("gedGsfElectrons"));
272  desc.add<edm::InputTag>("muons", edm::InputTag("muons"));
273  desc.add<edm::InputTag>("photons", edm::InputTag("gedPhotons"));
274  desc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
275  desc.add<std::string>("metSelection", "pt > 0");
276  desc.add<std::string>("jetSelection", "pt > 0");
277  desc.add<std::string>("jetId", "");
278  desc.add<std::string>("htjetSelection", "pt > 30");
279  desc.add<std::string>("eleSelection", "pt > 0");
280  desc.add<std::string>("muoSelection", "pt > 0");
281  desc.add<std::string>("phoSelection", "pt > 0");
282  desc.add<std::string>("trkSelection", "pt > 0");
283  desc.add<int>("njets", 0);
284  desc.add<int>("nelectrons", 0);
285  desc.add<int>("nmuons", 0);
286  desc.add<int>("nphotons", 0);
287  desc.add<int>("nmesons", 0);
288 
290  genericTriggerEventPSet.add<bool>("andOr");
291  genericTriggerEventPSet.add<edm::InputTag>("dcsInputTag", edm::InputTag("scalersRawToDigi"));
292  genericTriggerEventPSet.add<std::vector<int> >("dcsPartitions", {});
293  genericTriggerEventPSet.add<bool>("andOrDcs", false);
294  genericTriggerEventPSet.add<bool>("errorReplyDcs", true);
295  genericTriggerEventPSet.add<std::string>("dbLabel", "");
296  genericTriggerEventPSet.add<bool>("andOrHlt", true);
297  genericTriggerEventPSet.add<edm::InputTag>("hltInputTag", edm::InputTag("TriggerResults::HLT"));
298  genericTriggerEventPSet.add<std::vector<std::string> >("hltPaths", {});
299  genericTriggerEventPSet.add<std::string>("hltDBKey", "");
300  genericTriggerEventPSet.add<bool>("errorReplyHlt", false);
301  genericTriggerEventPSet.add<unsigned int>("verbosityLevel", 1);
302 
303  desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
304  desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
305 
306  desc.add<bool>("doMETHistos", true);
308  METDQM::fillMetDescription(histoPSet);
309  desc.add<bool>("doJetHistos", true);
310  JetDQM::fillJetDescription(histoPSet);
311  desc.add<bool>("doHTHistos", true);
312  HTDQM::fillHtDescription(histoPSet);
313  desc.add<bool>("doHMesonGammaHistos", true);
315 
316  desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
317 
318  descriptions.add("objMonitoring", desc);
319 }
320 
321 // Define this as a plug-in
void fillHistograms(const double &met, const double &phi, const int &ls, const bool passCond)
Definition: METDQM.cc:50
edm::EDGetTokenT< reco::TrackCollection > trkToken_
Definition: ObjMonitor.h:89
unsigned nelectrons_
Definition: ObjMonitor.h:115
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:131
void bookHistograms(DQMStore::IBooker &)
Definition: METDQM.cc:18
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: ObjMonitor.cc:265
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:229
StringCutObjectSelector< reco::GsfElectron, true > eleSelection_
Definition: ObjMonitor.h:109
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
StringCutObjectSelector< reco::MET, true > metSelection_
Definition: ObjMonitor.h:105
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
static void fillJetDescription(edm::ParameterSetDescription &histoPSet)
Definition: JetDQM.cc:221
void bookHistograms(DQMStore::IBooker &)
Definition: HTDQM.cc:18
static void fillMetDescription(edm::ParameterSetDescription &histoPSet)
Definition: METDQM.cc:68
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::GsfElectronCollection > eleToken_
Definition: ObjMonitor.h:86
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
double pt() const final
transverse momentum
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:247
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:39
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
static void fillHmgDescription(edm::ParameterSetDescription &histoPSet)
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
unsigned nmuons_
Definition: ObjMonitor.h:116
unsigned nphotons_
Definition: ObjMonitor.h:117
StringCutObjectSelector< reco::PFJet, true > jetSelection_
Definition: ObjMonitor.h:106
StringCutObjectSelector< reco::Muon, true > muoSelection_
Definition: ObjMonitor.h:110
void initialise(const edm::ParameterSet &iConfig)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::string jetId_
Definition: ObjMonitor.h:107
std::unique_ptr< GenericTriggerEventFlag > num_genTriggerEventFlag_
Definition: ObjMonitor.h:102
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:617
edm::EDGetTokenT< reco::PFMETCollection > metToken_
Definition: ObjMonitor.h:84
void initialise(const edm::ParameterSet &iConfig)
Definition: METDQM.cc:7
void bookHistograms(DQMStore::IBooker &)
Definition: JetDQM.cc:30
double pt() const
track transverse momentum
Definition: TrackBase.h:602
~ObjMonitor() override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< GenericTriggerEventFlag > den_genTriggerEventFlag_
Definition: ObjMonitor.h:103
std::string folderName_
Definition: ObjMonitor.h:81
bool do_met_
Definition: ObjMonitor.h:93
static void fillHtDescription(edm::ParameterSetDescription &histoPSet)
Definition: HTDQM.cc:69
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool do_jet_
Definition: ObjMonitor.h:95
edm::EDGetTokenT< reco::PhotonCollection > phoToken_
Definition: ObjMonitor.h:88
void bookHistograms(DQMStore::IBooker &)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:596
def ls(path, rec=False)
Definition: eostools.py:349
bool do_ht_
Definition: ObjMonitor.h:97
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
Definition: ObjMonitor.cc:84
JetDQM jetDQM_
Definition: ObjMonitor.h:96
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
unsigned nmesons_
Definition: ObjMonitor.h:118
HMesonGammaDQM hmgDQM_
Definition: ObjMonitor.h:100
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool do_hmg_
Definition: ObjMonitor.h:99
StringCutObjectSelector< reco::PFJet, true > htjetSelection_
Definition: ObjMonitor.h:108
edm::EventID id() const
Definition: EventBase.h:59
std::vector< PFJet > PFJetCollection
collection of PFJet objects
edm::EDGetTokenT< reco::PFJetCollection > jetToken_
Definition: ObjMonitor.h:85
fixed size matrix
HLT enums.
StringCutObjectSelector< reco::Photon, true > phoSelection_
Definition: ObjMonitor.h:111
void fillHistograms(const reco::PhotonCollection &photons, std::vector< TLorentzVector > mesons, const int &ls, const bool passCond)
void fillHistograms(const std::vector< reco::PFJet > &htjets, const double &met, const int &ls, const bool passCond)
Definition: HTDQM.cc:45
int charge() const
track electric charge
Definition: TrackBase.h:575
edm::EDGetTokenT< reco::MuonCollection > muoToken_
Definition: ObjMonitor.h:87
unsigned njets_
Definition: ObjMonitor.h:114
double phi() const final
momentum azimuthal angle
METDQM metDQM_
Definition: ObjMonitor.h:94
ObjMonitor(const edm::ParameterSet &)
Definition: ObjMonitor.cc:19
StringCutObjectSelector< reco::Track, true > trkSelection_
Definition: ObjMonitor.h:112
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ObjMonitor.cc:64
Definition: Run.h:45
HTDQM htDQM_
Definition: ObjMonitor.h:98
Collection of PF MET.