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