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.
5 
7 
8 #include <TROOT.h>
9 #include <TFile.h>
10 #include <TH1.h>
11 #include <TH2.h>
12 //
13 // -- Constructor
14 //
15 PFCandidateMonitor::PFCandidateMonitor( float dRMax, bool matchCharge, Benchmark::Mode mode) :
16  Benchmark(mode),
17  candBench_(mode),
18  matchCandBench_(mode),
19  dRMax_(dRMax),
20  matchCharge_(matchCharge) {
21 
22  setRange( 0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
23 
24  pt_gen_ = 0;
25  eta_gen_ = 0;
26  phi_gen_ = 0;
27 
28  pt_ref_ = 0;
29  eta_ref_ = 0;
30  phi_ref_ = 0;
31 
32  deltaR_ = 0;
33 
34  createReferenceHistos_ = false;
35  histogramBooked_ = false;
36 }
37 //
38 // -- Destructor
39 //
41 
42 //
43 // -- Set Parameters accessing them from ParameterSet
44 //
46 
47  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
48  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
49  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
50  createReferenceHistos_ = parameterSet.getParameter<bool>( "CreateReferenceHistos" );
51  createEfficiencyHistos_ = parameterSet.getParameter<bool>( "CreateEfficiencyHistos" );
52 
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 
62 }
63 //
64 // -- Set Parameters
65 //
66 void PFCandidateMonitor::setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
67  float ptmin, float ptmax, float etamin, float etamax,
68  float phimin, float phimax, bool refHistoFlag) {
69  dRMax_ = dRMax;
70  matchCharge_ = matchCharge;
71  mode_ = mode;
72  createReferenceHistos_ = refHistoFlag;
73 
74  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
75 
78 }
79 //
80 // -- Create histograms accessing parameters from ParameterSet
81 //
83  candBench_.setup(parameterSet);
84  matchCandBench_.setup(parameterSet);
85 
87  edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
88  edm::ParameterSet etaPS = parameterSet.getParameter<edm::ParameterSet>("EtaHistoParameter");
89  edm::ParameterSet phiPS = parameterSet.getParameter<edm::ParameterSet>("PhiHistoParameter");
90 
91  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
92 
93  if (ptPS.getParameter<bool>("switchOn")) {
94  pt_ref_ = book1D("pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
95  ptPS.getParameter<double>("xMin"),
96  ptPS.getParameter<double>("xMax"));
97  if (createEfficiencyHistos_) pt_gen_ = book1D( "pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
98  ptPS.getParameter<double>("xMin"),
99  ptPS.getParameter<double>("xMax") ) ;
100  }
101 
102  if (etaPS.getParameter<bool>("switchOn")) {
103  eta_ref_ = book1D("eta_ref_", "#eta_ref;#eta", etaPS.getParameter<int32_t>("nBin"),
104  etaPS.getParameter<double>("xMin"),
105  etaPS.getParameter<double>("xMax"));
106  if (createEfficiencyHistos_) eta_gen_ = book1D( "eta_gen_", "#eta_gen;#eta", etaPS.getParameter<int32_t>("nBin"),
107  etaPS.getParameter<double>("xMin"),
108  etaPS.getParameter<double>("xMax") ) ;
109  }
110  if (phiPS.getParameter<bool>("switchOn")) {
111  phi_ref_ = book1D("phi_ref_", "#phi_ref;#phi", phiPS.getParameter<int32_t>("nBin"),
112  phiPS.getParameter<double>("xMin"),
113  phiPS.getParameter<double>("xMax"));
114  if (createEfficiencyHistos_) phi_gen_ = book1D( "phi_gen_", "#phi_gen;#phi", phiPS.getParameter<int32_t>("nBin"),
115  phiPS.getParameter<double>("xMin"),
116  phiPS.getParameter<double>("xMax") ) ;
117  }
118  if ( createEfficiencyHistos_ && dR.getParameter<bool>("switchOn") )
119  deltaR_ = book1D("deltaR_", "#DeltaR;#DeltaR",
120  dR.getParameter<int32_t>("nBin"),
121  dR.getParameter<double>("xMin"),
122  dR.getParameter<double>("xMax"));
123 
124  histogramBooked_ = true;
125  }
126 }
127 //
128 // -- Create histograms using local parameters
129 //
131  candBench_.setup();
133 
135  PhaseSpace ptPS(100,0,100);
136  PhaseSpace phiPS(360, -3.1416, 3.1416);
137  PhaseSpace etaPS(100, -5,5);
138 
139  pt_ref_ = book1D("pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
140  if (createEfficiencyHistos_) pt_gen_ = book1D("pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
141 
142  eta_ref_ = book1D("eta_ref_", "#eta_ref;#eta", etaPS.n, etaPS.m, etaPS.M);
143  if (createEfficiencyHistos_) eta_gen_ = book1D("eta_gen_", "#eta_gen;#eta", etaPS.n, etaPS.m, etaPS.M);
144 
145  phi_ref_ = book1D("phi_ref_", "#phi_ref;#phi", phiPS.n, phiPS.m, phiPS.M);
146  if (createEfficiencyHistos_) phi_gen_ = book1D("phi_gen_", "#phi_gen;#phi", phiPS.n, phiPS.m, phiPS.M);
147 
148  histogramBooked_ = true;
149  }
150 }
151 //
152 // -- Set directory to book histograms using ROOT
153 //
156 
159 }
160 //
161 // -- fill histograms for a single collection
162 //
164 
165  if (matching_done_) {
167  if (pt_ref_) pt_ref_->Fill(cand.pt());
168  if (eta_ref_) eta_ref_->Fill(cand.eta() );
169  if (phi_ref_) phi_ref_->Fill(cand.phi() );
170  }
172  if (pt_gen_) pt_gen_->Fill(cand.pt());
173  if (eta_gen_) eta_gen_->Fill(cand.eta() );
174  if (phi_gen_) phi_gen_->Fill(cand.phi() );
175  }
176 
177 }
178 
T getParameter(std::string const &) const
void setup()
book histograms
virtual float eta() const =0
momentum pseudorapidity
void setup()
book histograms
void setup()
book histograms
abstract base class
Definition: Benchmark.h:21
virtual float phi() const =0
momentum azimuthal angle
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:19
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
virtual float pt() const =0
transverse momentum
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
void setParameters(Mode mode)
Definition: Benchmark.h:50
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
void fillOne(const reco::Candidate &cand)
double ptmin
Definition: HydjetWrapper.h:85
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:96
dbl *** dir
Definition: mlp_gen.cc:35
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
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
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)