CMS 3D CMS Logo

PFJetMonitor.cc
Go to the documentation of this file.
6 
7 #include <TFile.h>
8 #include <TH1.h>
9 #include <TH2.h>
10 #include <TROOT.h>
11 
12 //
13 // -- Constructor
14 //
16  : Benchmark(mode), candBench_(mode), matchCandBench_(mode), dRMax_(dRMax), matchCharge_(matchCharge) {
17  setRange(0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
18 
19  delta_frac_VS_frac_muon_ = nullptr;
24 
25  deltaR_ = nullptr;
26 
27  createPFractionHistos_ = false;
28  histogramBooked_ = false;
29 }
30 
31 //
32 // -- Destructor
33 //
35 
36 //
37 // -- Set Parameters accessing them from ParameterSet
38 //
40  dRMax_ = parameterSet.getParameter<double>("deltaRMax");
41  onlyTwoJets_ = parameterSet.getParameter<bool>("onlyTwoJets");
42  matchCharge_ = parameterSet.getParameter<bool>("matchCharge");
43  mode_ = (Benchmark::Mode)parameterSet.getParameter<int>("mode");
44  createPFractionHistos_ = parameterSet.getParameter<bool>("CreatePFractionHistos");
45 
46  setRange(parameterSet.getParameter<double>("ptMin"),
47  parameterSet.getParameter<double>("ptMax"),
48  parameterSet.getParameter<double>("etaMin"),
49  parameterSet.getParameter<double>("etaMax"),
50  parameterSet.getParameter<double>("phiMin"),
51  parameterSet.getParameter<double>("phiMax"));
52 
55 }
56 
57 //
58 // -- Set Parameters
59 //
61  bool matchCharge,
63  float ptmin,
64  float ptmax,
65  float etamin,
66  float etamax,
67  float phimin,
68  float phimax,
69  bool fracHistoFlag) {
70  dRMax_ = dRMax;
71  matchCharge_ = matchCharge;
72  mode_ = mode;
73  createPFractionHistos_ = fracHistoFlag;
74 
75  setRange(ptmin, ptmax, etamin, etamax, phimin, phimax);
76 
79 }
80 
82  bool onlyTwoJets,
83  bool matchCharge,
85  float ptmin,
86  float ptmax,
87  float etamin,
88  float etamax,
89  float phimin,
90  float phimax,
91  bool fracHistoFlag) {
92  dRMax_ = dRMax;
94  matchCharge_ = matchCharge;
95  mode_ = mode;
96  createPFractionHistos_ = fracHistoFlag;
97 
98  setRange(ptmin, ptmax, etamin, etamax, phimin, phimax);
99 
102 }
103 
104 //
105 // -- Create histograms accessing parameters from ParameterSet
106 //
108  candBench_.setup(b, parameterSet);
109  matchCandBench_.setup(b, parameterSet);
110 
111  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
112  if (dR.getParameter<bool>("switchOn")) {
113  deltaR_ = book1D(b,
114  "deltaR_",
115  "#DeltaR;#DeltaR",
116  dR.getParameter<int32_t>("nBin"),
117  dR.getParameter<double>("xMin"),
118  dR.getParameter<double>("xMax"));
119  }
122  book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
124  book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
126  b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
128  "delta_frac_VS_frac_charged_hadron_",
129  "#DeltaFraction_Vs_Fraction(charged hadron)",
130  100,
131  0.0,
132  1.0,
133  100,
134  -1.0,
135  1.0);
137  "delta_frac_VS_frac_neutral_hadron_",
138  "#DeltaFraction_Vs_Fraction(neutral hadron)",
139  100,
140  0.0,
141  1.0,
142  100,
143  -1.0,
144  1.0);
145 
146  histogramBooked_ = true;
147  }
148 }
149 
150 //
151 // -- Create histograms using local parameters
152 //
154  candBench_.setup(b);
156 
159  book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
161  book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
163  b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
165  "delta_frac_VS_frac_charged_hadron_",
166  "#DeltaFraction_Vs_Fraction(charged hadron)",
167  100,
168  0.0,
169  1.0,
170  100,
171  -1.0,
172  1.0);
174  "delta_frac_VS_frac_neutral_hadron_",
175  "#DeltaFraction_Vs_Fraction(neutral hadron)",
176  100,
177  0.0,
178  1.0,
179  100,
180  -1.0,
181  1.0);
182 
183  histogramBooked_ = true;
184  }
185 }
186 
187 //
188 // -- Set directory to book histograms using ROOT
189 //
190 
191 void PFJetMonitor::setDirectory(TDirectory *dir) {
193 
196 }
197 
198 //
199 // -- fill histograms for a given Jet pair
200 //
201 void PFJetMonitor::fillOne(const reco::Jet &jet, const reco::Jet &matchedJet) {
202  const reco::PFJet *pfJet = dynamic_cast<const reco::PFJet *>(&jet);
203  const reco::PFJet *pfMatchedJet = dynamic_cast<const reco::PFJet *>(&matchedJet);
204  if (pfJet && pfMatchedJet && createPFractionHistos_) {
205  float del_frac_muon = -99.9;
206  float del_frac_elec = -99.9;
207  float del_frac_phot = -99.9;
208  float del_frac_ch_had = -99.9;
209  float del_frac_neu_had = -99.9;
210 
211  int mult_muon = pfMatchedJet->muonMultiplicity();
212  int mult_elec = pfMatchedJet->electronMultiplicity();
213  int mult_phot = pfMatchedJet->photonMultiplicity();
214  int mult_ch_had = pfMatchedJet->chargedHadronMultiplicity();
215  int mult_neu_had = pfMatchedJet->neutralHadronMultiplicity();
216 
217  if (mult_muon > 0)
218  del_frac_muon = (pfJet->muonMultiplicity() - mult_muon) * 1.0 / mult_muon;
219  if (mult_elec > 0)
220  del_frac_elec = (pfJet->electronMultiplicity() - mult_elec) * 1.0 / mult_elec;
221  if (mult_phot > 0)
222  del_frac_phot = (pfJet->photonMultiplicity() - mult_phot) * 1.0 / mult_phot;
223  if (mult_ch_had > 0)
224  del_frac_ch_had = (pfJet->chargedHadronMultiplicity() - mult_ch_had) * 1.0 / mult_ch_had;
225  if (mult_neu_had > 0)
226  del_frac_neu_had = (pfJet->neutralHadronMultiplicity() - mult_neu_had) * 1.0 / mult_neu_had;
227 
228  delta_frac_VS_frac_muon_->Fill(mult_muon, del_frac_muon);
229  delta_frac_VS_frac_electron_->Fill(mult_elec, del_frac_elec);
230  delta_frac_VS_frac_photon_->Fill(mult_phot, del_frac_phot);
231  delta_frac_VS_frac_charged_hadron_->Fill(mult_ch_had, del_frac_ch_had);
232  delta_frac_VS_frac_neutral_hadron_->Fill(mult_neu_had, del_frac_neu_had);
233  }
234 }
void setDirectory(TDirectory *dir) override
set directory (to use in ROOT)
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 through IBooker or plain root
Definition: Benchmark.cc:16
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:131
bool histogramBooked_
Definition: PFJetMonitor.h:84
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:70
bool onlyTwoJets
Base class for all types of Jets.
Definition: Jet.h:20
abstract base class
Definition: Benchmark.h:21
void setup(DQMStore::IBooker &b)
book histograms
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:14
bool createPFractionHistos_
Definition: PFJetMonitor.h:83
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:73
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:71
Jets made from PFObjects.
Definition: PFJet.h:21
void setup(DQMStore::IBooker &b)
book histograms
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:77
void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool fracHistoFlag=true)
set the parameters locally
Definition: PFJetMonitor.cc:60
bool onlyTwoJets_
Definition: PFJetMonitor.h:81
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:129
void setParameters(Mode mode)
Definition: Benchmark.h:39
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.cc:15
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 through IBooker or plain root
Definition: Benchmark.cc:23
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:74
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:75
TH1F * deltaR_
Definition: PFJetMonitor.h:79
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
double b
Definition: hdecay.h:120
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:127
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:135
double ptmin
Definition: HydjetWrapper.h:90
Mode mode_
Definition: Benchmark.h:124
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:76
void setup(DQMStore::IBooker &b)
book histograms
dbl *** dir
Definition: mlp_gen.cc:35
~PFJetMonitor() override
Definition: PFJetMonitor.cc:34
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:133
bool matchCharge_
Definition: PFJetMonitor.h:82