CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes
PFMETMonitor Class Reference

#include <PFMETMonitor.h>

Inheritance diagram for PFMETMonitor:
Benchmark

Public Member Functions

void fillOne (const reco::MET &met, const reco::MET &matchedMet, float &minVal, float &maxVal)
 
void fillOne (const reco::MET &met, const reco::MET &matchedMet, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet)
 
 PFMETMonitor (Benchmark::Mode mode=Benchmark::DEFAULT)
 
void setDirectory (TDirectory *dir)
 set directory (to use in ROOT) More...
 
void setParameters (Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool metSpHistos)
 set the parameters locally More...
 
void setParameters (const edm::ParameterSet &parameterSet)
 set the parameters accessing them from ParameterSet More...
 
void setup ()
 book histograms More...
 
void setup (const edm::ParameterSet &parameterSet)
 book histograms More...
 
virtual ~PFMETMonitor ()
 
- Public Member Functions inherited from Benchmark
 Benchmark (Mode mode=DEFAULT)
 
bool isInRange (float pt, float eta, float phi) const
 
void setParameters (Mode mode)
 
void setRange (float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
 
void write ()
 write to the TFile, in plain ROOT mode. No need to call this function in DQM mode More...
 
virtual ~Benchmark ()
 

Protected Attributes

CandidateBenchmark candBench_
 
bool createMETSpecificHistos_
 
TH1F * delta_ex_
 
TH2F * delta_ex_VS_set_
 
TH2F * delta_set_Over_set_VS_set_
 
TH2F * delta_set_VS_set_
 
bool histogramBooked_
 
MatchCandidateBenchmark matchCandBench_
 
TProfile * profile_delta_ex_VS_set_
 
TProfile * profile_delta_set_Over_set_VS_set_
 
TProfile * profile_delta_set_VS_set_
 
TProfile * profileRMS_delta_ex_VS_set_
 
TProfile * profileRMS_delta_set_Over_set_VS_set_
 
TProfile * profileRMS_delta_set_VS_set_
 
TH1F * px_
 
TH1F * sumEt_
 
- Protected Attributes inherited from Benchmark
TDirectory * dir_
 
float etaMax_
 
float etaMin_
 
Mode mode_
 
float phiMax_
 
float phiMin_
 
float ptMax_
 
float ptMin_
 

Additional Inherited Members

- Public Types inherited from Benchmark
enum  Mode { DEFAULT, DQMOFFLINE, VALIDATION }
 
- Static Public Attributes inherited from Benchmark
static DQMStoreDQM_ = 0
 
- Protected Member Functions inherited from Benchmark
TH1F * book1D (const char *histname, const char *title, int nbins, float xmin, float xmax)
 book a 1D histogram, either with DQM or plain root. More...
 
TH2F * book2D (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. More...
 
TH2F * book2D (const char *histname, const char *title, int nbinsx, float *xbins, int nbinsy, float ymin, float ymax)
 book a 2D histogram, either with DQM or plain root. More...
 
TProfile * bookProfile (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. More...
 
TProfile * bookProfile (const char *histname, const char *title, int nbinsx, float *xbins, float ymin, float ymax, const char *option)
 book a TProfile histogram, either with DQM or plain root. More...
 

Detailed Description

Definition at line 11 of file PFMETMonitor.h.

Constructor & Destructor Documentation

PFMETMonitor::PFMETMonitor ( Benchmark::Mode  mode = Benchmark::DEFAULT)

Definition at line 17 of file PFMETMonitor.cc.

References createMETSpecificHistos_, delta_ex_, delta_ex_VS_set_, delta_set_Over_set_VS_set_, delta_set_VS_set_, histogramBooked_, profile_delta_ex_VS_set_, profile_delta_set_Over_set_VS_set_, profile_delta_set_VS_set_, profileRMS_delta_ex_VS_set_, profileRMS_delta_set_Over_set_VS_set_, profileRMS_delta_set_VS_set_, px_, Benchmark::setRange(), and sumEt_.

17  :
18  Benchmark(mode),
19  candBench_(mode),
21 
22  setRange( 0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
23 
24  px_ = 0;
25  sumEt_ = 0;
26  delta_ex_ = 0;
27  delta_ex_VS_set_ = 0;
30 
37 
39  histogramBooked_ = false;
40 }
TH1F * sumEt_
Definition: PFMETMonitor.h:44
TProfile * profile_delta_ex_VS_set_
Definition: PFMETMonitor.h:50
TProfile * profileRMS_delta_set_VS_set_
Definition: PFMETMonitor.h:55
TProfile * profile_delta_set_VS_set_
Definition: PFMETMonitor.h:51
TH2F * delta_ex_VS_set_
Definition: PFMETMonitor.h:46
TH2F * delta_set_VS_set_
Definition: PFMETMonitor.h:47
bool histogramBooked_
Definition: PFMETMonitor.h:63
TH1F * delta_ex_
Definition: PFMETMonitor.h:45
TProfile * profileRMS_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:56
TProfile * profileRMS_delta_ex_VS_set_
Definition: PFMETMonitor.h:54
Benchmark(Mode mode=DEFAULT)
Definition: Benchmark.h:42
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
TH2F * delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:48
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
TProfile * profile_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:52
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
PFMETMonitor::~PFMETMonitor ( )
virtual

Definition at line 44 of file PFMETMonitor.cc.

44 {}

Member Function Documentation

void PFMETMonitor::fillOne ( const reco::MET met,
const reco::MET matchedMet,
float &  minVal,
float &  maxVal 
)

Definition at line 255 of file PFMETMonitor.cc.

References candBench_, createMETSpecificHistos_, delta_ex_, delta_ex_VS_set_, delta_set_Over_set_VS_set_, delta_set_VS_set_, reco::LeafCandidate::eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), histogramBooked_, Benchmark::isInRange(), matchCandBench_, reco::LeafCandidate::phi(), profile_delta_ex_VS_set_, profile_delta_set_Over_set_VS_set_, profile_delta_set_VS_set_, profileRMS_delta_ex_VS_set_, profileRMS_delta_set_Over_set_VS_set_, profileRMS_delta_set_VS_set_, reco::LeafCandidate::pt(), reco::LeafCandidate::px(), px_, reco::LeafCandidate::py(), reco::MET::sumEt(), and sumEt_.

Referenced by PFMETDQMAnalyzer::analyze(), and PFRootEventManager::processEntry().

256  {
257  /*void PFMETMonitor::fillOne(const reco::MET& met,
258  const reco::MET& matchedMet, float& minVal, float& maxVal,
259  const edm::ParameterSet & parameterSet) {*/
260 
261  candBench_.fillOne(met); //std::cout <<"\nfillone MET candBench" <<std::endl;
262  matchCandBench_.fillOne(met, matchedMet); //std::cout <<"\nfillone MET MatchCandBench done" <<std::endl;
263  //matchCandBench_.fillOne(met, matchedMet, parameterSet); std::cout <<"\nfillone MET MatchCandBench done" <<std::endl;
264 
266  if( !isInRange(met.pt(), met.eta(), met.phi() ) ) return;
267 
268  if (px_) px_->Fill(met.px());
269  if (delta_ex_) {
270  delta_ex_->Fill( met.px() - matchedMet.px() );
271  delta_ex_->Fill( met.py() - matchedMet.py() );
272  }
273  if (sumEt_) sumEt_->Fill( met.sumEt() );
274 
275  if (delta_ex_VS_set_) {
276  delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
277  delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
278  profile_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
279  profile_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
280  profileRMS_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.px() - matchedMet.px() );
281  profileRMS_delta_ex_VS_set_->Fill( matchedMet.sumEt(), met.py() - matchedMet.py() );
282  }
283  if (delta_set_VS_set_) {
284  delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
285  profile_delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
286  profileRMS_delta_set_VS_set_->Fill( matchedMet.sumEt(), met.sumEt() - matchedMet.sumEt() );
287  }
288  if (delta_set_Over_set_VS_set_ && matchedMet.sumEt()>0.001 ) {
289  float setRes = (met.sumEt() - matchedMet.sumEt()) / matchedMet.sumEt();
290  if (setRes > maxVal) maxVal = setRes;
291  if (setRes < minVal) minVal = setRes;
292  delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
293  profile_delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
294  profileRMS_delta_set_Over_set_VS_set_->Fill( matchedMet.sumEt(), setRes );
295  }
296  }
297 }
TH1F * sumEt_
Definition: PFMETMonitor.h:44
TProfile * profile_delta_ex_VS_set_
Definition: PFMETMonitor.h:50
TProfile * profileRMS_delta_set_VS_set_
Definition: PFMETMonitor.h:55
TProfile * profile_delta_set_VS_set_
Definition: PFMETMonitor.h:51
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
TH2F * delta_ex_VS_set_
Definition: PFMETMonitor.h:46
TH2F * delta_set_VS_set_
Definition: PFMETMonitor.h:47
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
bool histogramBooked_
Definition: PFMETMonitor.h:63
TH1F * delta_ex_
Definition: PFMETMonitor.h:45
TProfile * profileRMS_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:56
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
TProfile * profileRMS_delta_ex_VS_set_
Definition: PFMETMonitor.h:54
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
double sumEt() const
Definition: MET.h:48
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
TH2F * delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:48
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
TProfile * profile_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:52
virtual float pt() const GCC11_FINAL
transverse momentum
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:59
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
void PFMETMonitor::fillOne ( const reco::MET met,
const reco::MET matchedMet,
float &  minVal,
float &  maxVal,
const edm::ParameterSet parameterSet 
)
void PFMETMonitor::setDirectory ( TDirectory *  dir)
virtual

