CMS 3D CMS Logo

PFCandidateChecker.cc
Go to the documentation of this file.
1 
17 
18 namespace edm {
19  class EventSetup;
20  class Run;
21 } // namespace edm
22 
24 public:
25  explicit PFCandidateChecker(const edm::ParameterSet&);
26 
27  void analyze(const edm::Event&, const edm::EventSetup&) override;
28 
29 private:
30  void printJets(const reco::PFJetCollection& pfJetsReco, const reco::PFJetCollection& pfJetsReReco) const;
31 
32  void printMet(const reco::PFCandidateCollection& pfReco, const reco::PFCandidateCollection& pfReReco) const;
33 
34  void printElementsInBlocks(const reco::PFCandidate& cand, std::ostream& out = std::cout) const;
35 
41 
43  double deltaEMax_;
44  double deltaEtaMax_;
45  double deltaPhiMax_;
46 
48  bool verbose_;
49 
52 
54  bool rankByPt_;
55 
57  unsigned entry_;
58 
59  static bool greaterPt(const reco::PFCandidate& a, const reco::PFCandidate& b) { return (a.pt() > b.pt()); }
60 };
61 
63 
64 using namespace std;
65 using namespace edm;
66 using namespace reco;
67 
69  inputTokenPFCandidatesReco_ =
70  consumes<reco::PFCandidateCollection>(iConfig.getParameter<InputTag>("pfCandidatesReco"));
71 
72  inputTokenPFCandidatesReReco_ =
73  consumes<reco::PFCandidateCollection>(iConfig.getParameter<InputTag>("pfCandidatesReReco"));
74 
75  inputTokenPFJetsReco_ = consumes<reco::PFJetCollection>(iConfig.getParameter<InputTag>("pfJetsReco"));
76 
77  inputTokenPFJetsReReco_ = consumes<reco::PFJetCollection>(iConfig.getParameter<InputTag>("pfJetsReReco"));
78 
79  deltaEMax_ = iConfig.getParameter<double>("deltaEMax");
80 
81  deltaEtaMax_ = iConfig.getParameter<double>("deltaEtaMax");
82 
83  deltaPhiMax_ = iConfig.getParameter<double>("deltaPhiMax");
84 
85  verbose_ = iConfig.getUntrackedParameter<bool>("verbose", false);
86 
87  printBlocks_ = iConfig.getUntrackedParameter<bool>("printBlocks", false);
88 
89  rankByPt_ = iConfig.getUntrackedParameter<bool>("rankByPt", false);
90 
91  entry_ = 0;
92 
93  LogDebug("PFCandidateChecker") << " input collections : " << iConfig.getParameter<InputTag>("pfCandidatesReco") << " "
94  << iConfig.getParameter<InputTag>("pfCandidatesReReco");
95 }
96 
97 void PFCandidateChecker::analyze(const Event& iEvent, const EventSetup& iSetup) {
98  LogDebug("PFCandidateChecker") << "START event: " << iEvent.id().event() << " in run " << iEvent.id().run() << endl;
99 
100  // get PFCandidates
101 
102  const auto& pfCandidatesReco = iEvent.getHandle(inputTokenPFCandidatesReco_);
103  const auto& pfCandidatesReReco = iEvent.getHandle(inputTokenPFCandidatesReReco_);
104  const auto& pfJetsReco = iEvent.getHandle(inputTokenPFJetsReco_);
105  const auto& pfJetsReReco = iEvent.getHandle(inputTokenPFJetsReReco_);
106 
107  reco::PFCandidateCollection pfReco, pfReReco;
108 
109  // to sort, one needs to copy
110  if (rankByPt_) {
111  pfReco = *pfCandidatesReco;
112  pfReReco = *pfCandidatesReReco;
113  sort(pfReco.begin(), pfReco.end(), greaterPt);
114  sort(pfReReco.begin(), pfReReco.end(), greaterPt);
115  }
116 
117  unsigned minSize = pfReco.size() < pfReReco.size() ? pfReco.size() : pfReReco.size();
118  bool differentCand = false;
119  bool differentSize = pfReco.size() != pfReReco.size();
120  if (differentSize)
121  std::cout << "+++WARNING+++ PFCandidate size changed for entry " << entry_ << " !" << endl
122  << " - RECO size : " << pfReco.size() << endl
123  << " - Re-RECO size : " << pfReReco.size() << endl;
124 
125  unsigned npr = 0;
126  for (unsigned i = 0; i < minSize; i++) {
127  const reco::PFCandidate& candReco = (rankByPt_) ? pfReco[i] : (*pfCandidatesReco)[i];
128  const reco::PFCandidate& candReReco = (rankByPt_) ? pfReReco[i] : (*pfCandidatesReReco)[i];
129 
130  double deltaE = (candReReco.energy() - candReco.energy()) / (candReReco.energy() + candReco.energy());
131  double deltaEta = candReReco.eta() - candReco.eta();
132  double deltaPhi = candReReco.phi() - candReco.phi();
133  if (fabs(deltaE) > deltaEMax_ || fabs(deltaEta) > deltaEtaMax_ || fabs(deltaPhi) > deltaPhiMax_) {
134  differentCand = true;
135  std::cout << "+++WARNING+++ PFCandidate " << i << " changed for entry " << entry_ << " ! " << std::endl
136  << " - RECO : " << candReco << std::endl
137  << " - Re-RECO : " << candReReco << std::endl
138  << " DeltaE = : " << deltaE << std::endl
139  << " DeltaEta = : " << deltaEta << std::endl
140  << " DeltaPhi = : " << deltaPhi << std::endl
141  << std::endl;
142  if (printBlocks_) {
143  std::cout << "Elements in Block for RECO: " << std::endl;
144  printElementsInBlocks(candReco);
145  std::cout << "Elements in Block for Re-RECO: " << std::endl;
146  printElementsInBlocks(candReReco);
147  }
148  if (++npr == 5)
149  break;
150  }
151  }
152 
153  if (differentSize || differentCand) {
154  printJets(*pfJetsReco, *pfJetsReReco);
155  printMet(pfReco, pfReReco);
156  }
157 
158  ++entry_;
159  LogDebug("PFCandidateChecker") << "STOP event: " << iEvent.id().event() << " in run " << iEvent.id().run()
160  << std::endl;
161 }
162 
164  double metX = 0.;
165  double metY = 0.;
166  for (unsigned i = 0; i < pfReco.size(); i++) {
167  metX += pfReco[i].px();
168  metY += pfReco[i].py();
169  }
170  double met = std::sqrt(metX * metX + metY * metY);
171  std::cout << "MET RECO = " << metX << " " << metY << " " << met << std::endl;
172 
173  metX = 0.;
174  metY = 0.;
175  for (unsigned i = 0; i < pfReReco.size(); i++) {
176  metX += pfReReco[i].px();
177  metY += pfReReco[i].py();
178  }
179  met = std::sqrt(metX * metX + metY * metY);
180  std::cout << "MET Re-RECO = " << metX << " " << metY << " " << met << std::endl;
181 }
182 
183 void PFCandidateChecker::printJets(const PFJetCollection& pfJetsReco, const PFJetCollection& pfJetsReReco) const {
184  bool differentSize = pfJetsReco.size() != pfJetsReReco.size();
185  if (differentSize)
186  std::cout << "+++WARNING+++ PFJet size changed for entry " << entry_ << " !" << endl
187  << " - RECO size : " << pfJetsReco.size() << endl
188  << " - Re-RECO size : " << pfJetsReReco.size() << endl;
189  unsigned minSize = pfJetsReco.size() < pfJetsReReco.size() ? pfJetsReco.size() : pfJetsReReco.size();
190  unsigned npr = 0;
191  for (unsigned i = 0; i < minSize; ++i) {
192  const reco::PFJet& candReco = pfJetsReco[i];
193  const reco::PFJet& candReReco = pfJetsReReco[i];
194  if (candReco.et() < 20. && candReReco.et() < 20.)
195  break;
196  double deltaE = (candReReco.et() - candReco.et()) / (candReReco.et() + candReco.et());
197  double deltaEta = candReReco.eta() - candReco.eta();
198  double deltaPhi = candReReco.phi() - candReco.phi();
199  if (fabs(deltaE) > deltaEMax_ || fabs(deltaEta) > deltaEtaMax_ || fabs(deltaPhi) > deltaPhiMax_) {
200  std::cout << "+++WARNING+++ PFJet " << i << " changed for entry " << entry_ << " ! " << std::endl
201  << " - RECO : " << candReco.et() << " " << candReco.eta() << " " << candReco.phi() << std::endl
202  << " - Re-RECO : " << candReReco.et() << " " << candReReco.eta() << " " << candReReco.phi()
203  << std::endl
204  << " DeltaE = : " << deltaE << std::endl
205  << " DeltaEta = : " << deltaEta << std::endl
206  << " DeltaPhi = : " << deltaPhi << std::endl
207  << std::endl;
208  if (++npr == 5)
209  break;
210  } else {
211  std::cout << "Jet " << i << " " << candReco.et() << std::endl;
212  }
213  }
214 }
215 
217  if (!out)
218  return;
219 
220  PFBlockRef firstRef;
221 
222  assert(!cand.elementsInBlocks().empty());
223  for (unsigned i = 0; i < cand.elementsInBlocks().size(); i++) {
224  PFBlockRef blockRef = cand.elementsInBlocks()[i].first;
225 
226  if (blockRef.isNull()) {
227  cerr << "ERROR! no block ref!";
228  continue;
229  }
230 
231  if (!i) {
232  out << (*blockRef);
233  firstRef = blockRef;
234  } else if (blockRef != firstRef) {
235  cerr << "WARNING! This PFCandidate is not made from a single block" << endl;
236  }
237 
238  out << "\t" << cand.elementsInBlocks()[i].second << endl;
239  }
240 }
void analyze(const edm::Event &, const edm::EventSetup &) override
static bool greaterPt(const reco::PFCandidate &a, const reco::PFCandidate &b)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void printMet(const reco::PFCandidateCollection &pfReco, const reco::PFCandidateCollection &pfReReco) const
void printJets(const reco::PFJetCollection &pfJetsReco, const reco::PFJetCollection &pfJetsReReco) const
unsigned entry_
Counter.
edm::EDGetTokenT< reco::PFCandidateCollection > inputTokenPFCandidatesReReco_
edm::EDGetTokenT< reco::PFCandidateCollection > inputTokenPFCandidatesReco_
PFCandidates in which we&#39;ll look for pile up particles.
void printElementsInBlocks(const reco::PFCandidate &cand, std::ostream &out=std::cout) const
double deltaEMax_
Cuts for comparison.
assert(be >=bs)
Jets made from PFObjects.
Definition: PFJet.h:20
static const double deltaEta
Definition: CaloConstants.h:8
T getUntrackedParameter(std::string const &, T const &) const
Checks what a re-reco changes in PFCandidates.
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:23
bool verbose_
verbose ?
edm::EDGetTokenT< reco::PFJetCollection > inputTokenPFJetsReco_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool isNull() const
Checks for null.
Definition: Ref.h:235
bool rankByPt_
rank the candidates by Pt
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
double b
Definition: hdecay.h:120
std::vector< PFJet > PFJetCollection
collection of PFJet objects
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
double et() const final
transverse energy
HLT enums.
double a
Definition: hdecay.h:121
bool printBlocks_
print the blocks associated to a given candidate ?
PFCandidateChecker(const edm::ParameterSet &)
double phi() const final
momentum azimuthal angle
Definition: Run.h:45
edm::EDGetTokenT< reco::PFJetCollection > inputTokenPFJetsReReco_
#define LogDebug(id)
double energy() const final
energy
double eta() const final
momentum pseudorapidity