CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFMETMonitor.cc
Go to the documentation of this file.
4 
6 
7 #include <TROOT.h>
8 #include <TFile.h>
9 #include <TH1.h>
10 #include <TH2.h>
11 
12 #include <TProfile.h>
13 
14 
15 //
16 // -- Constructor
17 //
19  Benchmark(mode),
20  candBench_(mode),
21  matchCandBench_(mode) {
22 
23  setRange( 0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
24 
25  px_ = 0;
26  sumEt_ = 0;
27  delta_ex_ = 0;
28  delta_ex_VS_set_ = 0;
31 
38 
40  histogramBooked_ = false;
41 }
42 
43 
44 //
45 // -- Destructor
46 //
48 
49 
50 //
51 // -- Set Parameters accessing them from ParameterSet
52 //
54 
55  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
56  createMETSpecificHistos_ = parameterSet.getParameter<bool>( "CreateMETSpecificHistos" );
57  setRange( parameterSet.getParameter<double>("ptMin"),
58  parameterSet.getParameter<double>("ptMax"),
59  parameterSet.getParameter<double>("etaMin"),
60  parameterSet.getParameter<double>("etaMax"),
61  parameterSet.getParameter<double>("phiMin"),
62  parameterSet.getParameter<double>("phiMax") );
63 
64 
67 }
68 
69 
70 //
71 // -- Set Parameters
72 //
74  float etamin, float etamax, float phimin,
75  float phimax, bool metSpHistos) {
76  mode_ = mode;
77  createMETSpecificHistos_ = metSpHistos;
78 
79  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax);
80 
83 }
84 
85 
86 //
87 // -- Create histograms accessing parameters from ParameterSet
88 //
90  candBench_.setup(b, parameterSet);
91  matchCandBench_.setup(b, parameterSet);
92 
94 
95  edm::ParameterSet pxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter");
96  edm::ParameterSet dpxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter");
97  edm::ParameterSet dptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtHistoParameter");
98  edm::ParameterSet setPS = parameterSet.getParameter<edm::ParameterSet>("SumEtHistoParameter");
99  edm::ParameterSet dsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtHistoParameter");
100  edm::ParameterSet setOvsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtOvSumEtHistoParameter");
101 
102  if (pxPS.getParameter<bool>("switchOn")) {
103  px_ = book1D(b, "px_", "px_;p_{X} (GeV)",
104  pxPS.getParameter<int32_t>("nBin"),
105  pxPS.getParameter<double>("xMin"),
106  pxPS.getParameter<double>("xMax"));
107  }
108  if (setPS.getParameter<bool>("switchOn")) {
109  sumEt_ = book1D(b, "sumEt_", "sumEt_;#sumE_{T}",
110  setPS.getParameter<int32_t>("nBin"),
111  setPS.getParameter<double>("xMin"),
112  setPS.getParameter<double>("xMax"));
113  }
114  if (dpxPS.getParameter<bool>("switchOn")) {
115  delta_ex_ = book1D(b, "delta_ex_", "#DeltaME_{X}",
116  dpxPS.getParameter<int32_t>("nBin"),
117  dpxPS.getParameter<double>("xMin"),
118  dpxPS.getParameter<double>("xMax"));
119  }
120 
121  if (dpxPS.getParameter<bool>("switchOn")) {
122  delta_ex_VS_set_ = book2D(b, "delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
123  setPS.getParameter<int32_t>("nBin"),
124  setPS.getParameter<double>("xMin"),
125  setPS.getParameter<double>("xMax"),
126  dptPS.getParameter<int32_t>("nBin"),
127  dptPS.getParameter<double>("xMin"),
128  dptPS.getParameter<double>("xMax"));
129  }
130  if (dsetPS.getParameter<bool>("switchOn")) {
131  delta_set_VS_set_ = book2D(b, "delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
132  setPS.getParameter<int32_t>("nBin"),
133  setPS.getParameter<double>("xMin"),
134  setPS.getParameter<double>("xMax"),
135  dsetPS.getParameter<int32_t>("nBin"),
136  dsetPS.getParameter<double>("xMin"),
137  dsetPS.getParameter<double>("xMax"));
138  }
139  if (setOvsetPS.getParameter<bool>("switchOn")) {
140  delta_set_Over_set_VS_set_ = book2D(b, "delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
141  setPS.getParameter<int32_t>("nBin"),
142  setPS.getParameter<double>("xMin"),
143  setPS.getParameter<double>("xMax"),
144  setOvsetPS.getParameter<int32_t>("nBin"),
145  setOvsetPS.getParameter<double>("xMin"),
146  setOvsetPS.getParameter<double>("xMax"));
147  }
148 
149  // TProfile
150  if (dpxPS.getParameter<bool>("switchOn")) {
151  profile_delta_ex_VS_set_ = bookProfile(b, "profile_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
152  setPS.getParameter<int32_t>("nBin"),
153  setPS.getParameter<double>("xMin"),
154  setPS.getParameter<double>("xMax"),
155  dptPS.getParameter<double>("xMin"),
156  dptPS.getParameter<double>("xMax"), "" );
157  profileRMS_delta_ex_VS_set_ = bookProfile(b, "profileRMS_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
158  setPS.getParameter<int32_t>("nBin"),
159  setPS.getParameter<double>("xMin"),
160  setPS.getParameter<double>("xMax"),
161  dptPS.getParameter<double>("xMin"),
162  dptPS.getParameter<double>("xMax"), "s" );
163  }
164  if (dsetPS.getParameter<bool>("switchOn")) {
165  profile_delta_set_VS_set_ = bookProfile(b, "profile_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
166  setPS.getParameter<int32_t>("nBin"),
167  setPS.getParameter<double>("xMin"),
168  setPS.getParameter<double>("xMax"),
169  dsetPS.getParameter<double>("xMin"),
170  dsetPS.getParameter<double>("xMax"), "" );
171  profileRMS_delta_set_VS_set_ = bookProfile(b, "profileRMS_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
172  setPS.getParameter<int32_t>("nBin"),
173  setPS.getParameter<double>("xMin"),
174  setPS.getParameter<double>("xMax"),
175  dsetPS.getParameter<double>("xMin"),
176  dsetPS.getParameter<double>("xMax"), "s" );
177  }
178  if (setOvsetPS.getParameter<bool>("switchOn")) {
179  profile_delta_set_Over_set_VS_set_ = bookProfile(b, "profile_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
180  setPS.getParameter<int32_t>("nBin"),
181  setPS.getParameter<double>("xMin"),
182  setPS.getParameter<double>("xMax"),
183  setOvsetPS.getParameter<double>("xMin"),
184  setOvsetPS.getParameter<double>("xMax"), "" );
185  profileRMS_delta_set_Over_set_VS_set_ = bookProfile(b, "profileRMS_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
186  setPS.getParameter<int32_t>("nBin"),
187  setPS.getParameter<double>("xMin"),
188  setPS.getParameter<double>("xMax"),
189  setOvsetPS.getParameter<double>("xMin"),
190  setOvsetPS.getParameter<double>("xMax"), "s" );
191  }
192  histogramBooked_ = true;
193  }
194 }
195 
196 
197 //
198 // -- Create histograms using local parameters
199 //
201  candBench_.setup(b);
203 
205 
206  PhaseSpace pxPS = PhaseSpace( 50, 0, 200);
207  PhaseSpace dpxPS = PhaseSpace( 50, -500, 500);
208  PhaseSpace setPS = PhaseSpace( 50, 0.0, 3000);
209  PhaseSpace dsetPS = PhaseSpace( 50, -1000.0, 1000);
210  PhaseSpace setOvsetPS = PhaseSpace( 100,0., 2.);
211 
212  px_ = book1D(b, "px_", "px_;p_{X} (GeV)", pxPS.n, pxPS.m, pxPS.M);
213  sumEt_ = book1D(b, "sumEt_", "sumEt_;#sumE_{T}", setPS.n, setPS.m, setPS.M);
214  delta_ex_ = book1D(b, "delta_ex_", "#DeltaME_{X}", dpxPS.n, dpxPS.m, dpxPS.M);
215  delta_ex_VS_set_ = book2D(b, "delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
216  setPS.n, setPS.m, setPS.M,
217  dpxPS.n, dpxPS.m, dpxPS.M );
218  delta_set_VS_set_ = book2D(b, "delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
219  setPS.n, setPS.m, setPS.M,
220  dsetPS.n, dsetPS.m, dsetPS.M );
221 
222  delta_set_Over_set_VS_set_ = book2D(b, "delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
223  setPS.n, setPS.m, setPS.M,
224  setOvsetPS.n, setOvsetPS.m, setOvsetPS.M );
225 
226  // TProfile
227  profile_delta_ex_VS_set_ = bookProfile(b, "profile_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
228  setPS.n, setPS.m, setPS.M,
229  setOvsetPS.m, setOvsetPS.M, "" );
230 
231  profile_delta_set_VS_set_ = bookProfile(b, "profile_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
232  setPS.n, setPS.m, setPS.M,
233  setOvsetPS.m, setOvsetPS.M, "" );
234 
235  profile_delta_set_Over_set_VS_set_ = bookProfile(b, "profile_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
236  setPS.n, setPS.m, setPS.M,
237  setOvsetPS.m, setOvsetPS.M, "" );
238 
239  // TProfile RMS
240  profileRMS_delta_ex_VS_set_ = bookProfile(b, "profileRMS_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
241  setPS.n, setPS.m, setPS.M,
242  setOvsetPS.m, setOvsetPS.M, "s" );
243 
244  profileRMS_delta_set_VS_set_ = bookProfile(b, "profileRMS_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
245  setPS.n, setPS.m, setPS.M,
246  setOvsetPS.m, setOvsetPS.M, "s" );
247 
248 
249  profileRMS_delta_set_Over_set_VS_set_ = bookProfile(b, "profileRMS_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
250  setPS.n, setPS.m, setPS.M,
251  setOvsetPS.m, setOvsetPS.M, "s" );
252  histogramBooked_ = true;
253  }
254 }
255 
256 
257 void PFMETMonitor::setDirectory(TDirectory* dir) {
259 
262 }
263 
264 
266  const reco::MET& matchedMet, float& minVal, float& maxVal) {
267 
268  candBench_.fillOne(met); //std::cout <<"\nfillone MET candBench" <<std::endl;
269  matchCandBench_.fillOne(met, matchedMet); //std::cout <<"\nfillone MET MatchCandBench done" <<std::endl;
270 
272  if( !isInRange(met.pt(), met.eta(), met.phi() ) ) return;
273 
274  if (px_) px_->Fill(met.px());
275  if (delta_ex_) {
276  delta_ex_->Fill( met.px() - matchedMet.px() );
277  delta_ex_->Fill( met.py() - matchedMet.py() );
278  }
279  if (sumEt_) sumEt_->Fill( met.sumEt() );
280 
281  if (delta_ex_VS_set_) {
282  delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
283  delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
284  profile_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
285  profile_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
286  profileRMS_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
287  profileRMS_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
288  }
289  if (delta_set_VS_set_) {
290  delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
291  profile_delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
292  profileRMS_delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
293  }
294  if (delta_set_Over_set_VS_set_ && matchedMet.sumEt()>0.001 ) {
295  float setRes = (met.sumEt() - matchedMet.sumEt()) / matchedMet.sumEt();
296  if (setRes > maxVal) maxVal = setRes;
297  if (setRes < minVal) minVal = setRes;
298  delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
299  profile_delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
300  profileRMS_delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
301  }
302  }
303 }
TH1F * sumEt_
Definition: PFMETMonitor.h:42
TProfile * profile_delta_ex_VS_set_
Definition: PFMETMonitor.h:48
T getParameter(std::string const &) const
void setParameters(Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool metSpHistos)
set the parameters locally
Definition: PFMETMonitor.cc:73
TH1F * book1D(DQMStore::IBooker &b, const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either with DQM or plain root depending if DQM_ has been initialized in a child ...
Definition: Benchmark.cc:23
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
TProfile * profileRMS_delta_set_VS_set_
Definition: PFMETMonitor.h:53
TProfile * profile_delta_set_VS_set_
Definition: PFMETMonitor.h:49
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
TH2F * delta_ex_VS_set_
Definition: PFMETMonitor.h:44
TH2F * delta_set_VS_set_
Definition: PFMETMonitor.h:45
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
abstract base class
Definition: Benchmark.h:22
void setup(DQMStore::IBooker &b)
book histograms
bool histogramBooked_
Definition: PFMETMonitor.h:61
TH1F * delta_ex_
Definition: PFMETMonitor.h:43
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:18
TProfile * profileRMS_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:54
TProfile * profileRMS_delta_ex_VS_set_
Definition: PFMETMonitor.h:52
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
double sumEt() const
Definition: MET.h:56
Definition: MET.h:42
PFMETMonitor(Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFMETMonitor.cc:18
void setParameters(Mode mode)
Definition: Benchmark.h:49
TH2F * book2D(DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float xmin, float xmax, int nbinsy, float ymin, float ymax)
book a 2D histogram, either with DQM or plain root depending if DQM_ has been initialized in a child ...
Definition: Benchmark.cc:32
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:57
void fillOne(const reco::MET &met, const reco::MET &matchedMet, float &minVal, float &maxVal)
virtual ~PFMETMonitor()
Definition: PFMETMonitor.cc:47
TH2F * delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:46
virtual double px() const
x coordinate of momentum vector
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:60
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:51
double b
Definition: hdecay.h:120
TProfile * profile_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:50
double ptmin
Definition: HydjetWrapper.h:85
Mode mode_
Definition: Benchmark.h:105
void setup(DQMStore::IBooker &b)
book histograms
dbl *** dir
Definition: mlp_gen.cc:35
void setup(DQMStore::IBooker &b)
book histograms
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:58
virtual double phi() const
momentum azimuthal angle
TProfile * bookProfile(DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float xmin, float xmax, float ymin, float ymax, const char *option)
book a TProfile histogram, either with DQM or plain root depending if DQM_ has been initialized in a ...
Definition: Benchmark.cc:60
virtual double py() const
y coordinate of momentum vector
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11