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
PFCandidateMonitor Class Reference

#include <PFCandidateMonitor.h>

Inheritance diagram for PFCandidateMonitor:
Benchmark

Public Member Functions

template<class T , class C >
void fill (const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet)
 fill histograms with all particle More...
 
void fillOne (const reco::Candidate &cand)
 
 PFCandidateMonitor (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 refHistoFlag)
 set the parameters locally More...
 
void setParameters (const edm::ParameterSet &parameterSet)
 set the parameters accessing them from ParameterSet More...
 
void setup ()
 book histograms More...
 
void setup (const edm::ParameterSet &parameterSet)
 book histograms More...
 
virtual ~PFCandidateMonitor ()
 
- 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 createEfficiencyHistos_
 
bool createReferenceHistos_
 
TH1F * deltaR_
 
float dRMax_
 
TH1F * eta_gen_
 
TH1F * eta_ref_
 
bool histogramBooked_
 
MatchCandidateBenchmark matchCandBench_
 
bool matchCharge_
 
bool matching_done_
 
TH1F * phi_gen_
 
TH1F * phi_ref_
 
TH1F * pt_gen_
 
TH1F * pt_ref_
 
- 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 }
 
- Static Public Attributes inherited from Benchmark
static DQMStoreDQM_ = 0
 
- Protected Member Functions inherited from Benchmark
TH1F * book1D (const char *histname, const char *title, int nbins, float xmin, float xmax)
 book a 1D histogram, either with DQM or plain root. More...
 
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. More...
 
