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");
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  parameterSet.getParameter<double>("ptMax"),
56  parameterSet.getParameter<double>("etaMin"),
57  parameterSet.getParameter<double>("etaMax"),
58  parameterSet.getParameter<double>("phiMin"),
59  parameterSet.getParameter<double>("phiMax"));
60 
63  parameterSet.getParameter<double>("ptMax"),
64  parameterSet.getParameter<double>("etaMin"),
65  parameterSet.getParameter<double>("etaMax"),
66  parameterSet.getParameter<double>("phiMin"),
67  parameterSet.getParameter<double>("phiMax"));
68 }
69 
70 //
71 // -- Create histograms accessing parameters from ParameterSet
72 //
76 
78  if (dR.getParameter<bool>("switchOn")) {
79  deltaR_ = book1D(b,
80  "deltaR_",
81  "#DeltaR;#DeltaR",
82  dR.getParameter<int32_t>("nBin"),
83  dR.getParameter<double>("xMin"),
84  dR.getParameter<double>("xMax"));
85  }
88  book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
90  book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
92  b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
94  "delta_frac_VS_frac_charged_hadron_",
95  "#DeltaFraction_Vs_Fraction(charged hadron)",
96  100,
97  0.0,
98  1.0,
99  100,
100  -1.0,
101  1.0);
103  "delta_frac_VS_frac_neutral_hadron_",
104  "#DeltaFraction_Vs_Fraction(neutral hadron)",
105  100,
106  0.0,
107  1.0,
108  100,
109  -1.0,
110  1.0);
111 
112  histogramBooked_ = true;
113  }
114 }
115 
116 //
117 // -- Create histograms using local parameters
118 //
120  candBench_.setup(b);
122 
125  book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
127  book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
129  b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
131  "delta_frac_VS_frac_charged_hadron_",
132  "#DeltaFraction_Vs_Fraction(charged hadron)",
133  100,
134  0.0,
135  1.0,
136  100,
137  -1.0,
138  1.0);
140  "delta_frac_VS_frac_neutral_hadron_",
141  "#DeltaFraction_Vs_Fraction(neutral hadron)",
142  100,
143  0.0,
144  1.0,
145  100,
146  -1.0,
147  1.0);
148 
149  histogramBooked_ = true;
150  }
151 }
152 
153 //
154 // -- Set directory to book histograms using ROOT
155 //
156 
157 void PFJetMonitor::setDirectory(TDirectory *dir) {
159 
162 }
163 
164 //
165 // -- fill histograms for a given Jet pair
166 //
167 void PFJetMonitor::fillOne(const reco::Jet &jet, const reco::Jet &matchedJet) {
168  const reco::PFJet *pfJet = dynamic_cast<const reco::PFJet *>(&jet);
169  const reco::PFJet *pfMatchedJet = dynamic_cast<const reco::PFJet *>(&matchedJet);
170  if (pfJet && pfMatchedJet && createPFractionHistos_) {
171  float del_frac_muon = -99.9;
172  float del_frac_elec = -99.9;
173  float del_frac_phot = -99.9;
174  float del_frac_ch_had = -99.9;
175  float del_frac_neu_had = -99.9;
176 
177  int mult_muon = pfMatchedJet->muonMultiplicity();
178  int mult_elec = pfMatchedJet->electronMultiplicity();
179  int mult_phot = pfMatchedJet->photonMultiplicity();
180  int mult_ch_had = pfMatchedJet->chargedHadronMultiplicity();
181  int mult_neu_had = pfMatchedJet->neutralHadronMultiplicity();
182 
183  if (mult_muon > 0)
184  del_frac_muon = (pfJet->muonMultiplicity() - mult_muon) * 1.0 / mult_muon;
185  if (mult_elec > 0)
186  del_frac_elec = (pfJet->electronMultiplicity() - mult_elec) * 1.0 / mult_elec;
187  if (mult_phot > 0)
188  del_frac_phot = (pfJet->photonMultiplicity() - mult_phot) * 1.0 / mult_phot;
189  if (mult_ch_had > 0)
190  del_frac_ch_had = (pfJet->chargedHadronMultiplicity() - mult_ch_had) * 1.0 / mult_ch_had;
191  if (mult_neu_had > 0)
192  del_frac_neu_had = (pfJet->neutralHadronMultiplicity() - mult_neu_had) * 1.0 / mult_neu_had;
193 
194  delta_frac_VS_frac_muon_->Fill(mult_muon, del_frac_muon);
195  delta_frac_VS_frac_electron_->Fill(mult_elec, del_frac_elec);
196  delta_frac_VS_frac_photon_->Fill(mult_phot, del_frac_phot);
197  delta_frac_VS_frac_charged_hadron_->Fill(mult_ch_had, del_frac_ch_had);
198  delta_frac_VS_frac_neutral_hadron_->Fill(mult_neu_had, del_frac_neu_had);
199  }
200 }
Benchmark
abstract base class
Definition: Benchmark.h:19
PFJetMonitor::delta_frac_VS_frac_neutral_hadron_
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:52
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
reco::PFJet::neutralHadronMultiplicity
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:126
MatchCandidateBenchmark::setup
void setup(DQMStore::IBooker &b)
book histograms
Definition: MatchCandidateBenchmark.cc:28
Benchmark::book2D
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:22
Matchers.h
PFJet.h
ALCARECOPromptCalibProdSiPixelAli0T_cff.mode
mode
Definition: ALCARECOPromptCalibProdSiPixelAli0T_cff.py:96
PFJetMonitor::histogramBooked_
bool histogramBooked_
Definition: PFJetMonitor.h:59
PFJetMonitor::setup
void setup(DQMStore::IBooker &b)
book histograms
Definition: PFJetMonitor.cc:119
Jet.h
Benchmark::setParameters
void setParameters(Mode mode)
Definition: Benchmark.h:39
Benchmark::setRange
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
CandidateBenchmark::setup
void setup(DQMStore::IBooker &b)
book histograms
Definition: CandidateBenchmark.cc:24
PFJetMonitor::dRMax_
float dRMax_
Definition: PFJetMonitor.h:55
pfCandidateManager_cfi.matchCharge
matchCharge
Definition: pfCandidateManager_cfi.py:15
PFJetMonitor::createPFractionHistos_
bool createPFractionHistos_
Definition: PFJetMonitor.h:58
PFJetMonitor::matchCharge_
bool matchCharge_
Definition: PFJetMonitor.h:57
reco::PFJet::chargedHadronMultiplicity
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:124
PFJetMonitor.h
PFJetMonitor::delta_frac_VS_frac_muon_
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:48
PFJetMonitor::setParameters
void setParameters(const edm::ParameterSet &parameterSet)
set the parameters accessing them from ParameterSet
Definition: PFJetMonitor.cc:39
Benchmark::mode_
Mode mode_
Definition: Benchmark.h:118
b
double b
Definition: hdecay.h:118
PFJetMonitor::delta_frac_VS_frac_photon_
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:49
PFJetMonitor::setDirectory
void setDirectory(TDirectory *dir) override
set directory (to use in ROOT)
Definition: PFJetMonitor.cc:157
edm::ParameterSet
Definition: ParameterSet.h:36
Benchmark::Mode
Mode
Definition: Benchmark.h:32
reco::PFJet::electronMultiplicity
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:130
reco::PFJet::photonMultiplicity
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:128
PFJetMonitor::~PFJetMonitor
~PFJetMonitor() override
Definition: PFJetMonitor.cc:34
PFJetMonitor::candBench_
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:45
Benchmark::setDirectory
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:13
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Benchmark::book1D
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:15
reco::PFJet::muonMultiplicity
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:132
metBenchmark_cfi.dRMax
dRMax
Definition: metBenchmark_cfi.py:18
PFJetMonitor::fillOne
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
Definition: PFJetMonitor.cc:167
metsig::jet
Definition: SignAlgoResolutions.h:47
PFJetMonitor::deltaR_
TH1F * deltaR_
Definition: PFJetMonitor.h:54
edm::parameterSet
ParameterSet const & parameterSet(Provenance const &provenance, ProcessHistory const &history)
Definition: Provenance.cc:11
reco::PFJet
Jets made from PFObjects.
Definition: PFJet.h:20
PFJetMonitor::matchCandBench_
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:46
dqm::implementation::IBooker
Definition: DQMStore.h:43
PFJetMonitor::onlyTwoJets_
bool onlyTwoJets_
Definition: PFJetMonitor.h:56
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
ParameterSet.h
PFJetMonitor::delta_frac_VS_frac_charged_hadron_
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:51
PFJetMonitor::PFJetMonitor
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.cc:15
PFJetMonitor::delta_frac_VS_frac_electron_
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:50
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23