CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
MatchMETBenchmark.cc
Go to the documentation of this file.
5 
6 #include <TFile.h>
7 #include <TH1.h>
8 #include <TH2.h>
9 #include <TROOT.h>
10 
11 using namespace std;
12 
14 
16  // std::cout << "FL: MatchMETBenchmark.cc: start setup()" << std::endl;
17  PhaseSpace ptPS;
18  PhaseSpace dptOvptPS;
19  PhaseSpace dptPS;
20  PhaseSpace dphiPS;
21  PhaseSpace setPS;
22  PhaseSpace dsetPS;
23  PhaseSpace setOvsetPS;
24 
25  switch (mode_) {
26  case VALIDATION:
27  ptPS = PhaseSpace(100, 0, 1000);
28  dptOvptPS = PhaseSpace(200, -1, 1);
29  dphiPS = PhaseSpace(100, -3.2, 3.2);
30  dptPS = PhaseSpace(200, -500, 500);
31  setPS = PhaseSpace(300, 0.0, 3000);
32  dsetPS = PhaseSpace(200, 0. - 1000, 1000);
33  setOvsetPS = PhaseSpace(500, 0., 2.);
34  break;
35  case DQMOFFLINE:
36  default:
37  ptPS = PhaseSpace(50, 0, 200);
38  dptOvptPS = PhaseSpace(50, -1, 1);
39  dphiPS = PhaseSpace(50, -3.2, 3.2);
40  dptPS = PhaseSpace(50, -500, 500);
41  setPS = PhaseSpace(50, 0.0, 3000);
42  dsetPS = PhaseSpace(50, -1000.0, 1000);
43  setOvsetPS = PhaseSpace(100, 0., 2.);
44  break;
45  }
46 
47  // variable bins to be done here, as they will save a lot of memory.
48 
49  // float ptBins[11] = {0, 1, 2, 5, 10, 20, 50, 100, 200, 400, 1000};
50 
51  delta_et_Over_et_VS_et_ = book2D(b,
52  "delta_et_Over_et_VS_et_",
53  ";ME_{T, true} (GeV);#DeltaME_{T}/ME_{T}",
54  ptPS.n,
55  ptPS.m,
56  ptPS.M,
57  dptOvptPS.n,
58  dptOvptPS.m,
59  dptOvptPS.M);
60 
61  delta_et_VS_et_ = book2D(
62  b, "delta_et_VS_et_", ";ME_{T, true} (GeV);#DeltaME_{T}", ptPS.n, ptPS.m, ptPS.M, dptPS.n, dptPS.m, dptPS.M);
63 
64  delta_phi_VS_et_ = book2D(
65  b, "delta_phi_VS_et_", ";ME_{T, true} (GeV);#Delta#phi", ptPS.n, ptPS.m, ptPS.M, dphiPS.n, dphiPS.m, dphiPS.M);
66 
67  delta_ex_ = book1D(b, "delta_ex_", "#DeltaME_{X}", dptPS.n, dptPS.m, dptPS.M);
68 
69  RecEt_VS_TrueEt_ =
70  book2D(b, "RecEt_VS_TrueEt_", ";ME_{T, true} (GeV);ME_{T}", ptPS.n, ptPS.m, ptPS.M, ptPS.n, ptPS.m, ptPS.M);
71 
72  delta_set_VS_set_ = book2D(b,
73  "delta_set_VS_set_",
74  ";SE_{T, true} (GeV);#DeltaSE_{T}",
75  setPS.n,
76  setPS.m,
77  setPS.M,
78  dsetPS.n,
79  dsetPS.m,
80  dsetPS.M);
81 
82  delta_set_Over_set_VS_set_ = book2D(b,
83  "delta_set_Over_set_VS_set_",
84  ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
85  setPS.n,
86  setPS.m,
87  setPS.M,
88  dptOvptPS.n,
89  dptOvptPS.m,
90  dptOvptPS.M);
91 
92  delta_ex_VS_set_ = book2D(
93  b, "delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}", setPS.n, setPS.m, setPS.M, ptPS.n, -ptPS.M, ptPS.M);
94 
95  RecSet_Over_TrueSet_VS_TrueSet_ = book2D(b,
96  "RecSet_Over_TrueSet_VS_TrueSet_",
97  ";SE_{T, true} (GeV);SE_{T}/SE_{T}",
98  setPS.n,
99  setPS.m,
100  setPS.M,
101  setOvsetPS.n,
102  setOvsetPS.m,
103  setOvsetPS.M);
104 }
105 
106 void MatchMETBenchmark::fillOne(const reco::MET &cand, const reco::MET &matchedCand) {
107  if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
108  return;
109 
110  if (matchedCand.pt() > 0.001)
111  delta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
112  else
113  edm::LogWarning("MatchMETBenchmark") << " matchedCand.pt()<0.001";
114  delta_et_VS_et_->Fill(matchedCand.pt(), cand.pt() - matchedCand.pt());
115  delta_phi_VS_et_->Fill(matchedCand.pt(), cand.phi() - matchedCand.phi());
116  delta_ex_->Fill(cand.px() - matchedCand.px());
117  delta_ex_->Fill(cand.py() - matchedCand.py());
118  RecEt_VS_TrueEt_->Fill(matchedCand.pt(), cand.pt());
119  delta_set_VS_set_->Fill(matchedCand.sumEt(), cand.sumEt() - matchedCand.sumEt());
120  if (matchedCand.sumEt() > 0.001)
121  delta_set_Over_set_VS_set_->Fill(matchedCand.sumEt(), (cand.sumEt() - matchedCand.sumEt()) / matchedCand.sumEt());
122  else
123  edm::LogWarning("MatchMETBenchmark") << " matchedCand.sumEt()<0.001";
124  delta_ex_VS_set_->Fill(matchedCand.sumEt(), cand.px() - matchedCand.px());
125  delta_ex_VS_set_->Fill(matchedCand.sumEt(), cand.py() - matchedCand.py());
126  if (matchedCand.sumEt() > 0.001)
127  RecSet_Over_TrueSet_VS_TrueSet_->Fill(matchedCand.sumEt(), cand.sumEt() / matchedCand.sumEt());
128 }
double pt() const final
transverse momentum
void setup(DQMStore::IBooker &b)
book histograms
~MatchMETBenchmark() override
double sumEt() const
Definition: MET.h:56
void fillOne(const reco::MET &candidate, const reco::MET &matchedCandidate)
fill histograms with a given particle
double px() const final
x coordinate of momentum vector
Definition: MET.h:41
double py() const final
y coordinate of momentum vector
double b
Definition: hdecay.h:120
Log< level::Warning, false > LogWarning
double phi() const final
momentum azimuthal angle