set directory (to use in ROOT)

Reimplemented from Benchmark.

Definition at line 248 of file PFMETMonitor.cc.

References candBench_, matchCandBench_, and Benchmark::setDirectory().

Referenced by PFRootEventManager::readOptions().

248  {
250 
253 }
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:19
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
dbl *** dir
Definition: mlp_gen.cc:35
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
void PFMETMonitor::setParameters ( Benchmark::Mode  mode,
float  ptmin,
float  ptmax,
float  etamin,
float  etamax,
float  phimin,
float  phimax,
bool  metSpHistos 
)

set the parameters locally

Definition at line 67 of file PFMETMonitor.cc.

References candBench_, createMETSpecificHistos_, matchCandBench_, alignBH_cfg::mode, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

Referenced by PFMETDQMAnalyzer::PFMETDQMAnalyzer(), and PFRootEventManager::readOptions().

69  {
70  mode_ = mode;
71  createMETSpecificHistos_ = metSpHistos;
72 
73  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax);
74 
77 }
void setParameters(Mode mode)
Definition: Benchmark.h:50
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
double ptmin
Definition: HydjetWrapper.h:85
Mode mode_
Definition: Benchmark.h:96
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
void PFMETMonitor::setParameters ( const edm::ParameterSet parameterSet)