TH2F * book2D (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. More...
 
TProfile * bookProfile (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. More...
 
TProfile * bookProfile (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. More...
 

Detailed Description

Definition at line 16 of file PFCandidateMonitor.h.

Constructor & Destructor Documentation

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

Definition at line 15 of file PFCandidateMonitor.cc.

References createReferenceHistos_, deltaR_, eta_gen_, eta_ref_, histogramBooked_, phi_gen_, phi_ref_, pt_gen_, pt_ref_, and Benchmark::setRange().

15  :
16  Benchmark(mode),
17  candBench_(mode),
19  dRMax_(dRMax),
20  matchCharge_(matchCharge) {
21 
22  setRange( 0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
23 
24  pt_gen_ = 0;
25  eta_gen_ = 0;
26  phi_gen_ = 0;
27 
28  pt_ref_ = 0;
29  eta_ref_ = 0;
30  phi_ref_ = 0;
31 
32  deltaR_ = 0;
33 
34  createReferenceHistos_ = false;
35  histogramBooked_ = false;
36 }
Benchmark(Mode mode=DEFAULT)
Definition: Benchmark.h:42
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
MatchCandidateBenchmark matchCandBench_
PFCandidateMonitor::~PFCandidateMonitor ( )
virtual

Definition at line 40 of file PFCandidateMonitor.cc.

40 {}

Member Function Documentation

template<class T , class C >
void PFCandidateMonitor::fill ( const T candidateCollection,
const C &  matchedCandCollection,
float &  minVal,
float &  maxVal,
const edm::ParameterSet parameterSet 
)

fill histograms with all particle

Definition at line 79 of file PFCandidateMonitor.h.

References candBench_, createEfficiencyHistos_, createReferenceHistos_, reco::deltaR(), deltaR_, dRMax_, reco::Candidate::eta(), eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), fillOne(), i, Benchmark::isInRange(), j, PFB::match(), matchCandBench_, matchCharge_, matching_done_, phi, reco::Candidate::phi(), RecoTauCleanerPlugins::pt, reco::Candidate::pt(), and findQualityFiles::size.

Referenced by PFCandidateDQMAnalyzer::analyze().

81  {
82 
83  matching_done_ = false;
85  for( unsigned i=0; i<candCollection.size(); ++i) {
86  if( !isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi() ) ) continue;
87  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for UNMATCHED generated candidate
88 
89  for( unsigned j=0; j<matchedCandCollection.size(); ++j) // for DeltaR spectrum
90  if (deltaR_) deltaR_->Fill( reco::deltaR( candCollection[i], matchedCandCollection[j] ) ) ;
91  }
92  }
93 
94  std::vector<int> matchIndices;
95  PFB::match( candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_ );
96  // now matchIndices[i] stores the j-th closest matched jet
97  matching_done_ = true;
98 
99  for (unsigned int i = 0; i < (candCollection).size(); i++) {
100  const reco::Candidate& cand = candCollection[i];
101 
102  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) continue;
103 
104  int iMatch = matchIndices[i];
105  assert(iMatch< static_cast<int>(matchedCandCollection.size()));
106 
107  if( iMatch!=-1 ) {
108  const reco::Candidate& matchedCand = matchedCandCollection[ iMatch ];
109  if(!isInRange(matchedCand.pt(),matchedCand.eta(),matchedCand.phi() ) ) continue;
110  //std::cout <<"PFJet pT " <<cand.pt() <<" eta " <<cand.eta() <<" phi " <<cand.phi() ;
111  //std::cout <<"\nmatched genJet pT " <<matchedCand.pt() <<" eta " <<matchedCand.eta() <<" phi " <<matchedCand.phi() <<"\n" <<std::endl ;
112  float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
113 
114  if (ptRes > maxVal) maxVal = ptRes;
115  if (ptRes < minVal) minVal = ptRes;
116 
117  if ( !createEfficiencyHistos_ ) {
118  candBench_.fillOne(cand); // fill pt, eta phi and charge histos for MATCHED candidate
119  //matchCandBench_.fillOne(cand, matchedCand); // fill delta_x_VS_y histos for matched couple
120  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
121  if (createReferenceHistos_) fillOne(matchedCand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
122  }
123  else {
124  candBench_.fillOne(matchedCand); // fill pt, eta phi and charge histos for MATCHED candidate
125  //matchCandBench_.fillOne(matchedCand, cand); // fill delta_x_VS_y histos for matched couple
126  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
127  if (createReferenceHistos_) fillOne(cand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
128  }
129 
130  }
131  }
132 }
int i
Definition: DBlmapReader.cc:9
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:14
virtual float eta() const =0
momentum pseudorapidity
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
virtual float phi() const =0
momentum azimuthal angle
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
T eta() const
virtual float pt() const =0
transverse momentum
int j
Definition: DBlmapReader.cc:9
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
MatchCandidateBenchmark matchCandBench_
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:59
tuple size
Write out results.
Definition: DDAxes.h:10
void PFCandidateMonitor::fillOne ( const reco::Candidate cand)

Definition at line 163 of file PFCandidateMonitor.cc.

References createEfficiencyHistos_, createReferenceHistos_, reco::Candidate::eta(), eta_gen_, eta_ref_, histogramBooked_, matching_done_, reco::Candidate::phi(), phi_gen_, phi_ref_, reco::Candidate::pt(), pt_gen_, and pt_ref_.

Referenced by fill().

163  {
164 
165  if (matching_done_) {
167  if (pt_ref_) pt_ref_->Fill(cand.pt());
168  if (eta_ref_) eta_ref_->Fill(cand.eta() );
169  if (phi_ref_) phi_ref_->Fill(cand.phi() );
170  }
172  if (pt_gen_) pt_gen_->Fill(cand.pt());
173  if (eta_gen_) eta_gen_->Fill(cand.eta() );
174  if (phi_gen_) phi_gen_->Fill(cand.phi() );
175  }
176 
177 }
virtual float eta() const =0
momentum pseudorapidity
virtual float phi() const =0
momentum azimuthal angle
virtual float pt() const =0
transverse momentum
void PFCandidateMonitor::setDirectory ( TDirectory *  dir)
virtual

set directory (to use in ROOT)

Reimplemented from Benchmark.

Definition at line 154 of file PFCandidateMonitor.cc.

References candBench_, matchCandBench_, and Benchmark::setDirectory().

154  {
156 
159 }
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:19
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
dbl *** dir
Definition: mlp_gen.cc:35
void PFCandidateMonitor::setParameters ( float  dRMax,
bool  matchCharge,
Benchmark::Mode  mode,
float  ptmin,
float  ptmax,
float  etamin,
float  etamax,
float  phimin,
float  phimax,
bool  refHistoFlag 
)

set the parameters locally

Definition at line 66 of file PFCandidateMonitor.cc.

References candBench_, createReferenceHistos_, dRMax_, matchCandBench_, matchCharge_, alignBH_cfg::mode, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

Referenced by PFCandidateDQMAnalyzer::PFCandidateDQMAnalyzer().

68  {
69  dRMax_ = dRMax;
70  matchCharge_ = matchCharge;
71  mode_ = mode;
72  createReferenceHistos_ = refHistoFlag;
73 
74  setRange( ptmin, ptmax, etamin, etamax, phimin, phimax );
75 
78 }
void setParameters(Mode mode)
Definition: Benchmark.h:50
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
double ptmin
Definition: HydjetWrapper.h:85
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:96
void PFCandidateMonitor::setParameters ( const edm::ParameterSet parameterSet)

set the parameters accessing them from ParameterSet

Definition at line 45 of file PFCandidateMonitor.cc.

References candBench_, createEfficiencyHistos_, createReferenceHistos_, dRMax_, edm::ParameterSet::getParameter(), matchCandBench_, matchCharge_, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

45  {
46 
47  dRMax_ = parameterSet.getParameter<double>( "deltaRMax" );
48  matchCharge_ = parameterSet.getParameter<bool>( "matchCharge" );
49  mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" );
50  createReferenceHistos_ = parameterSet.getParameter<bool>( "CreateReferenceHistos" );
51  createEfficiencyHistos_ = parameterSet.getParameter<bool>( "CreateEfficiencyHistos" );
52 
53  setRange( parameterSet.getParameter<double>("ptMin"),
54  parameterSet.getParameter<double>("ptMax"),
55  parameterSet.getParameter<double>("etaMin"),
56  parameterSet.getParameter<double>("etaMax"),
57  parameterSet.getParameter<double>("phiMin"),
58  parameterSet.getParameter<double>("phiMax") );
59 
62 }
T getParameter(std::string const &) const
void setParameters(Mode mode)
Definition: Benchmark.h:50
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:52
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:96
void PFCandidateMonitor::setup ( void  )

book histograms

Definition at line 130 of file PFCandidateMonitor.cc.

References Benchmark::book1D(), candBench_, createEfficiencyHistos_, createReferenceHistos_, eta_gen_, eta_ref_, histogramBooked_, Benchmark::PhaseSpace::m, Benchmark::PhaseSpace::M, matchCandBench_, Benchmark::PhaseSpace::n, phi_gen_, phi_ref_, pt_gen_, pt_ref_, CandidateBenchmark::setup(), and MatchCandidateBenchmark::setup().

Referenced by PFCandidateDQMAnalyzer::beginJob().

130  {
131  candBench_.setup();
133 
135  PhaseSpace ptPS(100,0,100);
136  PhaseSpace phiPS(360, -3.1416, 3.1416);
137  PhaseSpace etaPS(100, -5,5);
138 
139  pt_ref_ = book1D("pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
140  if (createEfficiencyHistos_) pt_gen_ = book1D("pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
141 
142  eta_ref_ = book1D("eta_ref_", "#eta_ref;#eta", etaPS.n, etaPS.m, etaPS.M);
143  if (createEfficiencyHistos_) eta_gen_ = book1D("eta_gen_", "#eta_gen;#eta", etaPS.n, etaPS.m, etaPS.M);
144 
145  phi_ref_ = book1D("phi_ref_", "#phi_ref;#phi", phiPS.n, phiPS.m, phiPS.M);
146  if (createEfficiencyHistos_) phi_gen_ = book1D("phi_gen_", "#phi_gen;#phi", phiPS.n, phiPS.m, phiPS.M);
147 
148  histogramBooked_ = true;
149  }
150 }
void setup()
book histograms
void setup()
book histograms
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
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
void PFCandidateMonitor::setup ( const edm::ParameterSet parameterSet)

book histograms

Definition at line 82 of file PFCandidateMonitor.cc.

References Benchmark::book1D(), candBench_, createEfficiencyHistos_, createReferenceHistos_, deltaR_, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eta_gen_, eta_ref_, edm::ParameterSet::getParameter(), histogramBooked_, matchCandBench_, phi_gen_, phi_ref_, pt_gen_, pt_ref_, CandidateBenchmark::setup(), and MatchCandidateBenchmark::setup().

82  {
83  candBench_.setup(parameterSet);
84  matchCandBench_.setup(parameterSet);
85 
87  edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
88  edm::ParameterSet etaPS = parameterSet.getParameter<edm::ParameterSet>("EtaHistoParameter");
89  edm::ParameterSet phiPS = parameterSet.getParameter<edm::ParameterSet>("PhiHistoParameter");
90 
91  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
92 
93  if (ptPS.getParameter<bool>("switchOn")) {
94  pt_ref_ = book1D("pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
95  ptPS.getParameter<double>("xMin"),
96  ptPS.getParameter<double>("xMax"));
97  if (createEfficiencyHistos_) pt_gen_ = book1D( "pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
98  ptPS.getParameter<double>("xMin"),
99  ptPS.getParameter<double>("xMax") ) ;
100  }
101 
102  if (etaPS.getParameter<bool>("switchOn")) {
103  eta_ref_ = book1D("eta_ref_", "#eta_ref;#eta", etaPS.getParameter<int32_t>("nBin"),
104  etaPS.getParameter<double>("xMin"),
105  etaPS.getParameter<double>("xMax"));
106  if (createEfficiencyHistos_) eta_gen_ = book1D( "eta_gen_", "#eta_gen;#eta", etaPS.getParameter<int32_t>("nBin"),
107  etaPS.getParameter<double>("xMin"),
108  etaPS.getParameter<double>("xMax") ) ;
109  }
110  if (phiPS.getParameter<bool>("switchOn")) {
111  phi_ref_ = book1D("phi_ref_", "#phi_ref;#phi", phiPS.getParameter<int32_t>("nBin"),
112  phiPS.getParameter<double>("xMin"),
113  phiPS.getParameter<double>("xMax"));
114  if (createEfficiencyHistos_) phi_gen_ = book1D( "phi_gen_", "#phi_gen;#phi", phiPS.getParameter<int32_t>("nBin"),
115  phiPS.getParameter<double>("xMin"),
116  phiPS.getParameter<double>("xMax") ) ;
117  }
118  if ( createEfficiencyHistos_ && dR.getParameter<bool>("switchOn") )
119  deltaR_ = book1D("deltaR_", "#DeltaR;#DeltaR",
120  dR.getParameter<int32_t>("nBin"),
121  dR.getParameter<double>("xMin"),
122  dR.getParameter<double>("xMax"));
123 
124  histogramBooked_ = true;
125  }
126 }
T getParameter(std::string const &) const
void setup()
book histograms
void setup()
book histograms
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
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

Member Data Documentation

CandidateBenchmark PFCandidateMonitor::candBench_
protected

Definition at line 54 of file PFCandidateMonitor.h.

Referenced by fill(), setDirectory(), setParameters(), and setup().

bool PFCandidateMonitor::createEfficiencyHistos_
protected

Definition at line 72 of file PFCandidateMonitor.h.

Referenced by fill(), fillOne(), setParameters(), and setup().

bool PFCandidateMonitor::createReferenceHistos_
protected

Definition at line 68 of file PFCandidateMonitor.h.

Referenced by fill(), fillOne(), PFCandidateMonitor(), setParameters(), and setup().

TH1F* PFCandidateMonitor::deltaR_
protected

Definition at line 65 of file PFCandidateMonitor.h.

Referenced by fill(), PFCandidateMonitor(), and setup().

float PFCandidateMonitor::dRMax_
protected

Definition at line 66 of file PFCandidateMonitor.h.

Referenced by fill(), and setParameters().

TH1F* PFCandidateMonitor::eta_gen_
protected

Definition at line 58 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::eta_ref_
protected

Definition at line 62 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

bool PFCandidateMonitor::histogramBooked_
protected

Definition at line 69 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

MatchCandidateBenchmark PFCandidateMonitor::matchCandBench_
protected

Definition at line 55 of file PFCandidateMonitor.h.

Referenced by fill(), setDirectory(), setParameters(), and setup().

bool PFCandidateMonitor::matchCharge_
protected

Definition at line 67 of file PFCandidateMonitor.h.

Referenced by fill(), and setParameters().

bool PFCandidateMonitor::matching_done_
protected

Definition at line 71 of file PFCandidateMonitor.h.

Referenced by fill(), and fillOne().

TH1F* PFCandidateMonitor::phi_gen_
protected

Definition at line 59 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::phi_ref_
protected

Definition at line 63 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::pt_gen_
protected

Definition at line 57 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::pt_ref_
protected

Definition at line 61 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().