CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateMonitor.cc
Go to the documentation of this file.
6 
7 #include <TROOT.h>
8 #include <TFile.h>
9 #include <TH1.h>
10 #include <TH2.h>
11 
12 
13 //
14 // -- Constructor
15 //
17  Benchmark(mode),
18  candBench_(mode),
19  matchCandBench_(mode),
20  dRMax_(dRMax),
21  matchCharge_(matchCharge) {
22 
23  setRange( 0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
24 
25  pt_gen_ = 0;
26  eta_gen_ = 0;
27  phi_gen_ = 0;
28 
29  pt_ref_ = 0;
30  eta_ref_ = 0;
31  phi_ref_ = 0;
32 
33  deltaR_ = 0;
34 
35  createReferenceHistos_ = false;
36  histogramBooked_ = false;
37 }
38 
39 
40 //
41 // -- Destructor
42 //
44 
45 
46 //
47 // -- Set Parameters accessing them from ParameterSet
48 //
50 
51  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
52  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
53  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
54  createReferenceHistos_ = parameterSet.getParameter<bool>( "CreateReferenceHistos" );
55  createEfficiencyHistos_ = parameterSet.getParameter<bool>( "CreateEfficiencyHistos" );
56 
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 
66 }
67 
68 
69 //
70 // -- Set Parameters
71 //
73  float ptmin, float ptmax, float etamin, float etamax,
74  float phimin, float phimax, bool refHistoFlag) {
75  dRMax_ = dRMax;
76  matchCharge_ = matchCharge;
77  mode_ = mode;
78  createReferenceHistos_ = refHistoFlag;
79 
80  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
81 
84 }
85 
86 
87 //
88 // -- Create histograms accessing parameters from ParameterSet
89 //
91  candBench_.setup(b, parameterSet);
92  matchCandBench_.setup(b, parameterSet);
93 
95  edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
96  edm::ParameterSet etaPS = parameterSet.getParameter<edm::ParameterSet>("EtaHistoParameter");
97  edm::ParameterSet phiPS = parameterSet.getParameter<edm::ParameterSet>("PhiHistoParameter");
98 
99  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
100 
101  if (ptPS.getParameter<bool>("switchOn")) {
102  pt_ref_ = book1D(b, "pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
103  ptPS.getParameter<double>("xMin"),
104  ptPS.getParameter<double>("xMax"));
106  pt_gen_ = book1D(b, "pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
107  ptPS.getParameter<double>("xMin"),
108  ptPS.getParameter<double>("xMax") ) ;
109  }
110  }
111 
112  if (etaPS.getParameter<bool>("switchOn")) {
113  eta_ref_ = book1D(b, "eta_ref_", "#eta_ref;#eta", etaPS.getParameter<int32_t>("nBin"),
114  etaPS.getParameter<double>("xMin"),
115  etaPS.getParameter<double>("xMax"));
117  eta_gen_ = book1D(b, "eta_gen_", "#eta_gen;#eta", etaPS.getParameter<int32_t>("nBin"),
118  etaPS.getParameter<double>("xMin"),
119  etaPS.getParameter<double>("xMax") ) ;
120  }
121  }
122 
123  if (phiPS.getParameter<bool>("switchOn")) {
124  phi_ref_ = book1D(b, "phi_ref_", "#phi_ref;#phi", phiPS.getParameter<int32_t>("nBin"),
125  phiPS.getParameter<double>("xMin"),
126  phiPS.getParameter<double>("xMax"));
128  phi_gen_ = book1D(b, "phi_gen_", "#phi_gen;#phi", phiPS.getParameter<int32_t>("nBin"),
129  phiPS.getParameter<double>("xMin"),
130  phiPS.getParameter<double>("xMax") ) ;
131  }
132  }
133 
134  if ( createEfficiencyHistos_ && dR.getParameter<bool>("switchOn") ) {
135  deltaR_ = book1D(b, "deltaR_", "#DeltaR;#DeltaR",
136  dR.getParameter<int32_t>("nBin"),
137  dR.getParameter<double>("xMin"),
138  dR.getParameter<double>("xMax"));
139  }
140 
141  histogramBooked_ = true;
142  }
143 }
144 
145 
146 //
147 // -- Create histograms using local parameters
148 //
150  candBench_.setup(b);
152 
154  PhaseSpace ptPS(100,0,100);
155  PhaseSpace phiPS(360, -3.1416, 3.1416);
156  PhaseSpace etaPS(100, -5,5);
157 
158  pt_ref_ = book1D(b, "pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
160  pt_gen_ = book1D(b, "pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
161  }
162 
163  eta_ref_ = book1D(b, "eta_ref_", "#eta_ref;#eta", etaPS.n, etaPS.m, etaPS.M);
165  eta_gen_ = book1D(b, "eta_gen_", "#eta_gen;#eta", etaPS.n, etaPS.m, etaPS.M);
166  }
167 
168  phi_ref_ = book1D(b, "phi_ref_", "#phi_ref;#phi", phiPS.n, phiPS.m, phiPS.M);
170  phi_gen_ = book1D(b, "phi_gen_", "#phi_gen;#phi", phiPS.n, phiPS.m, phiPS.M);
171  }
172 
173  histogramBooked_ = true;
174  }
175 }
176 
177 
178 //
179 // -- Set directory to book histograms using ROOT
180 //
183 
186 }
187 
188 
189 //
190 // -- fill histograms for a single collection
191 //
193 
194  if (matching_done_) {
196  if (pt_ref_) pt_ref_->Fill( cand.pt() );
197  if (eta_ref_) eta_ref_->Fill( cand.eta() );
198  if (phi_ref_) phi_ref_->Fill( cand.phi() );
199  }
201  if (pt_gen_) pt_gen_->Fill( cand.pt() );
202  if (eta_gen_) eta_gen_->Fill( cand.eta() );
203  if (phi_gen_) phi_gen_->Fill( cand.phi() );
204  }
205 
206 }
T getParameter(std::string const &) const
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
virtual double pt() const =0
transverse momentum
abstract base class
Definition: Benchmark.h:22
void setup(DQMStore::IBooker &b)
book histograms
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:18
void setup(DQMStore::IBooker &b)
book histograms
void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool refHistoFlag)
set the parameters locally
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
void setParameters(Mode mode)
Definition: Benchmark.h:49
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:51
void fillOne(const reco::Candidate &cand)
double b
Definition: hdecay.h:120
double ptmin
Definition: HydjetWrapper.h:90
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:105
void setup(DQMStore::IBooker &b)
book histograms
dbl *** dir
Definition: mlp_gen.cc:35
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)