CMS 3D CMS Logo

SingleTopTChannelLeptonDQM.cc
Go to the documentation of this file.
11 #include <iostream>
12 #include <memory>
13 
18 
19 using namespace std;
21 
22  // maximal number of leading jets
23  // to be used for top mass estimate
24  static const unsigned int MAXJETS = 4;
25  // nominal mass of the W boson to
26  // be used for the top mass estimate
27  static const double WMASS = 80.4;
28 
29  MonitorEnsemble::MonitorEnsemble(const char* label,
30  const edm::ParameterSet& cfg,
31  const edm::VParameterSet& vcfg,
33  : label_(label),
34  elecSelect_(nullptr),
35  pvSelect_(nullptr),
36  muonIso_(nullptr),
37  muonSelect_(nullptr),
38  jetIDSelect_(nullptr),
39  jetlooseSelection_(nullptr),
40  jetSelection_(nullptr),
41  includeBTag_(false),
42  lowerEdge_(-1.),
43  upperEdge_(-1.),
44  logged_(0) {
45  // sources have to be given; this PSet is not optional
46  edm::ParameterSet sources = cfg.getParameter<edm::ParameterSet>("sources");
47  muons_ = iC.consumes<edm::View<reco::PFCandidate>>(sources.getParameter<edm::InputTag>("muons"));
48  elecs_ = iC.consumes<edm::View<reco::PFCandidate>>(sources.getParameter<edm::InputTag>("elecs"));
49  jets_ = iC.consumes<edm::View<reco::Jet>>(sources.getParameter<edm::InputTag>("jets"));
50  for (edm::InputTag const& tag : sources.getParameter<std::vector<edm::InputTag>>("mets"))
51  mets_.push_back(iC.consumes<edm::View<reco::MET>>(tag));
52  pvs_ = iC.consumes<edm::View<reco::Vertex>>(sources.getParameter<edm::InputTag>("pvs"));
53  // electronExtras are optional; they may be omitted or
54  // empty
55  if (cfg.existsAs<edm::ParameterSet>("elecExtras")) {
56  // rho for PF isolation with EA corrections
57  // eventrhoToken_ =
58  // iC.consumes<double>(edm::InputTag("fixedGridRhoFastjetAll"));
59 
60  edm::ParameterSet elecExtras = cfg.getParameter<edm::ParameterSet>("elecExtras");
61  // select is optional; in case it's not found no
62  // selection will be applied
63  if (elecExtras.existsAs<std::string>("select")) {
64  elecSelect_ = std::make_unique<StringCutObjectSelector<reco::PFCandidate>>(
65  elecExtras.getParameter<std::string>("select"));
66  }
67 
68  if (elecExtras.existsAs<std::string>("rho")) {
69  rhoTag = elecExtras.getParameter<edm::InputTag>("rho");
70  }
71  // electronId is optional; in case it's not found the
72  // InputTag will remain empty
73  if (elecExtras.existsAs<edm::ParameterSet>("electronId")) {
74  edm::ParameterSet elecId = elecExtras.getParameter<edm::ParameterSet>("electronId");
75  electronId_ = iC.consumes<edm::ValueMap<float>>(elecId.getParameter<edm::InputTag>("src"));
76  eidCutValue_ = elecId.getParameter<double>("cutValue");
77  }
78  }
79  // pvExtras are optional; they may be omitted or empty
80  if (cfg.existsAs<edm::ParameterSet>("pvExtras")) {
81  edm::ParameterSet pvExtras = cfg.getParameter<edm::ParameterSet>("pvExtras");
82  // select is optional; in case it's not found no
83  // selection will be applied
84  if (pvExtras.existsAs<std::string>("select")) {
85  pvSelect_ =
86  std::make_unique<StringCutObjectSelector<reco::Vertex>>(pvExtras.getParameter<std::string>("select"));
87  }
88  }
89  // muonExtras are optional; they may be omitted or empty
90  if (cfg.existsAs<edm::ParameterSet>("muonExtras")) {
91  edm::ParameterSet muonExtras = cfg.getParameter<edm::ParameterSet>("muonExtras");
92  // select is optional; in case it's not found no
93  // selection will be applied
94  if (muonExtras.existsAs<std::string>("select")) {
95  muonSelect_ = std::make_unique<StringCutObjectSelector<reco::PFCandidate>>(
96  muonExtras.getParameter<std::string>("select"));
97  }
98  // isolation is optional; in case it's not found no
99  // isolation will be applied
100  if (muonExtras.existsAs<std::string>("isolation")) {
101  muonIso_ = std::make_unique<StringCutObjectSelector<reco::PFCandidate>>(
102  muonExtras.getParameter<std::string>("isolation"));
103  }
104  }
105 
106  // jetExtras are optional; they may be omitted or
107  // empty
108  if (cfg.existsAs<edm::ParameterSet>("jetExtras")) {
109  edm::ParameterSet jetExtras = cfg.getParameter<edm::ParameterSet>("jetExtras");
110  // jetCorrector is optional; in case it's not found
111  // the InputTag will remain empty
112  if (jetExtras.existsAs<std::string>("jetCorrector")) {
113  jetCorrector_ =
114  iC.consumes<reco::JetCorrector>(edm::InputTag(jetExtras.getParameter<std::string>("jetCorrector")));
115  }
116  // read jetID information if it exists
117  if (jetExtras.existsAs<edm::ParameterSet>("jetID")) {
118  edm::ParameterSet jetID = jetExtras.getParameter<edm::ParameterSet>("jetID");
119  jetIDLabel_ = iC.consumes<reco::JetIDValueMap>(jetID.getParameter<edm::InputTag>("label"));
120  jetIDSelect_ =
121  std::make_unique<StringCutObjectSelector<reco::JetID>>(jetID.getParameter<std::string>("select"));
122  }
123  // select is optional; in case it's not found no
124  // selection will be applied (only implemented for
125  // CaloJets at the moment)
126  if (jetExtras.existsAs<std::string>("select")) {
127  jetSelect_ = jetExtras.getParameter<std::string>("select");
128  jetSelection_ = std::make_unique<StringCutObjectSelector<reco::PFJet>>(jetSelect_);
129  jetlooseSelection_ = std::make_unique<StringCutObjectSelector<reco::PFJet>>(
130  "chargedHadronEnergyFraction()>0 && chargedMultiplicity()>0 && chargedEmEnergyFraction()<0.99 && "
131  "neutralHadronEnergyFraction()<0.99 && neutralEmEnergyFraction()<0.99 && "
132  "(chargedMultiplicity()+neutralMultiplicity())>1");
133  }
134 
135  // jetBDiscriminators are optional; in case they are
136  // not found the InputTag will remain empty; they
137  // consist of pairs of edm::JetFlavorAssociation's &
138  // corresponding working points
139  includeBTag_ = jetExtras.existsAs<edm::ParameterSet>("jetBTaggers");
140  if (includeBTag_) {
141  edm::ParameterSet btagCSV =
142  jetExtras.getParameter<edm::ParameterSet>("jetBTaggers").getParameter<edm::ParameterSet>("cvsVertex");
143  btagCSV_ = iC.consumes<reco::JetTagCollection>(btagCSV.getParameter<edm::InputTag>("label"));
144  btagCSVWP_ = btagCSV.getParameter<double>("workingPoint");
145  }
146  }
147 
148  // triggerExtras are optional; they may be omitted or empty
149  if (cfg.existsAs<edm::ParameterSet>("triggerExtras")) {
150  edm::ParameterSet triggerExtras = cfg.getParameter<edm::ParameterSet>("triggerExtras");
151  triggerTable_ = iC.consumes<edm::TriggerResults>(triggerExtras.getParameter<edm::InputTag>("src"));
152  triggerPaths_ = triggerExtras.getParameter<std::vector<std::string>>("paths");
153  }
154 
155  // massExtras is optional; in case it's not found no mass
156  // window cuts are applied for the same flavor monitor
157  // histograms
158  if (cfg.existsAs<edm::ParameterSet>("massExtras")) {
159  edm::ParameterSet massExtras = cfg.getParameter<edm::ParameterSet>("massExtras");
160  lowerEdge_ = massExtras.getParameter<double>("lowerEdge");
161  upperEdge_ = massExtras.getParameter<double>("upperEdge");
162  }
163 
164  // setup the verbosity level for booking histograms;
165  // per default the verbosity level will be set to
166  // STANDARD. This will also be the chosen level in
167  // the case when the monitoring PSet is not found
169  if (cfg.existsAs<edm::ParameterSet>("monitoring")) {
170  edm::ParameterSet monitoring = cfg.getParameter<edm::ParameterSet>("monitoring");
171  if (monitoring.getParameter<std::string>("verbosity") == "DEBUG")
172  verbosity_ = DEBUG;
173  if (monitoring.getParameter<std::string>("verbosity") == "VERBOSE")
175  if (monitoring.getParameter<std::string>("verbosity") == "STANDARD")
177  }
178  // and don't forget to do the histogram booking
179  directory_ = cfg.getParameter<std::string>("directory");
180  }
181 
183  // set up the current directory path
184  std::string current(directory_);
185  current += label_;
186  ibooker.setCurrentFolder(current);
187 
188  // --- [STANDARD] --- //
189  // number of selected primary vertices
190  hists_["pvMult_"] = ibooker.book1D("PvMult", "N_{good pvs}", 50, 0., 50.);
191  // pt of the leading muon
192  hists_["muonPt_"] = ibooker.book1D("MuonPt", "pt(#mu TightId, TightIso)", 40, 0., 200.);
193  // muon multiplicity before std isolation
194  hists_["muonMult_"] = ibooker.book1D("MuonMult", "N_{loose}(#mu)", 10, 0., 10.);
195  // muon multiplicity tight Id and tight Iso
196  hists_["muonMultTight_"] = ibooker.book1D("MuonMultTight", "N_{TightIso,TightId}(#mu)", 10, 0., 10.);
197  // pt of the leading electron
198  hists_["elecPt_"] = ibooker.book1D("ElecPt", "pt(e TightId, TightIso)", 40, 0., 200.);
199  // multiplicity of jets with pt>30 (corrected to L1+L2+L3)
200  hists_["jetMult_"] = ibooker.book1D("JetMult", "N_{30}(jet)", 10, 0., 10.);
201  // multiplicity of loose Id jets with pt>30
202  hists_["jetLooseMult_"] = ibooker.book1D("JetLooseMult", "N_{30,loose}(jet)", 10, 0., 10.);
203 
204  // MET (pflow)
205  hists_["metPflow_"] = ibooker.book1D("METPflow", "MET_{Pflow}", 50, 0., 200.);
206  // W mass estimate
207  hists_["massW_"] = ibooker.book1D("MassW", "M(W)", 60, 0., 300.);
208  // Top mass estimate
209  hists_["massTop_"] = ibooker.book1D("MassTop", "M(Top)", 50, 0., 500.);
210  // W mass transverse estimate mu
211  hists_["MTWm_"] = ibooker.book1D("MTWm", "M_{T}^{W}(#mu)", 60, 0., 300.);
212  // Top mass transverse estimate mu
213  hists_["mMTT_"] = ibooker.book1D("mMTT", "M_{T}^{t}(#mu)", 50, 0., 500.);
214 
215  // W mass transverse estimate e
216  hists_["MTWe_"] = ibooker.book1D("MTWe", "M_{T}^{W}(e)", 60, 0., 300.);
217  // Top mass transverse estimate e
218  hists_["eMTT_"] = ibooker.book1D("eMTT", "M_{T}^{t}(e)", 50, 0., 500.);
219 
220  // set bin labels for trigger monitoring
222 
223  if (verbosity_ == STANDARD)
224  return;
225 
226  // --- [VERBOSE] --- //
227 
228  // eta of the leading muon
229  hists_["muonEta_"] = ibooker.book1D("MuonEta", "#eta(#mu TightId, TightIso)", 30, -3., 3.);
230  // relative isolation of the candidate muon (depending on the decay channel)
231  hists_["muonRelIso_"] = ibooker.book1D("MuonRelIso", "Iso_{Rel}(#mu TightId) (#Delta#beta Corrected)", 50, 0., 1.);
232  // phi of the leading muon
233  hists_["muonPhi_"] = ibooker.book1D("MuonPhi", "#phi(#mu TightId, TightIso)", 40, -4., 4.);
234  // eta of the leading electron
235  hists_["elecEta_"] = ibooker.book1D("ElecEta", "#eta(e tightId, TightIso)", 30, -3., 3.);
236  // std isolation variable of the leading electron
237  hists_["elecRelIso_"] = ibooker.book1D("ElecRelIso", "Iso_{Rel}(e TightId)", 50, 0., 1.);
238  // phi of the leading electron
239  hists_["elecPhi_"] = ibooker.book1D("ElecPhi", "#phi(e tightId, TightIso)", 40, -4., 4.);
240  // multiplicity of tight Id, tight Iso electorns
241  hists_["elecMultTight_"] = ibooker.book1D("ElecMultTight", "N_{TightIso,TightId}(e)", 10, 0., 10.);
242 
243  hists_["jet1Eta_"] = ibooker.book1D("Jet1Eta", "#eta_{30,loose}(jet1)", 60, -3., 3.);
244  // pt of the 1. leading jet (corrected to L2+L3)
245  hists_["jet1Pt_"] = ibooker.book1D("Jet1Pt", "pt_{30,loose}(jet1)", 60, 0., 300.);
246  // eta of the 2. leading jet (corrected to L2+L3)
247  hists_["jet2Eta_"] = ibooker.book1D("Jet2Eta", "#eta_{30,loose}(jet2)", 60, -3., 3.);
248  // pt of the 2. leading jet (corrected to L2+L3)
249  hists_["jet2Pt_"] = ibooker.book1D("Jet2Pt", "pt_{30,loose}(jet2)", 60, 0., 300.);
250 
251  // dz for muons (to suppress cosmis)
252  hists_["muonDelZ_"] = ibooker.book1D("MuonDelZ", "d_{z}(#mu)", 50, -25., 25.);
253  // dxy for muons (to suppress cosmics)
254  hists_["muonDelXY_"] = ibooker.book2D("MuonDelXY", "d_{xy}(#mu)", 50, -0.1, 0.1, 50, -0.1, 0.1);
255 
256  // set axes titles for dxy for muons
257  hists_["muonDelXY_"]->setAxisTitle("x [cm]", 1);
258  hists_["muonDelXY_"]->setAxisTitle("y [cm]", 2);
259 
260  if (verbosity_ == VERBOSE)
261  return;
262 
263  // --- [DEBUG] --- //
264  // charged hadron isolation component of the candidate muon (depending on the
265  // decay channel)
266  hists_["muonChHadIso_"] = ibooker.book1D("MuonChHadIsoComp", "ChHad_{IsoComponent}(#mu TightId)", 50, 0., 5.);
267  // neutral hadron isolation component of the candidate muon (depending on the
268  // decay channel)
269  hists_["muonNeHadIso_"] = ibooker.book1D("MuonNeHadIsoComp", "NeHad_{IsoComponent}(#mu TightId)", 50, 0., 5.);
270  // photon isolation component of the candidate muon (depending on the decay
271  // channel)
272  hists_["muonPhIso_"] = ibooker.book1D("MuonPhIsoComp", "Photon_{IsoComponent}(#mu TightId)", 50, 0., 5.);
273  // charged hadron isolation component of the candidate electron (depending on
274  // the decay channel)
275  hists_["elecChHadIso_"] = ibooker.book1D("ElectronChHadIsoComp", "ChHad_{IsoComponent}(e tightId)", 50, 0., 5.);
276  // neutral hadron isolation component of the candidate electron (depending on
277  // the decay channel)
278  hists_["elecNeHadIso_"] = ibooker.book1D("ElectronNeHadIsoComp", "NeHad_{IsoComponent}(e tightId)", 50, 0., 5.);
279  // photon isolation component of the candidate electron (depending on the
280  // decay channel)
281  hists_["elecPhIso_"] = ibooker.book1D("ElectronPhIsoComp", "Photon_{IsoComponent}(e tightId)", 50, 0., 5.);
282 
283  // multiplicity for combined secondary vertex
284  hists_["jetMultBCSVM_"] = ibooker.book1D("JetMultBCSVM", "N_{30}(CSVM)", 10, 0., 10.);
285  // btag discriminator for combined secondary vertex
286  hists_["jetBCSV_"] = ibooker.book1D("JetDiscCSV", "BJet Disc_{CSV}(JET)", 100, -1., 2.);
287  // pt of the 1. leading jet (uncorrected)
288  // hists_["jet1PtRaw_"] = ibooker.book1D("Jet1PtRaw", "pt_{Raw}(jet1)", 60, 0., 300.);
289  // pt of the 2. leading jet (uncorrected)
290  // hists_["jet2PtRaw_"] = ibooker.book1D("Jet2PtRaw", "pt_{Raw}(jet2)", 60, 0., 300.);
291  // pt of the 3. leading jet (uncorrected)
292  // hists_["jet3PtRaw_"] = ibooker.book1D("Jet3PtRaw", "pt_{Raw}(jet3)", 60, 0., 300.);
293  // pt of the 4. leading jet (uncorrected)
294  // hists_["jet4PtRaw_"] = ibooker.book1D("Jet4PtRaw", "pt_{Raw}(jet4)", 60, 0., 300.);
295  // selected events
296  hists_["eventLogger_"] = ibooker.book2D("EventLogger", "Logged Events", 9, 0., 9., 10, 0., 10.);
297 
298  // set axes titles for selected events
299  hists_["eventLogger_"]->getTH1()->SetOption("TEXT");
300  hists_["eventLogger_"]->setBinLabel(1, "Run", 1);
301  hists_["eventLogger_"]->setBinLabel(2, "Block", 1);
302  hists_["eventLogger_"]->setBinLabel(3, "Event", 1);
303  hists_["eventLogger_"]->setBinLabel(4, "pt_{30,loose}(jet1)", 1);
304  hists_["eventLogger_"]->setBinLabel(5, "pt_{30,loose}(jet2)", 1);
305  hists_["eventLogger_"]->setBinLabel(6, "pt_{30,loose}(jet3)", 1);
306  hists_["eventLogger_"]->setBinLabel(7, "pt_{30,loose}(jet4)", 1);
307  hists_["eventLogger_"]->setBinLabel(8, "M_{W}", 1);
308  hists_["eventLogger_"]->setBinLabel(9, "M_{Top}", 1);
309  hists_["eventLogger_"]->setAxisTitle("logged evts", 2);
310  return;
311  }
312 
314  // fetch trigger event if configured such
316 
317  /*
318  ------------------------------------------------------------
319 
320  Primary Vertex Monitoring
321 
322  ------------------------------------------------------------
323  */
324 
325  // fill monitoring plots for primary verices
327 
328  if (!event.getByToken(pvs_, pvs))
329  return;
330  const reco::Vertex& Pvertex = pvs->front();
331  unsigned int pvMult = 0;
332  for (edm::View<reco::Vertex>::const_iterator pv = pvs->begin(); pv != pvs->end(); ++pv) {
333  if (!pvSelect_ || (*pvSelect_)(*pv))
334  pvMult++;
335  }
336  fill("pvMult_", pvMult);
337  /*
338  ------------------------------------------------------------
339 
340  Electron Monitoring
341 
342  ------------------------------------------------------------
343  */
344 
345  // fill monitoring plots for electrons
348  edm::Handle<double> _rhoHandle;
349  event.getByLabel(rhoTag, _rhoHandle);
350 
351  if (!event.getByToken(elecs_, elecs))
352  return;
353 
354  // check availability of electron id
356  if (!electronId_.isUninitialized()) {
357  if (!event.getByToken(electronId_, electronId)) {
358  return;
359  }
360  }
361  // loop electron collection
362  unsigned int eMult = 0, eMultIso = 0;
363  std::vector<const reco::PFCandidate*> isoElecs;
364 
365  for (edm::View<reco::PFCandidate>::const_iterator elec = elecs->begin(); elec != elecs->end(); ++elec) {
366  if (elec->gsfElectronRef().isNull()) {
367  continue;
368  }
369  reco::GsfElectronRef gsf_el = elec->gsfElectronRef();
370  // restrict to electrons with good electronId
372  ? true
373  : ((double)(*electronId)[gsf_el] >=
374  eidCutValue_)) { //This Electron Id is not currently used, but we can keep this for future needs
375  if (!elecSelect_ || (*elecSelect_)(*elec)) {
376  double el_ChHadIso = gsf_el->pfIsolationVariables().sumChargedHadronPt;
377  double el_NeHadIso = gsf_el->pfIsolationVariables().sumNeutralHadronEt;
378  double el_PhIso = gsf_el->pfIsolationVariables().sumPhotonEt;
379  double absEta = std::abs(gsf_el->superCluster()->eta());
380 
381  //Effective Area computation
382  double eA = 0;
383  if (absEta < 1.000)
384  eA = 0.1703;
385  else if (absEta < 1.479)
386  eA = 0.1715;
387  else if (absEta < 2.000)
388  eA = 0.1213;
389  else if (absEta < 2.200)
390  eA = 0.1230;
391  else if (absEta < 2.300)
392  eA = 0.1635;
393  else if (absEta < 2.400)
394  eA = 0.1937;
395  else if (absEta < 5.000)
396  eA = 0.2393;
397 
398  double rho = _rhoHandle.isValid() ? (float)(*_rhoHandle) : 0;
399  double el_pfRelIso = (el_ChHadIso + max(0., el_NeHadIso + el_PhIso - rho * eA)) / gsf_el->pt();
400 
401  //Only TightId
402  if (eMult == 0) { // Restricted to the leading tight electron
403  fill("elecRelIso_", el_pfRelIso);
404  fill("elecChHadIso_", el_ChHadIso);
405  fill("elecNeHadIso_", el_NeHadIso);
406  fill("elecPhIso_", el_PhIso);
407  }
408  ++eMult;
409 
410  if (!((el_pfRelIso < 0.0588 && absEta < 1.479) || (el_pfRelIso < 0.0571 && absEta > 1.479)))
411  continue; // PF Isolation with Effective Area Corrections according to https://twiki.cern.ch/twiki/bin/viewauth/CMS/CutBasedElectronIdentificationRun2
412 
413  // TightId and TightIso
414  if (eMultIso == 0) { //Only leading
415  e = *gsf_el;
416  fill("elecPt_", gsf_el->pt());
417  fill("elecEta_", gsf_el->eta());
418  fill("elecPhi_", gsf_el->phi());
419  }
420  ++eMultIso;
421  }
422  }
423  }
424  //fill("elecMult_", eMult);
425  fill("elecMultTight_", eMultIso);
426 
427  /*
428  ------------------------------------------------------------
429 
430  Muon Monitoring
431 
432  ------------------------------------------------------------
433  */
434 
435  // fill monitoring plots for muons
436  unsigned int mMult = 0, mTight = 0, mTightId = 0;
437 
441  reco::Muon mu;
442 
443  if (!event.getByToken(muons_, muons))
444  return;
445 
446  for (edm::View<reco::PFCandidate>::const_iterator muonit = muons->begin(); muonit != muons->end(); ++muonit) {
447  if (muonit->muonRef().isNull())
448  continue;
449  reco::MuonRef muon = muonit->muonRef();
450 
451  // restrict to globalMuons
452  if (muon->isGlobalMuon()) {
453  fill("muonDelZ_", muon->innerTrack()->vz()); // CB using inner track!
454  fill("muonDelXY_", muon->innerTrack()->vx(), muon->innerTrack()->vy());
455 
456  // apply preselection
457  if ((!muonSelect_ || (*muonSelect_)(*muonit))) {
458  mMult++;
459  double chHadPt = muon->pfIsolationR04().sumChargedHadronPt;
460  double neHadEt = muon->pfIsolationR04().sumNeutralHadronEt;
461  double phoEt = muon->pfIsolationR04().sumPhotonEt;
462  double pfRelIso = (chHadPt + max(0., neHadEt + phoEt - 0.5 * muon->pfIsolationR04().sumPUPt)) /
463  muon->pt(); // CB dBeta corrected iso!
464 
465  if (!(muon->isGlobalMuon() && muon->isPFMuon() && muon->globalTrack()->normalizedChi2() < 10. &&
466  muon->globalTrack()->hitPattern().numberOfValidMuonHits() > 0 && muon->numberOfMatchedStations() > 1 &&
467  muon->innerTrack()->hitPattern().numberOfValidPixelHits() > 0 &&
468  muon->innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
469  fabs(muon->muonBestTrack()->dxy(Pvertex.position())) < 0.2 &&
470  fabs(muon->muonBestTrack()->dz(Pvertex.position())) < 0.5))
471  continue; //Only tight muons
472 
473  if (mTightId == 0) {
474  fill("muonRelIso_", pfRelIso);
475  fill("muonChHadIso_", chHadPt);
476  fill("muonNeHadIso_", neHadEt);
477  fill("muonPhIso_", phoEt);
478  }
479  mTightId++;
480 
481  if (!(pfRelIso < 0.15))
482  continue; //Tight Iso
483 
484  if (mTight == 0) { //Leading tightId tightIso muon
485  mu = *(muon);
486  fill("muonPt_", muon->pt());
487  fill("muonEta_", muon->eta());
488  fill("muonPhi_", muon->phi());
489  }
490  mTight++;
491  }
492  }
493  }
494  fill("muonMult_", mMult);
495  fill("muonMultTight_", mTight);
496 
497  /*
498  ------------------------------------------------------------
499 
500  Jet Monitoring
501 
502  ------------------------------------------------------------
503  */
504 
505  // check availability of the btaggers
506  edm::Handle<reco::JetTagCollection> btagEff, btagPur, btagVtx, btagCSV;
507  if (includeBTag_) {
508  if (!event.getByToken(btagCSV_, btagCSV))
509  return;
510  }
511 
512  // load jet
513  // corrector if configured such
514  const reco::JetCorrector* corrector = nullptr;
516  // check whether a jet corrector is in the event or not
517  edm::Handle<reco::JetCorrector> correctorHandle = event.getHandle(jetCorrector_);
518  if (correctorHandle.isValid()) {
519  corrector = correctorHandle.product();
520  } else {
521  edm::LogVerbatim("SingleTopTChannelLeptonDQM")
522  << "\n"
523  << "-----------------------------------------------------------------"
524  "-------------------- \n"
525  << " No JetCorrector available from Event:\n"
526  << " - Jets will not be corrected.\n"
527  << "-----------------------------------------------------------------"
528  "-------------------- \n";
529  }
530  }
531  // loop jet collection
532  std::vector<reco::Jet> correctedJets;
533  std::vector<double> JetTagValues;
534  reco::Jet TaggedJetCand;
535  unsigned int mult = 0, multLoose = 0, multCSV = 0;
536  vector<double> bJetDiscVal;
537 
539  if (!event.getByToken(jets_, jets)) {
540  return;
541  }
542 
543  for (edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
544  bool isLoose = false;
545  // check jetID for calo jets, keep functionality if we ever want to go back to those
546  // unsigned int idx = jet - jets->begin();
547  // if (dynamic_cast<const reco::CaloJet*>(&*jet)) {
548  // if (jetIDSelect_ &&
549  // dynamic_cast<const reco::CaloJet*>(jets->refAt(idx).get())) {
550  // if (!(*jetIDSelect_)((*jetID)[jets->refAt(idx)])) continue;
551  // }
552  // }
553 
554  // check additional jet selection for pf jets
555  if (dynamic_cast<const reco::PFJet*>(&*jet)) {
556  reco::PFJet sel = dynamic_cast<const reco::PFJet&>(*jet);
557  if ((*jetlooseSelection_)(sel))
558  isLoose = true;
559  sel.scaleEnergy(corrector ? corrector->correction(*jet) : 1.);
560  if (!(*jetSelection_)(sel))
561  continue;
562  }
563  // prepare jet to fill monitor histograms
564  reco::Jet monitorJet = *jet;
565 
566  ++mult; // determine jet (no Id) multiplicity
567  monitorJet.scaleEnergy(corrector ? corrector->correction(*jet) : 1.);
568 
569  if (isLoose) { //Loose Id
570  unsigned int idx = jet - jets->begin();
571  correctedJets.push_back(monitorJet);
572  if (includeBTag_) {
573  // fill b-discriminators
574  edm::RefToBase<reco::Jet> jetRef = jets->refAt(idx);
575  fill("jetBCSV_", (*btagCSV)[jetRef]);
576  if ((*btagCSV)[jetRef] > btagCSVWP_) {
577  if (multCSV == 0) {
578  TaggedJetCand = monitorJet;
579  bJetDiscVal.push_back((*btagCSV)[jetRef]);
580  } else if (multCSV == 1) {
581  bJetDiscVal.push_back((*btagCSV)[jetRef]);
582  if (bJetDiscVal[1] > bJetDiscVal[0])
583  TaggedJetCand = monitorJet;
584  }
585  ++multCSV;
586  }
587  JetTagValues.push_back((*btagCSV)[jetRef]);
588  }
589 
590  // fill pt/eta for the leading four jets
591  if (multLoose == 0) {
592  fill("jet1Pt_", monitorJet.pt());
593  fill("jet1Eta_", monitorJet.eta());
594  };
595  if (multLoose == 1) {
596  fill("jet2Pt_", monitorJet.pt());
597  fill("jet2Eta_", monitorJet.eta());
598  }
599  multLoose++;
600  }
601  }
602  fill("jetMult_", mult);
603  fill("jetMultLoose_", multLoose);
604  fill("jetMultBCSVM_", multCSV);
605 
606  /*
607  ------------------------------------------------------------
608 
609  MET Monitoring
610 
611  ------------------------------------------------------------
612  */
613 
614  // fill monitoring histograms for met
615  reco::MET mET;
616  for (std::vector<edm::EDGetTokenT<edm::View<reco::MET>>>::const_iterator met_ = mets_.begin(); met_ != mets_.end();
617  ++met_) {
619  if (!event.getByToken(*met_, met))
620  continue;
621  if (met->begin() != met->end()) {
622  unsigned int idx = met_ - mets_.begin();
623  if (idx == 0) {
624  fill("metPflow_", met->begin()->et());
625  mET = *(met->begin());
626  }
627  }
628  }
629 
630  /*
631  ------------------------------------------------------------
632 
633  Event Monitoring
634 
635  ------------------------------------------------------------
636  */
637 
638  // fill W boson and top mass estimates
639  Calculate eventKinematics(MAXJETS, WMASS);
640  double wMass = eventKinematics.massWBoson(correctedJets);
641  double topMass = eventKinematics.massTopQuark(correctedJets);
642  if (wMass >= 0 && topMass >= 0) {
643  fill("massW_", wMass);
644  fill("massTop_", topMass);
645  }
646  // fill plots for trigger monitoring
647  if ((lowerEdge_ == -1. && upperEdge_ == -1.) || (lowerEdge_ < wMass && wMass < upperEdge_)) {
649  fill(event, *triggerTable, "trigger", triggerPaths_);
650  if (logged_ <= hists_.find("eventLogger_")->second->getNbinsY()) {
651  // log runnumber, lumi block, event number & some
652  // more pysics infomation for interesting events
653  fill("eventLogger_", 0.5, logged_ + 0.5, event.eventAuxiliary().run());
654  fill("eventLogger_", 1.5, logged_ + 0.5, event.eventAuxiliary().luminosityBlock());
655  fill("eventLogger_", 2.5, logged_ + 0.5, event.eventAuxiliary().event());
656  if (!correctedJets.empty())
657  fill("eventLogger_", 3.5, logged_ + 0.5, correctedJets[0].pt());
658  if (correctedJets.size() > 1)
659  fill("eventLogger_", 4.5, logged_ + 0.5, correctedJets[1].pt());
660  if (correctedJets.size() > 2)
661  fill("eventLogger_", 5.5, logged_ + 0.5, correctedJets[2].pt());
662  if (correctedJets.size() > 3)
663  fill("eventLogger_", 6.5, logged_ + 0.5, correctedJets[3].pt());
664  fill("eventLogger_", 7.5, logged_ + 0.5, wMass);
665  fill("eventLogger_", 8.5, logged_ + 0.5, topMass);
666  ++logged_;
667  }
668  }
669  if (multCSV != 0 && mTight == 1) {
670  double mtW = eventKinematics.tmassWBoson(&mu, mET, TaggedJetCand);
671  fill("MTWm_", mtW);
672  double MTT = eventKinematics.tmassTopQuark(&mu, mET, TaggedJetCand);
673  fill("mMTT_", MTT);
674  }
675 
676  if (multCSV != 0 && eMultIso == 1) {
677  double mtW = eventKinematics.tmassWBoson(&e, mET, TaggedJetCand);
678  fill("MTWe_", mtW);
679  double MTT = eventKinematics.tmassTopQuark(&e, mET, TaggedJetCand);
680  fill("eMTT_", MTT);
681  }
682  }
683 } // namespace SingleTopTChannelLepton
684 
686  : vertexSelect_(nullptr),
687  beamspot_(""),
688  beamspotSelect_(nullptr),
689  MuonStep(nullptr),
690  PFMuonStep(nullptr),
691  ElectronStep(nullptr),
692  PFElectronStep(nullptr),
693  PvStep(nullptr),
694  METStep(nullptr) {
695  JetSteps.clear();
696  CaloJetSteps.clear();
697  PFJetSteps.clear();
698  // configure preselection
699  edm::ParameterSet presel = cfg.getParameter<edm::ParameterSet>("preselection");
700  if (presel.existsAs<edm::ParameterSet>("trigger")) {
702  triggerTable__ = consumes<edm::TriggerResults>(trigger.getParameter<edm::InputTag>("src"));
703  triggerPaths_ = trigger.getParameter<std::vector<std::string>>("select");
704  }
705  if (presel.existsAs<edm::ParameterSet>("beamspot")) {
707  beamspot_ = beamspot.getParameter<edm::InputTag>("src");
708  beamspot__ = consumes<reco::BeamSpot>(beamspot.getParameter<edm::InputTag>("src"));
710  std::make_unique<StringCutObjectSelector<reco::BeamSpot>>(beamspot.getParameter<std::string>("select"));
711  }
712  // configure the selection
713  std::vector<edm::ParameterSet> sel = cfg.getParameter<std::vector<edm::ParameterSet>>("selection");
714 
715  for (unsigned int i = 0; i < sel.size(); ++i) {
716  selectionOrder_.push_back(sel.at(i).getParameter<std::string>("label"));
718  std::make_pair(sel.at(i),
719  std::make_unique<SingleTopTChannelLepton::MonitorEnsemble>(
720  selectionStep(selectionOrder_.back()).c_str(),
721  cfg.getParameter<edm::ParameterSet>("setup"),
722  cfg.getParameter<std::vector<edm::ParameterSet>>("selection"),
723  consumesCollector()));
724  }
725  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin(); selIt != selectionOrder_.end();
726  ++selIt) {
727  std::string key = selectionStep(*selIt), type = objectType(*selIt);
728  if (selection_.find(key) != selection_.end()) {
729  using std::unique_ptr;
730 
731  if (type == "muons") {
732  MuonStep = std::make_unique<SelectionStep<reco::Muon>>(selection_[key].first, consumesCollector());
733  }
734  if (type == "muons/pf") {
735  PFMuonStep = std::make_unique<SelectionStep<reco::PFCandidate>>(selection_[key].first, consumesCollector());
736  }
737  if (type == "elecs") {
738  ElectronStep = std::make_unique<SelectionStep<reco::GsfElectron>>(selection_[key].first, consumesCollector());
739  }
740  if (type == "elecs/pf") {
741  PFElectronStep = std::make_unique<SelectionStep<reco::PFCandidate>>(selection_[key].first, consumesCollector());
742  }
743  if (type == "pvs") {
744  PvStep = std::make_unique<SelectionStep<reco::Vertex>>(selection_[key].first, consumesCollector());
745  }
746  if (type == "jets") {
747  JetSteps.push_back(std::make_unique<SelectionStep<reco::Jet>>(selection_[key].first, consumesCollector()));
748  }
749  if (type == "jets/pf") {
751  }
752  if (type == "jets/calo") {
753  CaloJetSteps.push_back(
755  }
756  if (type == "met") {
757  METStep = std::make_unique<SelectionStep<reco::MET>>(selection_[key].first, consumesCollector());
758  }
759  }
760  }
761 }
763  for (auto selIt = selection_.begin(); selIt != selection_.end(); ++selIt) {
764  selIt->second.second->book(ibooker);
765  }
766 }
770  if (!event.getByToken(triggerTable__, triggerTable))
771  return;
772  if (!accept(event, *triggerTable, triggerPaths_))
773  return;
774  }
775  if (!beamspot__.isUninitialized()) {
777  if (!event.getByToken(beamspot__, beamspot))
778  return;
779  if (!(*beamspotSelect_)(*beamspot))
780  return;
781  }
782  // apply selection steps
783  unsigned int nJetSteps = -1;
784  unsigned int nPFJetSteps = -1;
785  unsigned int nCaloJetSteps = -1;
786  for (std::vector<std::string>::const_iterator selIt = selectionOrder_.begin(); selIt != selectionOrder_.end();
787  ++selIt) {
788  std::string key = selectionStep(*selIt), type = objectType(*selIt);
789  if (selection_.find(key) != selection_.end()) {
790  if (type == "empty") {
791  selection_[key].second->fill(event, setup);
792  }
793  if (type == "presel") {
794  selection_[key].second->fill(event, setup);
795  }
796  if (type == "elecs" && ElectronStep != nullptr) {
797  if (ElectronStep->select(event)) {
798  selection_[key].second->fill(event, setup);
799  } else
800  break;
801  }
802  if (type == "elecs/pf" && PFElectronStep != nullptr) {
803  if (PFElectronStep->select(event, "electron")) {
804  selection_[key].second->fill(event, setup);
805 
806  } else
807  break;
808  }
809  if (type == "muons" && MuonStep != nullptr) {
810  if (MuonStep->select(event)) {
811  selection_[key].second->fill(event, setup);
812  } else
813  break;
814  }
815  if (type == "muons/pf" && PFMuonStep != nullptr) {
816  if (PFMuonStep->select(event, "muon")) {
817  selection_[key].second->fill(event, setup);
818  } else
819  break;
820  }
821  if (type == "jets") {
822  nJetSteps++;
823  if (JetSteps[nJetSteps]) {
824  if (JetSteps[nJetSteps]->select(event, setup)) {
825  selection_[key].second->fill(event, setup);
826  } else
827  break;
828  }
829  }
830  if (type == "jets/pf") {
831  nPFJetSteps++;
832  if (PFJetSteps[nPFJetSteps]) {
833  if (PFJetSteps[nPFJetSteps]->select(event, setup)) {
834  selection_[key].second->fill(event, setup);
835  } else
836  break;
837  }
838  }
839  if (type == "jets/calo") {
840  nCaloJetSteps++;
841  if (CaloJetSteps[nCaloJetSteps]) {
842  if (CaloJetSteps[nCaloJetSteps]->select(event, setup)) {
843  selection_[key].second->fill(event, setup);
844  } else
845  break;
846  }
847  }
848  if (type == "met" && METStep != nullptr) {
849  if (METStep->select(event)) {
850  selection_[key].second->fill(event, setup);
851  } else
852  break;
853  }
854  }
855  }
856 }
std::vector< std::string > selectionOrder_
std::map< std::string, MonitorElement * > hists_
Log< level::Info, true > LogVerbatim
std::string objectType(const std::string &label)
edm::EDGetTokenT< reco::JetIDValueMap > jetIDLabel_
jetID as an extra selection type
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
virtual void scaleEnergy(double fScale)
scale energy of the jet
std::vector< std::string > triggerPaths_
trigger paths
Level verbosity_
verbosity level for booking
double pt() const final
transverse momentum
std::unique_ptr< StringCutObjectSelector< reco::PFJet > > jetlooseSelection_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::unique_ptr< SelectionStep< reco::Vertex > > PvStep
const Point & position() const
position
Definition: Vertex.h:128
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:35
Base class for all types of Jets.
Definition: Jet.h:20
int logged_
number of logged interesting events
std::vector< std::unique_ptr< SelectionStep< reco::CaloJet > > > CaloJetSteps
std::vector< std::unique_ptr< SelectionStep< reco::Jet > > > JetSteps
T const * product() const
Definition: Handle.h:70
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:104
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate t-quark mass estimate
edm::EDGetTokenT< edm::TriggerResults > triggerTable__
trigger table
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:172
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:212
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
edm::EDGetTokenT< edm::View< reco::PFCandidate > > muons_
Jets made from PFObjects.
Definition: PFJet.h:20
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
SingleTopTChannelLeptonDQM(const edm::ParameterSet &cfg)
default constructor
electronId
when omitted electron plots will be filled w/o cut on electronId
char const * label
std::unique_ptr< StringCutObjectSelector< reco::PFCandidate > > muonSelect_
extra selection on muons
Helper class for the calculation of a top and a W boson mass estime.
std::unique_ptr< StringCutObjectSelector< reco::PFCandidate > > muonIso_
extra isolation criterion on muon
double lowerEdge_
mass window upper and lower edge
void triggerBinLabels(std::string channel, const std::vector< std::string > labels)
set configurable labels for trigger monitoring histograms
std::vector< edm::ParameterSet > sel
Definition: MET.h:41
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< StringCutObjectSelector< reco::Vertex > > pvSelect_
key
prepare the HTCondor submission files and eventually submit them
double tmassTopQuark(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate top quark transverse mass estimate
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::unique_ptr< StringCutObjectSelector< reco::JetID > > jetIDSelect_
extra jetID selection on calo jets
std::unique_ptr< SelectionStep< reco::PFCandidate > > PFElectronStep
edm::EDGetTokenT< edm::TriggerResults > triggerTable_
trigger table
edm::EDGetTokenT< reco::JetCorrector > jetCorrector_
jetCorrector
edm::EDGetTokenT< edm::ValueMap< float > > electronId_
electronId label
edm::EDGetTokenT< reco::JetTagCollection > btagCSV_
std::string selectionStep(const std::string &label)
double tmassWBoson(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate W boson transverse mass estimate
std::unique_ptr< SelectionStep< reco::GsfElectron > > ElectronStep
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
edm::EDGetTokenT< edm::View< reco::PFCandidate > > elecs_
std::unique_ptr< SelectionStep< reco::MET > > METStep
static const unsigned int MAXJETS
bool isValid() const
Definition: HandleBase.h:70
std::vector< edm::EDGetTokenT< edm::View< reco::MET > > > mets_
considers a vector of METs
Templated helper class to allow a selection on a certain object collection.
void fill(const edm::Event &event, const edm::EventSetup &setup)
fill monitor histograms with electronId and jetCorrections
edm::EDGetTokenT< edm::View< reco::Jet > > jets_
input sources for monitoring
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
void book(DQMStore::IBooker &ibooker)
book histograms in subdirectory directory
edm::EDGetTokenT< reco::BeamSpot > beamspot__
std::unique_ptr< StringCutObjectSelector< reco::BeamSpot > > beamspotSelect_
string cut selector
std::unique_ptr< StringCutObjectSelector< reco::PFJet > > jetSelection_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
do this during the event loop
std::vector< std::unique_ptr< SelectionStep< reco::PFJet > > > PFJetSteps
std::map< std::string, std::pair< edm::ParameterSet, std::unique_ptr< SingleTopTChannelLepton::MonitorEnsemble > > > selection_
std::unique_ptr< SelectionStep< reco::Muon > > MuonStep
std::unique_ptr< StringCutObjectSelector< reco::PFCandidate > > elecSelect_
extra selection on electrons
Definition: event.py:1
Definition: Run.h:45
std::unique_ptr< SelectionStep< reco::PFCandidate > > PFMuonStep
double eta() const final
momentum pseudorapidity