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 
13 #include <TH1.h> //needed by the deltaR->Fill() call
14 
15 class PFJetMonitor : public Benchmark {
16 
17  public:
18 
19 
20  PFJetMonitor( float dRMax = 0.3, bool matchCharge = true,
22 
23  virtual ~PFJetMonitor();
24 
26  void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
27  float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax,
28  bool fracHistoFlag=true);
29 
30  void setParameters(float dRMax, bool onlyTwoJets, bool matchCharge, Benchmark::Mode mode,
31  float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax,
32  bool fracHistoFlag=true);
33 
36 
38  void setDirectory(TDirectory* dir);
39 
41  void setup(DQMStore::IBooker& b);
42  void setup(DQMStore::IBooker& b, const edm::ParameterSet & parameterSet);
43 
45  template< class T, class C>
46  void fill(const T& jetCollection, const C& matchedJetCollection,
47  float& minVal, float& maxVal);
48 
49  template< class T, class C>
50  void fill(const T& candidateCollection, const C& matchedCandCollection,
51  float& minVal, float& maxVal, float& jetpT,
52  const edm::ParameterSet & parameterSet);
53 
54  void fillOne(const reco::Jet& jet,
55  const reco::Jet& matchedJet);
56 
57  protected:
60 
66 
67  TH1F* deltaR_;
68  float dRMax_;
73 
74 };
75 
77 template< class T, class C>
78  void PFJetMonitor::fill(const T& jetCollection,
79  const C& matchedJetCollection, float& minVal, float& maxVal) {
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 }
107 
108 
109 template< class T, class C>
110  void PFJetMonitor::fill(const T& jetCollection,
111  const C& matchedJetCollection, float& minVal, float& maxVal, float& jetpT,
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 }
152 #endif
int i
Definition: DBlmapReader.cc:9
bool histogramBooked_
Definition: PFJetMonitor.h:72
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:58
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 setDirectory(TDirectory *dir)
set directory (to use in ROOT)
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())
void fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal)
fill histograms with all particle
Definition: PFJetMonitor.h:78
abstract base class
Definition: Benchmark.h:22
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
bool createPFractionHistos_
Definition: PFJetMonitor.h:71
TH2F * delta_frac_VS_frac_muon_
Definition: PFJetMonitor.h:61
MatchCandidateBenchmark matchCandBench_
Definition: PFJetMonitor.h:59
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
void setup(DQMStore::IBooker &b)
book histograms
virtual ~PFJetMonitor()
Definition: PFJetMonitor.cc:41
TH2F * delta_frac_VS_frac_neutral_hadron_
Definition: PFJetMonitor.h:65
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:70
int j
Definition: DBlmapReader.cc:9
bool onlyTwoJets_
Definition: PFJetMonitor.h:69
PFJetMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: PFJetMonitor.cc:16
void fillOne(const reco::Jet &jet, const reco::Jet &matchedJet)
TH2F * delta_frac_VS_frac_photon_
Definition: PFJetMonitor.h:62
TH2F * delta_frac_VS_frac_electron_
Definition: PFJetMonitor.h:63
TH1F * deltaR_
Definition: PFJetMonitor.h:67
double b
Definition: hdecay.h:120
double ptmin
Definition: HydjetWrapper.h:85
TH2F * delta_frac_VS_frac_charged_hadron_
Definition: PFJetMonitor.h:64
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
long double T
virtual double phi() const
momentum azimuthal angle
tuple size
Write out results.
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
virtual double phi() const =0
momentum azimuthal angle
bool matchCharge_
Definition: PFJetMonitor.h:70
virtual double eta() const =0
momentum pseudorapidity