CMS 3D CMS Logo

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 
13 #include <TH1.h> //needed by the deltaR->Fill() call
14 
15 class PFJetMonitor : public Benchmark {
16 public:
17  PFJetMonitor(float dRMax = 0.3, bool matchCharge = true, Benchmark::Mode mode = Benchmark::DEFAULT);
18 
19  ~PFJetMonitor() override;
20 
22  void setParameters(float dRMax,
23  bool matchCharge,
25  float ptmin,
26  float ptmax,
27  float etamin,
28  float etamax,
29  float phimin,
30  float phimax,
31  bool fracHistoFlag = true);
32 
33  void setParameters(float dRMax,
34  bool onlyTwoJets,
35  bool matchCharge,
36  Benchmark::Mode mode,
37  float ptmin,
38  float ptmax,
39  float etamin,
40  float etamax,
41  float phimin,
42  float phimax,
43  bool fracHistoFlag = true);
44 
47 
49  void setDirectory(TDirectory *dir) override;
50 
52  void setup(DQMStore::IBooker &b);
53  void setup(DQMStore::IBooker &b, const edm::ParameterSet &parameterSet);
54 
56  template <class T, class C>
57  void fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal);
58 
59  template <class T, class C>
60  void fill(const T &candidateCollection,
61  const C &matchedCandCollection,
62  float &minVal,
63  float &maxVal,
64  float &jetpT,
65  const edm::ParameterSet &parameterSet);
66 
67  void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet);
68 
69 protected:
72 
78 
79  TH1F *deltaR_;
80  float dRMax_;
85 };
86 
88 template <class T, class C>
89 void PFJetMonitor::fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal) {
90  std::vector<int> matchIndices;
91  PFB::match(jetCollection, matchedJetCollection, matchIndices, matchCharge_, dRMax_);
92 
93  for (unsigned int i = 0; i < (jetCollection).size(); i++) {
94  const reco::Jet &jet = jetCollection[i];
95 
96  if (!isInRange(jet.pt(), jet.eta(), jet.phi()))
97  continue;
98 
99  int iMatch = matchIndices[i];
100  assert(iMatch < static_cast<int>(matchedJetCollection.size()));
101 
102  if (iMatch != -1) {
103  const reco::Candidate &matchedJet = matchedJetCollection[iMatch];
104  if (!isInRange(matchedJet.pt(), matchedJet.eta(), matchedJet.phi()))
105  continue;
106  float ptRes = (jet.pt() - matchedJet.pt()) / matchedJet.pt();
107 
108  if (ptRes > maxVal)
109  maxVal = ptRes;
110  if (ptRes < minVal)
111  minVal = ptRes;
112 
113  candBench_.fillOne(jet);
114  matchCandBench_.fillOne(jet, matchedJetCollection[iMatch]);
116  fillOne(jet, matchedJetCollection[iMatch]);
117  }
118  }
119 }
120 
121 template <class T, class C>
122 void PFJetMonitor::fill(const T &jetCollection,
123  const C &matchedJetCollection,
124  float &minVal,
125  float &maxVal,
126  float &jetpT,
128  std::vector<int> matchIndices;
129  PFB::match(jetCollection, matchedJetCollection, matchIndices, matchCharge_, dRMax_);
130  // now matchIndices[i] stores the j-th closest matched jet
131 
132  for (unsigned i = 0; i < jetCollection.size(); ++i) {
133  // Count the number of jets with a larger energy = pT
134  unsigned int highJets = 0;
135  for (unsigned j = 0; j < jetCollection.size(); ++j) {
136  if (j != i && jetCollection[j].pt() > jetCollection[i].pt())
137  highJets++;
138  }
139  if (onlyTwoJets_ && highJets > 1)
140  continue;
141 
142  const reco::Jet &jet = jetCollection[i];
143 
144  if (!isInRange(jet.pt(), jet.eta(), jet.phi()))
145  continue;
146 
147  int iMatch = matchIndices[i];
148  assert(iMatch < static_cast<int>(matchedJetCollection.size()));
149 
150  if (iMatch != -1) {
151  const reco::Candidate &matchedJet = matchedJetCollection[iMatch];
152  if (!isInRange(matchedJet.pt(), matchedJet.eta(), matchedJet.phi()))
153  continue;
154 
155  float ptRes = (jet.pt() - matchedJet.pt()) / matchedJet.pt();
156 
157  jetpT = jet.pt();
158  if (ptRes > maxVal)
159  maxVal = ptRes;
160  if (ptRes < minVal)
161  minVal = ptRes;
162 
163  candBench_.fillOne(jet); // fill pt eta phi and charge histos for MATCHED candidate jet
165  matchedJetCollection[iMatch],
166  parameterSet); // fill delta_x_VS_y histos for matched couple
168  fillOne(jet,
169  matchedJetCollection[iMatch]); // book and fill delta_frac_VS_frac
170  // histos for matched couple
171  }
172 
173  for (unsigned j = 0; j < matchedJetCollection.size(); ++j) // for DeltaR spectrum
174  if (deltaR_)
175  deltaR_->Fill(reco::deltaR(jetCollection[i], matchedJetCollection[j]));
176  } // end loop on jetCollection
177 }
178 #endif
void setDirectory(TDirectory *dir) override
set directory (to use in ROOT)
size
Write out results.
bool histogramBooked_
Definition: PFJetMonitor.h:84
To plot Candidate quantities.
double eta() const final
momentum pseudorapidity
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:70
bool onlyTwoJets
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
Base class for all types of Jets.
Definition: Jet.h:20
void fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal)
fill histograms with all particle
Definition: PFJetMonitor.h:89
abstract base class
Definition: Benchmark.h:21
double pt() const final
transverse momentum
bool createPFractionHistos_
Definition: PFJetMonitor.h:83
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:73
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:71
void setup(DQMStore::IBooker &b)
book histograms
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:77
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:60
bool onlyTwoJets_
Definition: PFJetMonitor.h:81
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)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:74
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:75
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
TH1F * deltaR_
Definition: PFJetMonitor.h:79
double b
Definition: hdecay.h:120
double ptmin
Definition: HydjetWrapper.h:90
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:76
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:50
long double T
double phi() const final
momentum azimuthal angle
virtual double phi() const =0
momentum azimuthal angle
~PFJetMonitor() override
Definition: PFJetMonitor.cc:34
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
bool matchCharge_
Definition: PFJetMonitor.h:82