CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFJetMonitor.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 PFJetMonitor::PFJetMonitor( 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 
29 
30  deltaR_ = 0;
31 
32  createPFractionHistos_ = false;
33  histogramBooked_ = false;
34 }
35 //
36 // -- Destructor
37 //
39 
40 //
41 // -- Set Parameters accessing them from ParameterSet
42 //
44 
45  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
46  onlyTwoJets_ = parameterSet.getParameter<bool>( "onlyTwoJets" );
47  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
48  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
49  createPFractionHistos_ = parameterSet.getParameter<bool>( "CreatePFractionHistos" );
50 
51  setRange( parameterSet.getParameter<double>("ptMin"),
52  parameterSet.getParameter<double>("ptMax"),
53  parameterSet.getParameter<double>("etaMin"),
54  parameterSet.getParameter<double>("etaMax"),
55  parameterSet.getParameter<double>("phiMin"),
56  parameterSet.getParameter<double>("phiMax") );
57 
60 }
61 //
62 // -- Set Parameters
63 //
64 void PFJetMonitor::setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
65  float ptmin, float ptmax, float etamin, float etamax,
66  float phimin, float phimax, bool fracHistoFlag) {
67  dRMax_ = dRMax;
68  matchCharge_ = matchCharge;
69  mode_ = mode;
70  createPFractionHistos_ = fracHistoFlag;
71 
72  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
73 
76 }
77 
78 
79 void PFJetMonitor::setParameters(float dRMax, bool onlyTwoJets, bool matchCharge, Benchmark::Mode mode,
80  float ptmin, float ptmax, float etamin, float etamax,
81  float phimin, float phimax, bool fracHistoFlag) {
82  dRMax_ = dRMax;
84  matchCharge_ = matchCharge;
85  mode_ = mode;
86  createPFractionHistos_ = fracHistoFlag;
87 
88  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
89 
92 }
93 //
94 // -- Create histograms accessing parameters from ParameterSet
95 //
97  candBench_.setup(parameterSet);
98  matchCandBench_.setup(parameterSet);
99 
100  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
101  if ( dR.getParameter<bool>("switchOn") )
102  deltaR_ = book1D("deltaR_", "#DeltaR;#DeltaR",
103  dR.getParameter<int32_t>("nBin"),
104  dR.getParameter<double>("xMin"),
105  dR.getParameter<double>("xMax"));
106 
108  delta_frac_VS_frac_muon_ = book2D("delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
109  delta_frac_VS_frac_photon_ = book2D("delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
110  delta_frac_VS_frac_electron_ = book2D("delta_frac_VS_frac_electron_", "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
111  delta_frac_VS_frac_charged_hadron_ = book2D("delta_frac_VS_frac_charged_hadron_", "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
112  delta_frac_VS_frac_neutral_hadron_ = book2D("delta_frac_VS_frac_neutral_hadron_", "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
113 
114  histogramBooked_ = true;
115  }
116 }
117 //
118 // -- Create histograms using local parameters
119 //
121  candBench_.setup();
123 
125  delta_frac_VS_frac_muon_ = book2D("delta_frac_VS_frac_muon_", "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
126  delta_frac_VS_frac_photon_ = book2D("delta_frac_VS_frac_photon_", "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
127  delta_frac_VS_frac_electron_ = book2D("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_ = book2D("delta_frac_VS_frac_charged_hadron_", "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
129  delta_frac_VS_frac_neutral_hadron_ = book2D("delta_frac_VS_frac_neutral_hadron_", "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
130 
131  histogramBooked_ = true;
132  }
133 
134 }
135 //
136 // -- Set directory to book histograms using ROOT
137 //
138 
139 void PFJetMonitor::setDirectory(TDirectory* dir) {
141 
144 }
145 //
146 // -- fill histograms for a given Jet pair
147 //
149  const reco::Jet& matchedJet) {
150  std::cout <<"\nfillone Jet histos" <<std::endl;
151 
152  const reco::PFJet* pfJet = dynamic_cast<const reco::PFJet*>(&jet);
153  const reco::PFJet* pfMatchedJet = dynamic_cast<const reco::PFJet*>(&matchedJet);
154  if (pfJet && pfMatchedJet && createPFractionHistos_) {
155  float del_frac_muon = -99.9;
156  float del_frac_elec = -99.9;
157  float del_frac_phot = -99.9;
158  float del_frac_ch_had = -99.9;
159  float del_frac_neu_had = -99.9;
160 
161  int mult_muon = pfMatchedJet->muonMultiplicity();
162  int mult_elec = pfMatchedJet->electronMultiplicity();
163  int mult_phot = pfMatchedJet->photonMultiplicity();
164  int mult_ch_had = pfMatchedJet->chargedHadronMultiplicity();
165  int mult_neu_had = pfMatchedJet->neutralHadronMultiplicity();
166 
167  if (mult_muon > 0) del_frac_muon = (pfJet->muonMultiplicity() - mult_muon)*1.0/mult_muon;
168  if (mult_elec > 0) del_frac_elec = (pfJet->electronMultiplicity() - mult_elec)*1.0/mult_elec;
169  if (mult_phot > 0) del_frac_phot = (pfJet->photonMultiplicity() - mult_phot)*1.0/mult_phot;
170  if (mult_ch_had > 0) del_frac_ch_had = (pfJet->chargedHadronMultiplicity() - mult_ch_had)*1.0/mult_ch_had;
171  if (mult_neu_had > 0) del_frac_neu_had = (pfJet->neutralHadronMultiplicity() - mult_neu_had)*1.0/mult_neu_had;
172 
173  delta_frac_VS_frac_muon_->Fill(mult_muon, del_frac_muon);
174  delta_frac_VS_frac_electron_->Fill(mult_elec, del_frac_elec);
175  delta_frac_VS_frac_photon_->Fill(mult_phot, del_frac_phot);
176  delta_frac_VS_frac_charged_hadron_->Fill(mult_ch_had, del_frac_ch_had);
177  delta_frac_VS_frac_neutral_hadron_->Fill(mult_neu_had, del_frac_neu_had);
178  }
179 }
180 
T getParameter(std::string const &) const
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:127
bool histogramBooked_
Definition: PFJetMonitor.h:75
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:61
bool onlyTwoJets
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
void setup()
book histograms
Base class for all types of Jets.
Definition: Jet.h:20
void setup()
book histograms
abstract base class
Definition: Benchmark.h:21
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:19
bool createPFractionHistos_
Definition: PFJetMonitor.h:74
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:64
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:62
Jets made from PFObjects.
Definition: PFJet.h:21
virtual ~PFJetMonitor()
Definition: PFJetMonitor.cc:38
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:68
TH2F * book2D(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.
Definition: Benchmark.cc:43
void setup()
book histograms
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:64
bool onlyTwoJets_
Definition: PFJetMonitor.h:72
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:125
void setParameters(Mode mode)
Definition: Benchmark.h:50
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.cc:15
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:65
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:66
TH1F * deltaR_
Definition: PFJetMonitor.h:70
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:123
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:131
double ptmin
Definition: HydjetWrapper.h:85
Mode mode_
Definition: Benchmark.h:96
tuple cout
Definition: gather_cfg.py:121
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:67
dbl *** dir
Definition: mlp_gen.cc:35
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:129
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
bool matchCharge_
Definition: PFJetMonitor.h:73