CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes
PFJetMonitor Class Reference

#include <PFJetMonitor.h>

Inheritance diagram for PFJetMonitor:
Benchmark

Public Member Functions

template<class T , class C >
void fill (const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal)
 fill histograms with all particle More...
 
template<class T , class C >
void fill (const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, float &jetpT, const edm::ParameterSet &parameterSet)
 
void fillOne (const reco::Jet &jet, const reco::Jet &matchedJet)
 
 PFJetMonitor (float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
 
void setDirectory (TDirectory *dir)
 set directory (to use in ROOT) More...
 
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 More...
 
void setParameters (float dRMax, bool onlyTwoJets, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool fracHistoFlag=true)
 
void setParameters (const edm::ParameterSet &parameterSet)
 set the parameters accessing them from ParameterSet More...
 
void setup (DQMStore::IBooker &b)
 book histograms More...
 
void setup (DQMStore::IBooker &b, const edm::ParameterSet &parameterSet)
 
virtual ~PFJetMonitor ()
 
- Public Member Functions inherited from Benchmark
 Benchmark (Mode mode=DEFAULT)
 
bool isInRange (float pt, float eta, float phi) const
 
void setParameters (Mode mode)
 
void setRange (float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
 
void write ()
 write to the TFile, in plain ROOT mode. No need to call this function in DQM mode More...
 
virtual ~Benchmark ()
 

Protected Attributes

CandidateBenchmark candBench_
 
bool createPFractionHistos_
 
TH2F * delta_frac_VS_frac_charged_hadron_
 
TH2F * delta_frac_VS_frac_electron_
 
TH2F * delta_frac_VS_frac_muon_
 
TH2F * delta_frac_VS_frac_neutral_hadron_
 
TH2F * delta_frac_VS_frac_photon_
 
TH1F * deltaR_
 
float dRMax_
 
bool histogramBooked_
 
MatchCandidateBenchmark matchCandBench_
 
bool matchCharge_
 
bool onlyTwoJets_
 
- Protected Attributes inherited from Benchmark
TDirectory * dir_
 
float etaMax_
 
float etaMin_
 
Mode mode_
 
float phiMax_
 
float phiMin_
 
float ptMax_
 
float ptMin_
 

Additional Inherited Members

- Public Types inherited from Benchmark
enum  Mode { DEFAULT, DQMOFFLINE, VALIDATION }
 
- Protected Member Functions inherited from Benchmark
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 analyzer or not. More...
 
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 analyzer or not. More...
 
TH2F * book2D (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float *xbins, 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 analyzer or not. More...
 
TProfile * bookProfile (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float xmin, float xmax, float ymin, float ymax, const char *option)
 book a TProfile histogram, either with DQM or plain root depending if DQM_ has been initialized in a child analyzer or not. More...
 
TProfile * bookProfile (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float *xbins, float ymin, float ymax, const char *option)
 book a TProfile histogram, either with DQM or plain root depending if DQM_ has been initialized in a child analyzer or not. More...
 

Detailed Description

Definition at line 15 of file PFJetMonitor.h.

Constructor & Destructor Documentation

PFJetMonitor::PFJetMonitor ( float  dRMax = 0.3,
bool  matchCharge = true,
Benchmark::Mode  mode = Benchmark::DEFAULT 
)

Definition at line 16 of file PFJetMonitor.cc.

PFJetMonitor::~PFJetMonitor ( )
virtual

Definition at line 41 of file PFJetMonitor.cc.

Member Function Documentation

template<class T , class C >
void PFJetMonitor::fill ( const T jetCollection,
const C &  matchedJetCollection,
float &  minVal,
float &  maxVal 
)

fill histograms with all particle

Definition at line 78 of file PFJetMonitor.h.

References assert(), candBench_, createPFractionHistos_, dRMax_, reco::Candidate::eta(), reco::LeafCandidate::eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), fillOne(), histogramBooked_, i, Benchmark::isInRange(), metsig::jet, PFB::match(), matchCandBench_, matchCharge_, reco::Candidate::phi(), reco::LeafCandidate::phi(), reco::Candidate::pt(), reco::LeafCandidate::pt(), and findQualityFiles::size.

Referenced by PFJetDQMAnalyzer::analyze().

79  {
80 
81 
82  std::vector<int> matchIndices;
83  PFB::match( jetCollection, matchedJetCollection, matchIndices, matchCharge_, dRMax_ );
84 
85  for (unsigned int i = 0; i < (jetCollection).size(); i++) {
86  const reco::Jet& jet = jetCollection[i];
87 
88  if( !isInRange(jet.pt(), jet.eta(), jet.phi() ) ) continue;
89 
90  int iMatch = matchIndices[i];
91  assert(iMatch< static_cast<int>(matchedJetCollection.size()));
92 
93  if( iMatch != -1 ) {
94  const reco::Candidate& matchedJet = matchedJetCollection[ iMatch ];
95  if( !isInRange( matchedJet.pt(), matchedJet.eta(), matchedJet.phi() ) ) continue;
96  float ptRes = (jet.pt() - matchedJet.pt())/matchedJet.pt();
97 
98  if (ptRes > maxVal) maxVal = ptRes;
99  if (ptRes < minVal) minVal = ptRes;
100 
101  candBench_.fillOne(jet);
102  matchCandBench_.fillOne(jet, matchedJetCollection[ iMatch ]);
103  if (createPFractionHistos_ && histogramBooked_) fillOne(jet, matchedJetCollection[ iMatch ]);
104  }
105  }
106 }
int i
Definition: DBlmapReader.cc:9
bool histogramBooked_
Definition: PFJetMonitor.h:72
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:58
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
virtual double pt() const =0
transverse momentum
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
Base class for all types of Jets.
Definition: Jet.h:20
assert(m_qm.get())
virtual double phi() const final
momentum azimuthal angle
bool createPFractionHistos_
Definition: PFJetMonitor.h:71
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:59
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
virtual double eta() const final
momentum pseudorapidity
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
tuple size
Write out results.
virtual double phi() const =0
momentum azimuthal angle
bool matchCharge_
Definition: PFJetMonitor.h:70
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const final
transverse momentum
template<class T , class C >
void PFJetMonitor::fill ( const T candidateCollection,
const C &  matchedCandCollection,
float &  minVal,
float &  maxVal,
float &  jetpT,
const edm::ParameterSet parameterSet 
)

Definition at line 110 of file PFJetMonitor.h.

References assert(), candBench_, createPFractionHistos_, reco::deltaR(), deltaR_, dRMax_, reco::Candidate::eta(), reco::LeafCandidate::eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), fillOne(), histogramBooked_, i, Benchmark::isInRange(), j, metsig::jet, PFB::match(), matchCandBench_, matchCharge_, onlyTwoJets_, reco::Candidate::phi(), reco::LeafCandidate::phi(), EnergyCorrector::pt, reco::Candidate::pt(), and reco::LeafCandidate::pt().

112  {
113 
114  std::vector<int> matchIndices;
115  PFB::match( jetCollection, matchedJetCollection, matchIndices, matchCharge_, dRMax_ );
116  // now matchIndices[i] stores the j-th closest matched jet
117 
118  for( unsigned i=0; i<jetCollection.size(); ++i) {
119  // Count the number of jets with a larger energy = pT
120  unsigned int highJets = 0;
121  for( unsigned j=0; j<jetCollection.size(); ++j) {
122  if (j != i && jetCollection[j].pt() > jetCollection[i].pt()) highJets++;
123  }
124  if ( onlyTwoJets_ && highJets > 1 ) continue;
125 
126  const reco::Jet& jet = jetCollection[i];
127 
128  if( !isInRange(jet.pt(), jet.eta(), jet.phi() ) ) continue;
129 
130  int iMatch = matchIndices[i];
131  assert( iMatch < static_cast<int>(matchedJetCollection.size()) );
132 
133  if( iMatch != -1 ) {
134  const reco::Candidate& matchedJet = matchedJetCollection[ iMatch ];
135  if ( !isInRange( matchedJet.pt(), matchedJet.eta(), matchedJet.phi() ) ) continue;
136 
137  float ptRes = (jet.pt() - matchedJet.pt()) / matchedJet.pt();
138 
139  jetpT = jet.pt();
140  if (ptRes > maxVal) maxVal = ptRes;
141  if (ptRes < minVal) minVal = ptRes;
142 
143  candBench_.fillOne(jet); // fill pt eta phi and charge histos for MATCHED candidate jet
144  matchCandBench_.fillOne(jet, matchedJetCollection[iMatch], parameterSet); // fill delta_x_VS_y histos for matched couple
145  if (createPFractionHistos_ && histogramBooked_) fillOne(jet, matchedJetCollection[iMatch]); // book and fill delta_frac_VS_frac histos for matched couple
146  }
147 
148  for( unsigned j=0; j<matchedJetCollection.size(); ++j) // for DeltaR spectrum
149  if (deltaR_) deltaR_->Fill( reco::deltaR( jetCollection[i], matchedJetCollection[j] ) ) ;
150  } // end loop on jetCollection
151 }
int i
Definition: DBlmapReader.cc:9
bool histogramBooked_
Definition: PFJetMonitor.h:72
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:58
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
virtual double pt() const =0
transverse momentum
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
Base class for all types of Jets.
Definition: Jet.h:20
assert(m_qm.get())
virtual double phi() const final
momentum azimuthal angle
bool createPFractionHistos_
Definition: PFJetMonitor.h:71
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:59
int j
Definition: DBlmapReader.cc:9
bool onlyTwoJets_
Definition: PFJetMonitor.h:69
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
TH1F * deltaR_
Definition: PFJetMonitor.h:67
virtual double eta() const final
momentum pseudorapidity
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
virtual double phi() const =0
momentum azimuthal angle
bool matchCharge_
Definition: PFJetMonitor.h:70
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const final
transverse momentum
void PFJetMonitor::fillOne ( const reco::Jet jet,
const reco::Jet matchedJet 
)

Definition at line 162 of file PFJetMonitor.cc.

Referenced by fill().

void PFJetMonitor::setDirectory ( TDirectory *  dir)
virtual

set directory (to use in ROOT)

Reimplemented from Benchmark.

Definition at line 151 of file PFJetMonitor.cc.

void PFJetMonitor::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 at line 70 of file PFJetMonitor.cc.

Referenced by PFJetDQMAnalyzer::PFJetDQMAnalyzer().

void PFJetMonitor::setParameters ( float  dRMax,
bool  onlyTwoJets,
bool  matchCharge,
Benchmark::Mode  mode,
float  ptmin,
float  ptmax,
float  etamin,
float  etamax,
float  phimin,
float  phimax,
bool  fracHistoFlag = true 
)

Definition at line 85 of file PFJetMonitor.cc.

void PFJetMonitor::setParameters ( const edm::ParameterSet parameterSet)

set the parameters accessing them from ParameterSet

Definition at line 47 of file PFJetMonitor.cc.

void PFJetMonitor::setup ( DQMStore::IBooker b)

book histograms

Definition at line 130 of file PFJetMonitor.cc.

Referenced by PFJetDQMAnalyzer::bookHistograms().

void PFJetMonitor::setup ( DQMStore::IBooker b,
const edm::ParameterSet parameterSet 
)

Definition at line 104 of file PFJetMonitor.cc.

Member Data Documentation

CandidateBenchmark PFJetMonitor::candBench_
protected

Definition at line 58 of file PFJetMonitor.h.

Referenced by fill().

bool PFJetMonitor::createPFractionHistos_
protected

Definition at line 71 of file PFJetMonitor.h.

Referenced by fill().

TH2F* PFJetMonitor::delta_frac_VS_frac_charged_hadron_
protected

Definition at line 64 of file PFJetMonitor.h.

TH2F* PFJetMonitor::delta_frac_VS_frac_electron_
protected

Definition at line 63 of file PFJetMonitor.h.

TH2F* PFJetMonitor::delta_frac_VS_frac_muon_
protected

Definition at line 61 of file PFJetMonitor.h.

TH2F* PFJetMonitor::delta_frac_VS_frac_neutral_hadron_
protected

Definition at line 65 of file PFJetMonitor.h.

TH2F* PFJetMonitor::delta_frac_VS_frac_photon_
protected

Definition at line 62 of file PFJetMonitor.h.

TH1F* PFJetMonitor::deltaR_
protected

Definition at line 67 of file PFJetMonitor.h.

Referenced by fill().

float PFJetMonitor::dRMax_
protected

Definition at line 68 of file PFJetMonitor.h.

Referenced by fill().

bool PFJetMonitor::histogramBooked_
protected

Definition at line 72 of file PFJetMonitor.h.

Referenced by fill().

MatchCandidateBenchmark PFJetMonitor::matchCandBench_
protected

Definition at line 59 of file PFJetMonitor.h.

Referenced by fill().

bool PFJetMonitor::matchCharge_
protected

Definition at line 70 of file PFJetMonitor.h.

Referenced by fill().

bool PFJetMonitor::onlyTwoJets_
protected

Definition at line 69 of file PFJetMonitor.h.

Referenced by fill().