CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EwkDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Valentina Gori, University of Firenze
5  */
6 
8 
9 #include <vector>
10 #include <string>
11 #include <cmath>
12 
17 
21 
23 
26 
27 // Physics Objects
36 
37 #include "TLorentzVector.h"
38 
39 using namespace std;
40 using namespace edm;
41 using namespace reco;
42 
43 
44 
45 // EwkDQM::EwkDQM(const ParameterSet& parameters) {
47  eJetMin_ = parameters.getUntrackedParameter<double>("EJetMin", 999999.);
48 
49  // riguardare questa sintassi
50  // Get parameters from configuration file
51  thePFJetCollectionLabel_ =
52  parameters.getParameter<InputTag>("PFJetCollection");
53  theCaloMETCollectionLabel_ =
54  parameters.getParameter<InputTag>("caloMETCollection");
55  theTriggerResultsCollection_ =
56  parameters.getParameter<InputTag>("triggerResultsCollection");
57 
58  theElecTriggerPathToPass_ =
59  parameters.getParameter<std::vector<string> >("elecTriggerPathToPass");
60  theMuonTriggerPathToPass_ =
61  parameters.getParameter<std::vector<string> >("muonTriggerPathToPass");
62  theTriggerResultsToken_ = consumes<edm::TriggerResults>(
63  parameters.getParameter<InputTag>("triggerResultsCollection"));
64  theMuonCollectionLabel_ = consumes<reco::MuonCollection>(
65  parameters.getParameter<InputTag>("muonCollection"));
66  theElectronCollectionLabel_ = consumes<reco::GsfElectronCollection>(
67  parameters.getParameter<InputTag>("electronCollection"));
68  thePFJetCollectionToken_ = consumes<edm::View<reco::Jet> >(
69  parameters.getParameter<InputTag>("PFJetCollection"));
70  theCaloMETCollectionToken_ = consumes<edm::View<reco::MET> >(
71  parameters.getParameter<InputTag>("caloMETCollection"));
72  theVertexToken_ = consumes<reco::VertexCollection>(
73  parameters.getParameter<InputTag>("vertexCollection"));
74 
75  // just to initialize
76  isValidHltConfig_ = false;
77 
78  h_vertex_number = 0;
79  h_vertex_chi2 = 0;
80  h_vertex_numTrks = 0;
81  h_vertex_sumTrks = 0;
82  h_vertex_d0 = 0;
83 
84  h_jet_count = 0;
85  h_jet_et = 0;
86  h_jet_pt = 0;
87  h_jet_eta = 0;
88  h_jet_phi = 0;
89  h_jet2_et = 0;
90  // h_jet2_pt = 0;
91  h_jet2_eta = 0;
92  h_jet2_phi = 0;
93 
94  h_e1_et = 0;
95  h_e2_et = 0;
96  h_e1_eta = 0;
97  h_e2_eta = 0;
98  h_e1_phi = 0;
99  h_e2_phi = 0;
100 
101  h_m1_pt = 0;
102  h_m2_pt = 0;
103  h_m1_eta = 0;
104  h_m2_eta = 0;
105  h_m1_phi = 0;
106  h_m2_phi = 0;
107 
108  // h_t1_et = 0;
109  // h_t1_eta = 0;
110  // h_t1_phi = 0;
111 
112  h_met = 0;
113  h_met_phi = 0;
114 
115  h_e_invWMass = 0;
116  h_m_invWMass = 0;
117  h_mumu_invMass = 0;
118  h_ee_invMass = 0;
119 
120  theDbe = Service<DQMStore>().operator->();
121 }
122 
124 }
125 
126 
128  char chtitle[256] = "";
129  const size_t title_s = sizeof(chtitle);
130 
131  logTraceName = "EwkAnalyzer";
132 
133  LogTrace(logTraceName) << "Parameters initialization";
134  theDbe->setCurrentFolder("Physics/EwkDQM"); // Use folder with name of PAG
135 
136  const float pi = 4*atan(1);
137 
138  // Keep the number of plots and number of bins to a minimum!
139  h_vertex_number = theDbe->book1D("vertex_number",
140  "Number of event vertices in collection",
141  10, -0.5, 9.5);
142  h_vertex_chi2 = theDbe->book1D("vertex_chi2",
143  "Event Vertex #chi^{2}/n.d.o.f.",
144  20, 0.0, 2.0);
145  h_vertex_numTrks = theDbe->book1D("vertex_numTrks",
146  "Event Vertex, number of tracks",
147  20, -0.5, 59.5);
148  h_vertex_sumTrks = theDbe->book1D("vertex_sumTrks",
149  "Event Vertex, sum of track pt",
150  20, 0.0, 100.0);
151  h_vertex_d0 = theDbe->book1D("vertex_d0",
152  "Event Vertex d0", 20, 0.0, 0.05);
153 
154  snprintf(chtitle, title_s, "Number of %s (E_{T} > 15 GeV);Number of Jets",
155  thePFJetCollectionLabel_.label().data());
156  h_jet_count = theDbe->book1D("jet_count", chtitle, 8, -0.5, 7.5);
157 
158  snprintf(chtitle, title_s, "Leading jet E_{T} (from %s);E_{T}(1^{st} jet) (GeV)",
159  thePFJetCollectionLabel_.label().data());
160  h_jet_et = theDbe->book1D("jet_et", chtitle, 20, 0., 200.0);
161 
162  snprintf(chtitle, title_s, "Leading jet p_{T} (from %s);p_{T}(1^{st} jet) (GeV/c)",
163  thePFJetCollectionLabel_.label().data());
164  h_jet_pt = theDbe->book1D("jet_pt", chtitle, 20, 0., 200.0);
165 
166  snprintf(chtitle, title_s, "Leading jet #eta (from %s); #eta (1^{st} jet)",
167  thePFJetCollectionLabel_.label().data());
168  h_jet_eta = theDbe->book1D("jet_eta", chtitle, 20, -10., 10.0);
169  snprintf(chtitle, title_s, "Leading jet #phi (from %s); #phi(1^{st} jet)",
170  thePFJetCollectionLabel_.label().data());
171  h_jet_phi = theDbe->book1D("jet_phi", chtitle, 22, -1.1*pi, 1.1*pi);
172 
173  snprintf(chtitle, title_s, "2^{nd} leading jet E_{T} (from %s);E_{T}(2^{nd} jet) (GeV)",
174  thePFJetCollectionLabel_.label().data());
175  h_jet2_et = theDbe->book1D("jet2_et", chtitle, 20, 0., 200.0);
176  // snprintf(chtitle, title_s, "2^{nd} leading jet p_{T} (from %s);p_{T}(2^{nd} jet) (GeV/c)",
177  // thePFJetCollectionLabel_.label().data());
178  // h_jet2_pt = theDbe->book1D("jet2_pt", chtitle, 20, 0., 200.0);
179 
180  snprintf(chtitle, title_s, "2^{nd} leading jet #eta (from %s); #eta (2^{nd} jet)",
181  thePFJetCollectionLabel_.label().data());
182  h_jet2_eta = theDbe->book1D("jet2_eta", chtitle, 20, -10., 10.0);
183 
184  snprintf(chtitle, title_s, "2^{nd} leading jet #phi (from %s); #phi(2^{nd} jet)",
185  thePFJetCollectionLabel_.label().data());
186  h_jet2_phi = theDbe->book1D("jet2_phi", chtitle, 22, -1.1*pi, 1.1*pi);
187 
188  h_e1_et = theDbe->book1D("e1_et", "E_{T} of Leading Electron;E_{T} (GeV)",
189  20, 0.0, 100.0);
190  h_e2_et = theDbe->book1D("e2_et", "E_{T} of Second Electron;E_{T} (GeV)",
191  20, 0.0, 100.0);
192  h_e1_eta = theDbe->book1D("e1_eta", "#eta of Leading Electron;#eta",
193  20, -4.0, 4.0);
194  h_e2_eta = theDbe->book1D("e2_eta", "#eta of Second Electron;#eta",
195  20, -4.0, 4.0);
196  h_e1_phi = theDbe->book1D("e1_phi", "#phi of Leading Electron;#phi",
197  22, -1.1*pi, 1.1*pi);
198  h_e2_phi = theDbe->book1D("e2_phi", "#phi of Second Electron;#phi",
199  22, -1.1*pi, 1.1*pi);
200  h_m1_pt = theDbe->book1D("m1_pt", "p_{T} of Leading Muon;p_{T}(1^{st} #mu) (GeV)",
201  20, 0.0, 100.0);
202  h_m2_pt = theDbe->book1D("m2_pt", "p_{T} of Second Muon;p_{T}(2^{nd} #mu) (GeV)",
203  20, 0.0, 100.0);
204  h_m1_eta = theDbe->book1D("m1_eta", "#eta of Leading Muon;#eta(1^{st} #mu)",
205  20, -4.0, 4.0);
206  h_m2_eta = theDbe->book1D("m2_eta", "#eta of Second Muon;#eta(2^{nd} #mu)",
207  20, -4.0, 4.0);
208  h_m1_phi = theDbe->book1D("m1_phi", "#phi of Leading Muon;#phi(1^{st} #mu)",
209  20, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
210  h_m2_phi = theDbe->book1D("m2_phi", "#phi of Second Muon;#phi(2^{nd} #mu)",
211  20, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
212  // h_t1_et = theDbe->book1D("t1_et", "E_{T} of Leading Tau;E_{T} (GeV)",
213  // 20, 0.0 , 100.0);
214  // h_t1_eta = theDbe->book1D("t1_eta", "#eta of Leading Tau;#eta",
215  // 20, -4.0, 4.0);
216  // h_t1_phi = theDbe->book1D("t1_phi", "#phi of Leading Tau;#phi",
217  // 20, -4.0, 4.0);
218  snprintf(chtitle, title_s, "Missing E_{T} (%s); GeV",
219  theCaloMETCollectionLabel_.label().data());
220  h_met = theDbe->book1D("met", chtitle, 20, 0.0, 100);
221  h_met_phi = theDbe->book1D("met_phi", "Missing E_{T} #phi;#phi(MET)",
222  22, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
223 
224  h_e_invWMass = theDbe->book1D("we_invWMass", "W-> e #nu Transverse Mass;M_{T} (GeV)",
225  20, 0.0, 140.0);
226  h_m_invWMass = theDbe->book1D("wm_invWMass", "W-> #mu #nu Transverse Mass;M_{T} (GeV)",
227  20, 0.0, 140.0);
228  h_mumu_invMass = theDbe->book1D("z_mm_invMass", "#mu#mu Invariant Mass;InvMass (GeV)",
229  20, 40.0, 140.0);
230  h_ee_invMass = theDbe->book1D("z_ee_invMass", "ee Invariant Mass;InvMass (Gev)",
231  20, 40.0, 140.0);
232 }
233 
237 void EwkDQM::beginRun(const edm::Run& theRun, const edm::EventSetup& theSetup) {
238  // passed as parameter to HLTConfigProvider::init(), not yet used
239  bool isConfigChanged = false;
240 
241  // isValidHltConfig_ used to short-circuit analyze() in case of problems
242  const std::string hltProcessName(theTriggerResultsCollection_.process());
243  isValidHltConfig_ = hltConfigProvider_.init(theRun, theSetup,
244  hltProcessName, isConfigChanged);
245 }
246 
247 
248 void EwkDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
249  // short-circuit if hlt problems
250  if (!isValidHltConfig_)
251  return;
252 
253  LogTrace(logTraceName) << "Analysis of event # ";
254  // Did it pass certain HLT path?
255  Handle<TriggerResults> HLTresults;
256  iEvent.getByToken(theTriggerResultsToken_, HLTresults);
257  if (!HLTresults.isValid())
258  return;
259 
260  const edm::TriggerNames & trigNames = iEvent.triggerNames(*HLTresults);
261 
262  // a temporary, until we have a list of triggers of interest
263  std::vector<std::string> eleTrigPathNames;
264  std::vector<std::string> muTrigPathNames;
265 
266  // eleTrigPathNames.push_back(theElecTriggerPathToPass_);
267  // muTrigPathNames.push_back(theMuonTriggerPathToPass_);
268  // end of temporary
269 
270  bool passed_electron_HLT = false;
271  bool passed_muon_HLT = false;
272  for (unsigned int i = 0; i < HLTresults->size(); i++) {
273  const std::string trigName = trigNames.triggerName(i);
274  // check if triggerName matches electronPath
275  for (unsigned int index = 0;
276  index < theElecTriggerPathToPass_.size() && !passed_electron_HLT;
277  index++) {
278  // 0 if found, pos if not
279  size_t trigPath = trigName.find(theElecTriggerPathToPass_[index]);
280  if (trigPath == 0) {
281  // cout << "MuonTrigger passed (=trigName): " << trigName <<endl;
282  passed_electron_HLT = HLTresults->accept(i);
283  }
284  }
285  // check if triggerName matches muonPath
286  for (unsigned int index = 0;
287  index < theMuonTriggerPathToPass_.size() && !passed_muon_HLT;
288  index++) {
289  // 0 if found, pos if not
290  size_t trigPath = trigName.find(theMuonTriggerPathToPass_[index]);
291  if (trigPath == 0) {
292  // cout << "MuonTrigger passed (=trigName): " << trigName <<endl;
293  passed_muon_HLT = HLTresults->accept(i);
294  }
295  }
296  }
297 
298  // we are interested in events with a valid electron or muon
299  if (!(passed_electron_HLT || passed_muon_HLT))
300  return;
301 
303  // Vertex information
304  Handle<VertexCollection> vertexHandle;
305  iEvent.getByToken(theVertexToken_, vertexHandle);
306  if (!vertexHandle.isValid())
307  return;
308  VertexCollection vertexCollection = *(vertexHandle.product());
309  VertexCollection::const_iterator v = vertexCollection.begin();
310  int vertex_number = vertexCollection.size();
311  double vertex_chi2 = v->normalizedChi2(); // v->chi2();
312  double vertex_d0 = sqrt(v->x()*v->x()+v->y()*v->y());
313  double vertex_numTrks = v->tracksSize();
314  double vertex_sumTrks = 0.0;
315  // std::cout << "vertex_d0=" << vertex_d0 << "\n";
316  // double vertex_ndof = v->ndof();cout << "ndof="<<vertex_ndof<<endl;
317  for (Vertex::trackRef_iterator vertex_curTrack = v->tracks_begin();
318  vertex_curTrack != v->tracks_end(); vertex_curTrack++)
319  vertex_sumTrks += (*vertex_curTrack)->pt();
320 
322  // Missing ET
323  Handle< View<MET> > caloMETCollection;
324  iEvent.getByToken(theCaloMETCollectionToken_, caloMETCollection);
325  if (!caloMETCollection.isValid())
326  return;
327  float missing_et = caloMETCollection->begin()->et();
328  float met_phi = caloMETCollection->begin()->phi();
329 
330 
332  // grab "gaussian sum fitting" electrons
333  Handle<GsfElectronCollection> electronCollection;
334  iEvent.getByToken(theElectronCollectionLabel_, electronCollection);
335  if (!electronCollection.isValid())
336  return;
337 
338  // Find the highest and 2nd highest electron
339  float electron_et = -8.0;
340  float electron_eta = -8.0;
341  float electron_phi = -8.0;
342  float electron2_et = -9.0;
343  float electron2_eta = -9.0;
344  float electron2_phi = -9.0;
345  float ee_invMass = -9.0;
346  TLorentzVector e1, e2;
347 
348  // If it passed electron HLT and the collection was found, find electrons near Z mass
349  if (passed_electron_HLT) {
350  for (reco::GsfElectronCollection::const_iterator recoElectron = electronCollection->begin();
351  recoElectron != electronCollection->end(); recoElectron++) {
352  // Require electron to pass some basic cuts
353  if (recoElectron->et() < 20 || fabs(recoElectron->eta()) > 2.5)
354  continue;
355 
356  // Tighter electron cuts
357  if (recoElectron->deltaPhiSuperClusterTrackAtVtx() > 0.58 ||
358  recoElectron->deltaEtaSuperClusterTrackAtVtx() > 0.01 ||
359  recoElectron->sigmaIetaIeta() > 0.027)
360  continue;
361 
362  if (recoElectron->et() > electron_et) {
363  electron2_et = electron_et; // 2nd highest gets values from current highest
364  electron2_eta = electron_eta;
365  electron2_phi = electron_phi;
366  electron_et = recoElectron->et(); // 1st highest gets values from new highest
367  electron_eta = recoElectron->eta();
368  electron_phi = recoElectron->phi();
369  e1 = TLorentzVector(recoElectron->momentum().x(),
370  recoElectron->momentum().y(),
371  recoElectron->momentum().z(),
372  recoElectron->p());
373  } else if (recoElectron->et() > electron2_et) {
374  electron2_et = recoElectron->et();
375  electron2_eta = recoElectron->eta();
376  electron2_phi = recoElectron->phi();
377  e2 = TLorentzVector(recoElectron->momentum().x(),
378  recoElectron->momentum().y(),
379  recoElectron->momentum().z(),
380  recoElectron->p());
381  }
382  } // end of loop over electrons
383  if (electron2_et > 0.0) {
384  TLorentzVector pair = e1+e2;
385  ee_invMass = pair.M();
386  }
387  } // end of "are electrons valid"
389 
390 
391 
393  // Take the STA muon container
394  Handle<MuonCollection> muonCollection;
395  iEvent.getByToken(theMuonCollectionLabel_, muonCollection);
396  if (!muonCollection.isValid())
397  return;
398 
399  // Find the highest pt muons
400  float mm_invMass = -9.0;
401  float muon_pt = -9.0;
402  float muon_eta = -9.0;
403  float muon_phi = -9.0;
404  float muon2_pt = -9.0;
405  float muon2_eta = -9.0;
406  float muon2_phi = -9.0;
407  TLorentzVector m1, m2;
408 
409  if (passed_muon_HLT) {
410  for (reco::MuonCollection::const_iterator recoMuon = muonCollection->begin();
411  recoMuon != muonCollection->end(); recoMuon++) {
412  // Require muon to pass some basic cuts
413  if (recoMuon->pt() < 20 || !recoMuon->isGlobalMuon())
414  continue;
415  // Some tighter muon cuts
416  if (recoMuon->globalTrack()->normalizedChi2() > 10)
417  continue;
418 
419  if (recoMuon->pt() > muon_pt) {
420  muon2_pt = muon_pt; // 2nd highest gets values from current highest
421  muon2_eta = muon_eta;
422  muon2_phi = muon_phi;
423  muon_pt = recoMuon->pt(); // 1st highest gets values from new highest
424  muon_eta = recoMuon->eta();
425  muon_phi = recoMuon->phi();
426  m1 = TLorentzVector(recoMuon->momentum().x(),
427  recoMuon->momentum().y(),
428  recoMuon->momentum().z(),
429  recoMuon->p());
430  } else if (recoMuon->pt() > muon2_pt) {
431  muon2_pt = recoMuon->pt();
432  muon2_eta = recoMuon->eta();
433  muon2_phi = recoMuon->phi();
434  m2 = TLorentzVector(recoMuon->momentum().x(),
435  recoMuon->momentum().y(),
436  recoMuon->momentum().z(),
437  recoMuon->p());
438  }
439  }
440  }
441  if (muon2_pt > 0.0) {
442  TLorentzVector pair = m1+m2;
443  mm_invMass = pair.M();
444  }
446 
447 
449  // Find the highest et jet
450 
451  // Handle<CaloJetCollection> caloJetCollection;
453  iEvent.getByToken(thePFJetCollectionToken_, PFJetCollection);
454  if (!PFJetCollection.isValid())
455  return;
456 
457  unsigned int muonCollectionSize = muonCollection->size();
458  // unsigned int jetCollectionSize = jetCollection->size();
459  unsigned int PFJetCollectionSize = PFJetCollection->size();
460  int jet_count = 0;
461  // int LEADJET=-1; double max_pt=0;
462 
463 
464  float jet_et = -80.0;
465  float jet_pt = -80.0; // prova
466  float jet_eta = -80.0; // now USED
467  float jet_phi = -80.0; // now USED
468  float jet2_et = -90.0;
469  float jet2_eta = -90.0; // now USED
470  float jet2_phi = -90.0; // now USED
471  // for (CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin();
472  // i_calojet != caloJetCollection->end(); i_calojet++) {
473  // for (PFJetCollection::const_iterator i_pfjet = PFJetCollection->begin();
474  // i_pfjet != PFJetCollection->end(); i_pfjet++) {
475  // float jet_current_et = i_calojet->et();
476  // float jet_current_et = i_pfjet->et(); // e` identico a jet.et()
477  // jet_count++;
478 
479  // cleaning: va messo prima del riempimento dell'istogramma // This is in order to use PFJets
480  for (unsigned int i = 0; i < PFJetCollectionSize; i++) {
481  const Jet& jet = PFJetCollection->at(i);
482  // la classe "jet" viene definita qui!!!
483  double minDistance = 99999;
484  for (unsigned int j = 0; j < muonCollectionSize; j++) {
485  const Muon& mu = muonCollection->at(j);
486  double distance = sqrt((mu.eta() - jet.eta()) * (mu.eta() - jet.eta()) +
487  (mu.phi() - jet.phi()) * (mu.phi() - jet.phi()));
488  if (minDistance > distance)
489  minDistance = distance;
490  }
491  if (minDistance < 0.3)
492  continue; // 0.3 is the isolation cone around the muon
493  // se la distanza muone-cono del jet e` minore di 0.3, passo avanti e non conteggio il mio jet
494 
495  // If it overlaps with ELECTRON, it is not a jet
496  if (electron_et > 0.0 &&
497  fabs(jet.eta() - electron_eta) < 0.2 &&
498  calcDeltaPhi(jet.phi(), electron_phi) < 0.2)
499  continue;
500  if (electron2_et > 0.0 &&
501  fabs(jet.eta() - electron2_eta) < 0.2 &&
502  calcDeltaPhi(jet.phi(), electron2_phi) < 0.2)
503  continue;
504 
505  // provo a cambiare la parte degli elettroni in modo simmetrico alla parte per i muoni
506 
507  // ...
508  // ...
509 
510 
511  // if it has too low Et, throw away
512  if (jet.et() < eJetMin_)
513  continue;
514  jet_count++;
515 
516  // ovvero: incrementa jet_count se:
517  // - non c'e un muone entro 0.3 di distanza dal cono del jet;
518  // - se il jet non si sovrappone ad un elettrone;
519  // - se l'energia trasversa e` maggiore della soglia impostata (15?)
520 
521  // if(jet.et()>max_pt) { LEADJET=i; max_pt=jet.et();}
522  // se l'energia del jet e` maggiore di max_pt, diventa "i"
523  // l'indice del jet piu` energetico e max_pt la sua energia
524 
525  // riguardare questo!!!
526  // fino ad ora, jet_et era inizializzato a -8.0
527  if (jet.et() > jet_et) {
528  jet2_et = jet_et; // 2nd highest jet gets et from current highest
529  // perche` prende l'energia del primo jet??
530  jet2_eta = jet_eta; // now USED
531  jet2_phi = jet_phi; // now USED
532  // jet_et = i_calojet->et(); // current highest jet gets
533  // et from the new highest
534  jet_et = jet.et(); // current highest jet gets et from the new highest
535  // ah, ok! lo riaggiorna solo dopo!
536  jet_pt = jet.pt(); // e` il pT del leading jet
537  jet_eta = jet.eta(); // now USED
538  jet_phi = jet.phi() * (Geom::pi() / 180.); // now USED
539  } else if (jet.et() > jet2_et) {
540  // jet2_et = i_calojet->et();
541  jet2_et = jet.et();
542  // jet2_eta = i_calojet->eta(); // UNUSED
543  // jet2_phi = i_calojet->phi(); // UNUSED
544  jet2_eta = jet.eta(); // now USED
545  jet2_phi = jet.phi(); // now USED
546  }
547  // questo elseif funziona
548  }
550 
551 
552 
554  // Fill Histograms //
556 
557  bool fill_e1 = false;
558  bool fill_e2 = false;
559  bool fill_m1 = false;
560  bool fill_m2 = false;
561  bool fill_met = false;
562 
563  // Was Z->ee found?
564  if (ee_invMass > 0.0) {
565  h_ee_invMass->Fill(ee_invMass);
566  fill_e1 = true;
567  fill_e2 = true;
568  }
569 
570  // Was Z->mu mu found?
571  if (mm_invMass > 0.0) {
572  h_mumu_invMass->Fill(mm_invMass);
573  fill_m1 = true;
574  fill_m2 = true;
575  h_jet2_et ->Fill(jet2_et);
576  }
577 
578  // Was W->e nu found?
579  if (electron_et > 0.0 && missing_et > 20.0) {
580  float dphiW = fabs(met_phi-electron_phi);
581  float W_mt_e = sqrt(2 * missing_et * electron_et * (1 - cos(dphiW)));
582  h_e_invWMass->Fill(W_mt_e);
583  fill_e1 = true;
584  fill_met = true;
585  }
586 
587  // Was W->mu nu found?
588  if (muon_pt > 0.0 && missing_et > 20.0) {
589  float dphiW = fabs(met_phi - muon_phi);
590  float W_mt_m = sqrt(2 * missing_et * muon_pt * (1 - cos(dphiW)));
591  h_m_invWMass->Fill(W_mt_m);
592  fill_m1 = true;
593  fill_met = true;
594  }
595 
596  if (jet_et > -10.0) {
597  h_jet_et->Fill(jet_et);
598  h_jet_count->Fill(jet_count);
599  }
600 
601  if (jet_pt > 0.) {
602  h_jet_pt->Fill(jet_pt);
603  }
604 
605  if (jet_eta > -50.) {
606  h_jet_eta->Fill(jet_eta);
607  }
608 
609  if (jet_phi > -10.) {
610  h_jet_phi->Fill(jet_phi);
611  }
612 
613  if (jet2_et > -10.0) {
614  h_jet2_et->Fill(jet2_et);
615  }
616 
617  // if (jet2_pt>0.) {
618  // h_jet2_pt ->Fill(jet2_pt);
619  // }
620 
621  if (jet2_eta > -50.) {
622  h_jet2_eta->Fill(jet2_eta);
623  }
624 
625  if (jet2_phi > -10.) {
626  h_jet2_phi->Fill(jet2_phi);
627  }
628 
629 
630 
631  if (fill_e1 || fill_m1) {
632  h_vertex_number->Fill(vertex_number);
633  h_vertex_chi2->Fill(vertex_chi2);
634  h_vertex_d0->Fill(vertex_d0);
635  h_vertex_numTrks->Fill(vertex_numTrks);
636  h_vertex_sumTrks->Fill(vertex_sumTrks);
637  }
638 
639  if (fill_e1) {
640  h_e1_et->Fill(electron_et);
641  h_e1_eta->Fill(electron_eta);
642  h_e1_phi->Fill(electron_phi);
643  }
644  if (fill_e2) {
645  h_e2_et->Fill(electron2_et);
646  h_e2_eta->Fill(electron2_eta);
647  h_e2_phi->Fill(electron2_phi);
648  }
649  if (fill_m1) {
650  h_m1_pt->Fill(muon_pt);
651  h_m1_eta->Fill(muon_eta);
652  h_m1_phi->Fill(muon_phi);
653  }
654  if (fill_m2) {
655  h_m2_pt->Fill(muon2_pt);
656  h_m2_eta->Fill(muon2_eta);
657  h_m2_phi->Fill(muon2_phi);
658  }
659  if (fill_met) {
660  h_met->Fill(missing_et);
661  h_met_phi->Fill(met_phi);
662  }
664 }
665 
666 
667 void EwkDQM::endJob(void) {}
668 
669 
670 // This always returns only a positive deltaPhi
671 double EwkDQM::calcDeltaPhi(double phi1, double phi2) {
672  double deltaPhi = phi1 - phi2;
673 
674  if (deltaPhi < 0)
675  deltaPhi = -deltaPhi;
676 
677  if (deltaPhi > 3.1415926)
678  deltaPhi = 2 * 3.1415926 - deltaPhi;
679 
680  return deltaPhi;
681 }
682 
683 // Local Variables:
684 // show-trailing-whitespace: t
685 // truncate-lines: t
686 // End:
void beginRun(const edm::Run &, const edm::EventSetup &)
Definition: EwkDQM.cc:237
T getParameter(std::string const &) const
double calcDeltaPhi(double phi1, double phi2)
Definition: EwkTauDQM.cc:972
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:204
virtual float pt() const
transverse momentum
EwkDQM(const edm::ParameterSet &)
Constructor.
Definition: EwkDQM.cc:46
virtual double et() const
transverse energy
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual float phi() const
momentum azimuthal angle
Base class for all types of Jets.
Definition: Jet.h:20
void endJob(void)
Save the histos.
Definition: EwkDQM.cc:667
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
void analyze(const edm::Event &, const edm::EventSetup &)
Get the analysis.
Definition: EwkDQM.cc:248
tuple vertexCollection
const Double_t pi
int iEvent
Definition: GenABIO.cc:230
virtual float eta() const
momentum pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int j
Definition: DBlmapReader.cc:9
void beginJob()
Inizialize parameters for histo binning.
Definition: EwkDQM.cc:127
const int mu
Definition: Constants.h:22
bool isValid() const
Definition: HandleBase.h:76
#define LogTrace(id)
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
virtual ~EwkDQM()
Destructor.
Definition: EwkDQM.cc:123
T const * product() const
Definition: Handle.h:81
std::vector< PFJet > PFJetCollection
collection of PFJet objects
double pi()
Definition: Pi.h:31
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
double calcDeltaPhi(double phi1, double phi2)
Definition: EwkDQM.cc:671
Definition: Run.h:41