set the parameters accessing them from ParameterSet

Definition at line 49 of file PFMETMonitor.cc.

References candBench_, createMETSpecificHistos_, edm::ParameterSet::getParameter(), matchCandBench_, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

49  {
50 
51  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
52  createMETSpecificHistos_ = parameterSet.getParameter<bool>( "CreateMETSpecificHistos" );
53  setRange( parameterSet.getParameter<double>("ptMin"),
54  parameterSet.getParameter<double>("ptMax"),
55  parameterSet.getParameter<double>("etaMin"),
56  parameterSet.getParameter<double>("etaMax"),
57  parameterSet.getParameter<double>("phiMin"),
58  parameterSet.getParameter<double>("phiMax") );
59 
60 
63 }
T getParameter(std::string const &) const
void setParameters(Mode mode)
Definition: Benchmark.h:50
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
Mode mode_
Definition: Benchmark.h:96
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
void PFMETMonitor::setup ( void  )

book histograms

Definition at line 189 of file PFMETMonitor.cc.

References Benchmark::book1D(), Benchmark::book2D(), Benchmark::bookProfile(), candBench_, createMETSpecificHistos_, delta_ex_, delta_ex_VS_set_, delta_set_Over_set_VS_set_, delta_set_VS_set_, histogramBooked_, Benchmark::PhaseSpace::m, Benchmark::PhaseSpace::M, matchCandBench_, Benchmark::PhaseSpace::n, profile_delta_ex_VS_set_, profile_delta_set_Over_set_VS_set_, profile_delta_set_VS_set_, profileRMS_delta_ex_VS_set_, profileRMS_delta_set_Over_set_VS_set_, profileRMS_delta_set_VS_set_, px_, CandidateBenchmark::setup(), MatchCandidateBenchmark::setup(), and sumEt_.

Referenced by PFMETDQMAnalyzer::beginJob(), and PFRootEventManager::readOptions().

