CMS 3D CMS Logo

PFJetMonitor.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 
30 
31  deltaR_ = 0;
32 
33  createPFractionHistos_ = false;
34  histogramBooked_ = false;
35 }
36 
37 
38 //
39 // -- Destructor
40 //
42 
43 
44 //
45 // -- Set Parameters accessing them from ParameterSet
46 //
48 
49  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
50  onlyTwoJets_ = parameterSet.getParameter<bool>( "onlyTwoJets" );
51  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
52  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
53  createPFractionHistos_ = parameterSet.getParameter<bool>( "CreatePFractionHistos" );
54 
55  setRange( parameterSet.getParameter<double>("ptMin"),
56  parameterSet.getParameter<double>("ptMax"),
57  parameterSet.getParameter<double>("etaMin"),
58  parameterSet.getParameter<double>("etaMax"),
59  parameterSet.getParameter<double>("phiMin"),
60  parameterSet.getParameter<double>("phiMax") );
61 
64 }
65 
66 
67 //
68 // -- Set Parameters
69 //
70 void PFJetMonitor::setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
71  float ptmin, float ptmax, float etamin, float etamax,
72  float phimin, float phimax, bool fracHistoFlag) {
73  dRMax_ = dRMax;
74  matchCharge_ = matchCharge;
75  mode_ = mode;
76  createPFractionHistos_ = fracHistoFlag;
77 
78  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
79 
82 }
83 
84 
86  float ptmin, float ptmax, float etamin, float etamax,
87  float phimin, float phimax, bool fracHistoFlag) {
88  dRMax_ = dRMax;
90  matchCharge_ = matchCharge;
91  mode_ = mode;
92  createPFractionHistos_ = fracHistoFlag;
93 
94  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
95 
98 }
99 
100 
101 //
102 // -- Create histograms accessing parameters from ParameterSet
103 //
105  candBench_.setup(b, parameterSet);
106  matchCandBench_.setup(b, parameterSet);
107 
108  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
109  if ( dR.getParameter<bool>("switchOn") ) {
110  deltaR_ = book1D(b, "deltaR_", "#DeltaR;#DeltaR",
111  dR.getParameter<int32_t>("nBin"),
112  dR.getParameter<double>("xMin"),
113  dR.getParameter<double>("xMax"));
114  }
116  delta_frac_VS_frac_muon_ = book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
117  delta_frac_VS_frac_photon_ = book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
118  delta_frac_VS_frac_electron_ = book2D(b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
119  delta_frac_VS_frac_charged_hadron_ = book2D(b, "delta_frac_VS_frac_charged_hadron_", "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
120  delta_frac_VS_frac_neutral_hadron_ = book2D(b, "delta_frac_VS_frac_neutral_hadron_", "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
121 
122  histogramBooked_ = true;
123  }
124 }
125 
126 
127 //
128 // -- Create histograms using local parameters
129 //
131  candBench_.setup(b);
133 
135  delta_frac_VS_frac_muon_ = book2D(b, "delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
136  delta_frac_VS_frac_photon_ = book2D(b, "delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
137  delta_frac_VS_frac_electron_ = book2D(b, "delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
138  delta_frac_VS_frac_charged_hadron_ = book2D(b, "delta_frac_VS_frac_charged_hadron_", "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
139  delta_frac_VS_frac_neutral_hadron_ = book2D(b, "delta_frac_VS_frac_neutral_hadron_", "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
140 
141  histogramBooked_ = true;
142  }
143 
144 }
145 
146 
147 //
148 // -- Set directory to book histograms using ROOT
149 //
150 
151 void PFJetMonitor::setDirectory(TDirectory* dir) {
153 
156 }
157 
158 
159 //
160 // -- fill histograms for a given Jet pair
161 //
163  const reco::Jet& matchedJet) {
164 
165  const reco::PFJet* pfJet = dynamic_cast<const reco::PFJet*>(&jet);
166  const reco::PFJet* pfMatchedJet = dynamic_cast<const reco::PFJet*>(&matchedJet);
167  if (pfJet && pfMatchedJet && createPFractionHistos_) {
168  float del_frac_muon = -99.9;
169  float del_frac_elec = -99.9;
170  float del_frac_phot = -99.9;
171  float del_frac_ch_had = -99.9;
172  float del_frac_neu_had = -99.9;
173 
174  int mult_muon = pfMatchedJet->muonMultiplicity();
175  int mult_elec = pfMatchedJet->electronMultiplicity();
176  int mult_phot = pfMatchedJet->photonMultiplicity();
177  int mult_ch_had = pfMatchedJet->chargedHadronMultiplicity();
178  int mult_neu_had = pfMatchedJet->neutralHadronMultiplicity();
179 
180  if (mult_muon > 0) del_frac_muon = (pfJet->muonMultiplicity() - mult_muon)*1.0/mult_muon;
181  if (mult_elec > 0) del_frac_elec = (pfJet->electronMultiplicity() - mult_elec)*1.0/mult_elec;
182  if (mult_phot > 0) del_frac_phot = (pfJet->photonMultiplicity() - mult_phot)*1.0/mult_phot;
183  if (mult_ch_had > 0) del_frac_ch_had = (pfJet->chargedHadronMultiplicity() - mult_ch_had)*1.0/mult_ch_had;
184  if (mult_neu_had > 0) del_frac_neu_had = (pfJet->neutralHadronMultiplicity() - mult_neu_had)*1.0/mult_neu_had;
185 
186  delta_frac_VS_frac_muon_->Fill(mult_muon, del_frac_muon);
187  delta_frac_VS_frac_electron_->Fill(mult_elec, del_frac_elec);
188  delta_frac_VS_frac_photon_->Fill(mult_phot, del_frac_phot);
189  delta_frac_VS_frac_charged_hadron_->Fill(mult_ch_had, del_frac_ch_had);
190  delta_frac_VS_frac_neutral_hadron_->Fill(mult_neu_had, del_frac_neu_had);
191  }
192 }
193 
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
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:131
bool histogramBooked_
Definition: PFJetMonitor.h:72
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:58
bool onlyTwoJets
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
Base class for all types of Jets.
Definition: Jet.h:20
abstract base class
Definition: Benchmark.h:22
void setup(DQMStore::IBooker &b)
book histograms
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:18
bool createPFractionHistos_
Definition: PFJetMonitor.h:71
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:61
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:59
Jets made from PFObjects.
Definition: PFJet.h:21
void setup(DQMStore::IBooker &b)
book histograms
virtual ~PFJetMonitor()
Definition: PFJetMonitor.cc:41
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:65
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:70
bool onlyTwoJets_
Definition: PFJetMonitor.h:69
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:129
void setParameters(Mode mode)
Definition: Benchmark.h:49
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.cc:16
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:29
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:62
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:63
TH1F * deltaR_
Definition: PFJetMonitor.h:67
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:51
double b
Definition: hdecay.h:120
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:127
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:135
Mode mode_
Definition: Benchmark.h:105
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:64
void setup(DQMStore::IBooker &b)
book histograms
dbl *** dir
Definition: mlp_gen.cc:35
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:133
bool matchCharge_
Definition: PFJetMonitor.h:70