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  createPFractionHistos_ = false;
31  histogramBooked_ = false;
32 }
33 //
34 // -- Destructor
35 //
37 
38 //
39 // -- Set Parameters accessing them from ParameterSet
40 //
41 void PFJetMonitor::setParameters( const edm::ParameterSet & parameterSet) {
42 
43  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
44  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
45  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
46  createPFractionHistos_ = parameterSet.getParameter<bool>( "CreatePFractionHistos" );
47 
48 
49  setRange( parameterSet.getParameter<double>("ptMin"),
50  parameterSet.getParameter<double>("ptMax"),
51  parameterSet.getParameter<double>("etaMin"),
52  parameterSet.getParameter<double>("etaMax"),
53  parameterSet.getParameter<double>("phiMin"),
54  parameterSet.getParameter<double>("phiMax") );
55 
58 }
59 //
60 // -- Set Parameters
61 //
62 void PFJetMonitor::setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
63  float ptmin, float ptmax, float etamin, float etamax,
64  float phimin, float phimax, bool fracHistoFlag) {
65  dRMax_ = dRMax;
66  matchCharge_ = matchCharge;
67  mode_ = mode;
68  createPFractionHistos_ = fracHistoFlag;
69 
70  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
71 
74 
75 }
76 //
77 // -- Create histograms accessing parameters from ParameterSet
78 //
79 void PFJetMonitor::setup(const edm::ParameterSet & parameterSet) {
80  candBench_.setup(parameterSet);
81  matchCandBench_.setup(parameterSet);
82 
84  delta_frac_VS_frac_muon_ = book2D("delta_frac_VS_frac_muon_",
85  "#DeltaFraction_Vs_Fraction(muon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
86  delta_frac_VS_frac_photon_ = book2D("delta_frac_VS_frac_photon_",
87  "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
88  delta_frac_VS_frac_electron_ = book2D("delta_frac_VS_frac_electron_",
89  "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
90  delta_frac_VS_frac_charged_hadron_ = book2D("delta_frac_VS_frac_charged_hadron_",
91  "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
92  delta_frac_VS_frac_neutral_hadron_ = book2D("delta_frac_VS_frac_neutral_hadron_",
93  "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
94 
95  histogramBooked_ = true;
96  }
97 
98 }
99 //
100 // -- Create histograms using local parameters
101 //
103  candBench_.setup();
105 
107  delta_frac_VS_frac_muon_ = book2D("delta_frac_VS_frac_muon_",
108  "#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_",
110  "#DeltaFraction_Vs_Fraction(photon)", 100, 0.0, 1.0, 100, -1.0, 1.0);
111  delta_frac_VS_frac_electron_ = book2D("delta_frac_VS_frac_electron_",
112  "#DeltaFraction_Vs_Fraction(electron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
113  delta_frac_VS_frac_charged_hadron_ = book2D("delta_frac_VS_frac_charged_hadron_",
114  "#DeltaFraction_Vs_Fraction(charged hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
115  delta_frac_VS_frac_neutral_hadron_ = book2D("delta_frac_VS_frac_neutral_hadron_",
116  "#DeltaFraction_Vs_Fraction(neutral hadron)", 100, 0.0, 1.0, 100, -1.0, 1.0);
117 
118  histogramBooked_ = true;
119  }
120 
121 }
122 //
123 // -- Set directory to book histograms using ROOT
124 //
125 
126 void PFJetMonitor::setDirectory(TDirectory* dir) {
128 
131 }
132 //
133 // -- fill histograms for a given Jet pair
134 //
136  const reco::Jet& matchedJet) {
137  const reco::PFJet* pfJet = dynamic_cast<const reco::PFJet*>(&jet);
138  const reco::PFJet* pfMatchedJet = dynamic_cast<const reco::PFJet*>(&matchedJet);
139  if (pfJet && pfMatchedJet && createPFractionHistos_) {
140  float frac_muon = -99.9;
141  float frac_elec = -99.9;
142  float frac_phot = -99.9;
143  float frac_ch_had = -99.9;
144  float frac_neu_had = -99.9;
145 
146  if (pfMatchedJet->muonMultiplicity() > 0) frac_muon = (pfJet->muonMultiplicity() - pfMatchedJet->muonMultiplicity())*1.0/pfMatchedJet->muonMultiplicity();
147  if (pfMatchedJet->chargedHadronMultiplicity() > 0) frac_ch_had = (pfJet->chargedHadronMultiplicity() - pfMatchedJet->chargedHadronMultiplicity())*1.0/pfMatchedJet->chargedHadronMultiplicity();
148  if (pfMatchedJet->neutralHadronMultiplicity() > 0) frac_neu_had = (pfJet->neutralHadronMultiplicity() - pfMatchedJet->neutralHadronMultiplicity())*1.0/pfMatchedJet->neutralHadronMultiplicity();
149  if (pfMatchedJet->photonMultiplicity() > 0) frac_phot = (pfJet->photonMultiplicity() - pfMatchedJet->photonMultiplicity())*1.0/pfMatchedJet->photonMultiplicity();
150  if (pfMatchedJet->electronMultiplicity() > 0) frac_elec = (pfJet->electronMultiplicity() - pfMatchedJet->electronMultiplicity())*1.0/pfMatchedJet->electronMultiplicity();
151 
152  delta_frac_VS_frac_muon_->Fill(frac_muon);
153  delta_frac_VS_frac_electron_->Fill(frac_elec);
154  delta_frac_VS_frac_photon_->Fill(frac_phot);
155  delta_frac_VS_frac_charged_hadron_->Fill(frac_ch_had);
156  delta_frac_VS_frac_neutral_hadron_->Fill(frac_neu_had);
157  }
158 }
159 
T getParameter(std::string const &) const
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:124
bool histogramBooked_
Definition: PFJetMonitor.h:61
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:49
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
void setup()
book histograms
Base class for all types of Jets.
Definition: Jet.h:21
void setup()
book histograms
abstract base class
Definition: Benchmark.h:20
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:19
bool createPFractionHistos_
Definition: PFJetMonitor.h:60
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:52
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:50
Jets made from PFObjects.
Definition: PFJet.h:22
virtual ~PFJetMonitor()
Definition: PFJetMonitor.cc:36
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:56
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:62
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:122
void setParameters(Mode mode)
Definition: Benchmark.h:49
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:53
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:54
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:51
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:120
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:128
double ptmin
Definition: HydjetWrapper.h:86
Mode mode_
Definition: Benchmark.h:86
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:55
dbl *** dir
Definition: mlp_gen.cc:35
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:126
bool matchCharge_
Definition: PFJetMonitor.h:59