189  {
190  candBench_.setup();
192 
194 
195  PhaseSpace pxPS = PhaseSpace( 50, 0, 200);
196  PhaseSpace dpxPS = PhaseSpace( 50, -500, 500);
197  PhaseSpace setPS = PhaseSpace( 50, 0.0, 3000);
198  PhaseSpace dsetPS = PhaseSpace( 50, -1000.0, 1000);
199  PhaseSpace setOvsetPS = PhaseSpace( 100,0., 2.);
200 
201  px_ = book1D("px_", "px_;p_{X} (GeV)", pxPS.n, pxPS.m, pxPS.M);
202 
203  sumEt_ = book1D("sumEt_", "sumEt_;#sumE_{T}", setPS.n, setPS.m, setPS.M);
204 
205  delta_ex_ = book1D("delta_ex_", "#DeltaME_{X}", dpxPS.n, dpxPS.m, dpxPS.M);
206 
207  delta_ex_VS_set_ = book2D("delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
208  setPS.n, setPS.m, setPS.M,
209  dpxPS.n, dpxPS.m, dpxPS.M );
210 
211  delta_set_VS_set_ = book2D("delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
212  setPS.n, setPS.m, setPS.M,
213  dsetPS.n, dsetPS.m, dsetPS.M );
214 
215  delta_set_Over_set_VS_set_ = book2D("delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
216  setPS.n, setPS.m, setPS.M,
217  setOvsetPS.n, setOvsetPS.m, setOvsetPS.M );
218  // TProfile
219  profile_delta_ex_VS_set_ = bookProfile("profile_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
220  setPS.n, setPS.m, setPS.M,
221  setOvsetPS.m, setOvsetPS.M, "" );
222 
223  profile_delta_set_VS_set_ = bookProfile("profile_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
224  setPS.n, setPS.m, setPS.M,
225  setOvsetPS.m, setOvsetPS.M, "" );
226 
227 
228  profile_delta_set_Over_set_VS_set_ = bookProfile("profile_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
229  setPS.n, setPS.m, setPS.M,
230  setOvsetPS.m, setOvsetPS.M, "" );
231  // TProfile RMS
232  profileRMS_delta_ex_VS_set_ = bookProfile("profileRMS_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
233  setPS.n, setPS.m, setPS.M,
234  setOvsetPS.m, setOvsetPS.M, "s" );
235 
236  profileRMS_delta_set_VS_set_ = bookProfile("profileRMS_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
237  setPS.n, setPS.m, setPS.M,
238  setOvsetPS.m, setOvsetPS.M, "s" );
239 
240 
241  profileRMS_delta_set_Over_set_VS_set_ = bookProfile("profileRMS_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
242  setPS.n, setPS.m, setPS.M,
243  setOvsetPS.m, setOvsetPS.M, "s" );
244  histogramBooked_ = true;
245  }
246 }
TH1F * sumEt_
Definition: PFMETMonitor.h:44
TProfile * profile_delta_ex_VS_set_
Definition: PFMETMonitor.h:50
TProfile * profileRMS_delta_set_VS_set_
Definition: PFMETMonitor.h:55
TProfile * profile_delta_set_VS_set_
Definition: PFMETMonitor.h:51
TH2F * delta_ex_VS_set_
Definition: PFMETMonitor.h:46
TH2F * delta_set_VS_set_
Definition: PFMETMonitor.h:47
void setup()
book histograms
void setup()
book histograms
bool histogramBooked_
Definition: PFMETMonitor.h:63
TH1F * delta_ex_
Definition: PFMETMonitor.h:45
TProfile * profileRMS_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:56
TProfile * profileRMS_delta_ex_VS_set_
Definition: PFMETMonitor.h:54
TH2F * book2D(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.
Definition: Benchmark.cc:43
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
TH2F * delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:48
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
TProfile * profile_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:52
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
TProfile * bookProfile(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.
Definition: Benchmark.cc:95
TH1F * book1D(const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either with DQM or plain root.
Definition: Benchmark.cc:25
void PFMETMonitor::setup ( const edm::ParameterSet parameterSet)

book histograms

Definition at line 81 of file PFMETMonitor.cc.

References Benchmark::book1D(), Benchmark::book2D(), Benchmark::bookProfile(), candBench_, createMETSpecificHistos_, delta_ex_, delta_ex_VS_set_, delta_set_Over_set_VS_set_, delta_set_VS_set_, edm::ParameterSet::getParameter(), histogramBooked_, matchCandBench_, profile_delta_ex_VS_set_, profile_delta_set_Over_set_VS_set_, profile_delta_set_VS_set_, profileRMS_delta_ex_VS_set_, profileRMS_delta_set_Over_set_VS_set_, profileRMS_delta_set_VS_set_, px_, CandidateBenchmark::setup(), MatchCandidateBenchmark::setup(), and sumEt_.

81  {
82  candBench_.setup(parameterSet);
83  matchCandBench_.setup(parameterSet);
84 
86 
87  edm::ParameterSet pxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter");
88  edm::ParameterSet dpxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter");
89  edm::ParameterSet dptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtHistoParameter");
90  edm::ParameterSet setPS = parameterSet.getParameter<edm::ParameterSet>("SumEtHistoParameter");
91  edm::ParameterSet dsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtHistoParameter");
92  edm::ParameterSet setOvsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtOvSumEtHistoParameter");
93 
94  if (pxPS.getParameter<bool>("switchOn")) {
95  px_ = book1D("px_", "px_;p_{X} (GeV)",
96  pxPS.getParameter<int32_t>("nBin"),
97  pxPS.getParameter<double>("xMin"),
98  pxPS.getParameter<double>("xMax"));
99  }
100  if (setPS.getParameter<bool>("switchOn")) {
101  sumEt_ = book1D("sumEt_", "sumEt_;#sumE_{T}",
102  setPS.getParameter<int32_t>("nBin"),
103  setPS.getParameter<double>("xMin"),
104  setPS.getParameter<double>("xMax"));
105  }
106  if (dpxPS.getParameter<bool>("switchOn")) {
107  delta_ex_ = book1D("delta_ex_", "#DeltaME_{X}",
108  dpxPS.getParameter<int32_t>("nBin"),
109  dpxPS.getParameter<double>("xMin"),
110  dpxPS.getParameter<double>("xMax"));
111  }
112 
113  if (dpxPS.getParameter<bool>("switchOn")) {
114  delta_ex_VS_set_ = book2D("delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
115  setPS.getParameter<int32_t>("nBin"),
116  setPS.getParameter<double>("xMin"),
117  setPS.getParameter<double>("xMax"),
118  dptPS.getParameter<int32_t>("nBin"),
119  dptPS.getParameter<double>("xMin"),
120  dptPS.getParameter<double>("xMax"));
121  }
122  if (dsetPS.getParameter<bool>("switchOn")) {
123  delta_set_VS_set_ = book2D("delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
124  setPS.getParameter<int32_t>("nBin"),
125  setPS.getParameter<double>("xMin"),
126  setPS.getParameter<double>("xMax"),
127  dsetPS.getParameter<int32_t>("nBin"),
128  dsetPS.getParameter<double>("xMin"),
129  dsetPS.getParameter<double>("xMax"));
130  }
131  if (setOvsetPS.getParameter<bool>("switchOn")) {
132  delta_set_Over_set_VS_set_ = book2D("delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
133  setPS.getParameter<int32_t>("nBin"),
134  setPS.getParameter<double>("xMin"),
135  setPS.getParameter<double>("xMax"),
136  setOvsetPS.getParameter<int32_t>("nBin"),
137  setOvsetPS.getParameter<double>("xMin"),
138  setOvsetPS.getParameter<double>("xMax"));
139  }
140  // TProfile
141  if (dpxPS.getParameter<bool>("switchOn")) {
142  profile_delta_ex_VS_set_ = bookProfile("profile_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
143  setPS.getParameter<int32_t>("nBin"),
144  setPS.getParameter<double>("xMin"),
145  setPS.getParameter<double>("xMax"),
146  dptPS.getParameter<double>("xMin"),
147  dptPS.getParameter<double>("xMax"), "" );
148  profileRMS_delta_ex_VS_set_ = bookProfile("profileRMS_delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}",
149  setPS.getParameter<int32_t>("nBin"),
150  setPS.getParameter<double>("xMin"),
151  setPS.getParameter<double>("xMax"),
152  dptPS.getParameter<double>("xMin"),
153  dptPS.getParameter<double>("xMax"), "s" );
154  }
155  if (dsetPS.getParameter<bool>("switchOn")) {
156  profile_delta_set_VS_set_ = bookProfile("profile_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
157  setPS.getParameter<int32_t>("nBin"),
158  setPS.getParameter<double>("xMin"),
159  setPS.getParameter<double>("xMax"),
160  dsetPS.getParameter<double>("xMin"),
161  dsetPS.getParameter<double>("xMax"), "" );
162  profileRMS_delta_set_VS_set_ = bookProfile("profileRMS_delta_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}",
163  setPS.getParameter<int32_t>("nBin"),
164  setPS.getParameter<double>("xMin"),
165  setPS.getParameter<double>("xMax"),
166  dsetPS.getParameter<double>("xMin"),
167  dsetPS.getParameter<double>("xMax"), "s" );
168  }
169  if (setOvsetPS.getParameter<bool>("switchOn")) {
170  profile_delta_set_Over_set_VS_set_ = bookProfile("profile_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
171  setPS.getParameter<int32_t>("nBin"),
172  setPS.getParameter<double>("xMin"),
173  setPS.getParameter<double>("xMax"),
174  setOvsetPS.getParameter<double>("xMin"),
175  setOvsetPS.getParameter<double>("xMax"), "" );
176  profileRMS_delta_set_Over_set_VS_set_ = bookProfile("profileRMS_delta_set_Over_set_VS_set_", ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
177  setPS.getParameter<int32_t>("nBin"),
178  setPS.getParameter<double>("xMin"),
179  setPS.getParameter<double>("xMax"),
180  setOvsetPS.getParameter<double>("xMin"),
181  setOvsetPS.getParameter<double>("xMax"), "s" );
182  }
183  histogramBooked_ = true;
184  }
185 }
TH1F * sumEt_
Definition: PFMETMonitor.h:44
TProfile * profile_delta_ex_VS_set_
Definition: PFMETMonitor.h:50
T getParameter(std::string const &) const
TProfile * profileRMS_delta_set_VS_set_
Definition: PFMETMonitor.h:55
TProfile * profile_delta_set_VS_set_
Definition: PFMETMonitor.h:51
TH2F * delta_ex_VS_set_
Definition: PFMETMonitor.h:46
TH2F * delta_set_VS_set_
Definition: PFMETMonitor.h:47
void setup()
book histograms
void setup()
book histograms
bool histogramBooked_
Definition: PFMETMonitor.h:63
TH1F * delta_ex_
Definition: PFMETMonitor.h:45
TProfile * profileRMS_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:56
TProfile * profileRMS_delta_ex_VS_set_
Definition: PFMETMonitor.h:54
TH2F * book2D(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.
Definition: Benchmark.cc:43
CandidateBenchmark candBench_
Definition: PFMETMonitor.h:59
TH2F * delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:48
bool createMETSpecificHistos_
Definition: PFMETMonitor.h:62
TProfile * profile_delta_set_Over_set_VS_set_
Definition: PFMETMonitor.h:52
MatchCandidateBenchmark matchCandBench_
Definition: PFMETMonitor.h:60
TProfile * bookProfile(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.
Definition: Benchmark.cc:95
TH1F * book1D(const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either with DQM or plain root.
Definition: Benchmark.cc:25

Member Data Documentation

CandidateBenchmark PFMETMonitor::candBench_
protected

Definition at line 59 of file PFMETMonitor.h.

Referenced by fillOne(), setDirectory(), setParameters(), and setup().

bool PFMETMonitor::createMETSpecificHistos_
protected

Definition at line 62 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), setParameters(), and setup().

TH1F* PFMETMonitor::delta_ex_
protected

Definition at line 45 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TH2F* PFMETMonitor::delta_ex_VS_set_
protected

Definition at line 46 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TH2F* PFMETMonitor::delta_set_Over_set_VS_set_
protected

Definition at line 48 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TH2F* PFMETMonitor::delta_set_VS_set_
protected

Definition at line 47 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

bool PFMETMonitor::histogramBooked_
protected

Definition at line 63 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

MatchCandidateBenchmark PFMETMonitor::matchCandBench_
protected

Definition at line 60 of file PFMETMonitor.h.

Referenced by fillOne(), setDirectory(), setParameters(), and setup().

TProfile* PFMETMonitor::profile_delta_ex_VS_set_
protected

Definition at line 50 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TProfile* PFMETMonitor::profile_delta_set_Over_set_VS_set_
protected

Definition at line 52 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TProfile* PFMETMonitor::profile_delta_set_VS_set_
protected

Definition at line 51 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TProfile* PFMETMonitor::profileRMS_delta_ex_VS_set_
protected

Definition at line 54 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TProfile* PFMETMonitor::profileRMS_delta_set_Over_set_VS_set_
protected

Definition at line 56 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TProfile* PFMETMonitor::profileRMS_delta_set_VS_set_
protected

Definition at line 55 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TH1F* PFMETMonitor::px_
protected

Definition at line 43 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().

TH1F* PFMETMonitor::sumEt_
protected

Definition at line 44 of file PFMETMonitor.h.

Referenced by fillOne(), PFMETMonitor(), and setup().