CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFJetMonitor.h
Go to the documentation of this file.
1 #ifndef DQMOffline_PFTau_PFJetMonitor_h
2 #define DQMOffline_PFTau_PFJetMonitor_h
3 
8 
10 
11 #include <vector>
12 class PFJetMonitor : public Benchmark {
13 
14  public:
15 
16  PFJetMonitor( float dRMax = 0.3,
17  bool matchCharge = true,
19  :
20  Benchmark(mode),
22  dRMax_(dRMax), matchCharge_(matchCharge) {}
23 
24  virtual ~PFJetMonitor();
25 
27  void setParameters( const edm::ParameterSet& parameterSet);
28 
30  void setDirectory(TDirectory* dir);
31 
33  void setup();
34 
36  void setup(const edm::ParameterSet & parameterSet);
37 
39  template< class T, class C>
40  void fill(const T& jetCollection,
41  const C& matchedJetCollection, float& minVal, float& maxVal);
42 
43 
44  void fillOne(const reco::Jet& jet,
45  const reco::Jet& matchedJet);
46 
47  protected:
50 
56 
57  float dRMax_;
60 
61 };
62 
64 template< class T, class C>
65 void PFJetMonitor::fill(const T& jetCollection,
66  const C& matchedJetCollection, float& minVal, float& maxVal) {
67 
68 
69  std::vector<int> matchIndices;
70  PFB::match( jetCollection, matchedJetCollection, matchIndices,
72 
73  for (unsigned int i = 0; i < (jetCollection).size(); i++) {
74  const reco::Jet& jet = jetCollection[i];
75 
76  if( !isInRange(jet.pt(), jet.eta(), jet.phi() ) ) continue;
77 
78  int iMatch = matchIndices[i];
79  assert(iMatch< static_cast<int>(matchedJetCollection.size()));
80 
81  if( iMatch!=-1 ) {
82  const reco::Candidate& matchedJet = matchedJetCollection[ iMatch ];
83  if(!isInRange(matchedJet.pt(),matchedJet.eta(),matchedJet.phi() ) ) continue;
84  float ptRes = (jet.pt() - matchedJet.pt())/matchedJet.pt();
85 
86  if (ptRes > maxVal) maxVal = ptRes;
87  if (ptRes < minVal) minVal = ptRes;
88 
89  candBench_.fillOne(jet);
90  matchCandBench_.fillOne(jet, matchedJetCollection[ iMatch ]);
91  if (createPFractionHistos_) fillOne(jet, matchedJetCollection[ iMatch ]);
92  }
93  }
94 }
95 #endif
int i
Definition: DBlmapReader.cc:9
To plot Candidate quantities.
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
To plot Candidate quantities.
CandidateBenchmark candBench_
Definition: PFJetMonitor.h:48
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:13
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
Definition: PFJetMonitor.cc:71
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:21
void fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal)
fill histograms with all particle
Definition: PFJetMonitor.h:65
abstract base class
Definition: Benchmark.h:20
bool createPFractionHistos_
Definition: PFJetMonitor.h:59
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:51
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:49
virtual double eta() const
momentum pseudorapidity
virtual ~PFJetMonitor()
Definition: PFJetMonitor.cc:13
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:55
void setup()
book histograms
Definition: PFJetMonitor.cc:53
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.h:16
void setParameters(const edm::ParameterSet &parameterSet)
set the parameters
Definition: PFJetMonitor.cc:16
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
Definition: PFJetMonitor.cc:78
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:52
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:53
virtual double pt() const
transverse momentum
int mode
Definition: AMPTWrapper.h:139
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:54
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
virtual double phi() const
momentum azimuthal angle
tuple size
Write out results.
virtual double phi() const =0
momentum azimuthal angle
bool matchCharge_
Definition: PFJetMonitor.h:58
virtual double eta() const =0
momentum pseudorapidity