CMS 3D CMS Logo

PFCandidateRecalibrator.cc
Go to the documentation of this file.
1 // - a PF candidate collection (which uses bugged HCAL respcorrs)
3 // - respCorrs values from fixed GT and bugged GT
4 // Produce:
5 // - a new PFCandidate collection containing the recalibrated PFCandidates in HF and where the neutral had pointing to problematic cells are removed
6 // - a second PFCandidate collection with just those discarded hadrons
7 // - a ValueMap<reco::PFCandidateRef> that maps the old to the new, and vice-versa
8 
17 
19 
28 
31 
34 #include <iostream>
35 
36 struct HEChannel {
37  float eta;
38  float phi;
39  float ratio;
40  HEChannel(float eta, float phi, float ratio) : eta(eta), phi(phi), ratio(ratio) {}
41 };
42 struct HFChannel {
43  int ieta;
44  int iphi;
45  int depth;
46  float ratio;
47  HFChannel(int ieta, int iphi, int depth, float ratio) : ieta(ieta), iphi(iphi), depth(depth), ratio(ratio) {}
48 };
49 
51 public:
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56 
57 private:
58  void beginRun(const edm::Run& iRun, edm::EventSetup const& iSetup) override;
59  void endRun(const edm::Run& iRun, edm::EventSetup const& iSetup) override;
60  void produce(edm::Event&, const edm::EventSetup&) override;
61 
64 
66 
72 
73  std::vector<HEChannel> badChHE_;
74  std::vector<HFChannel> badChHF_;
75 
78 };
79 
81  : pfcandidates_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfcandidates"))),
82  gtCondToken_(esConsumes<edm::Transition::BeginRun>()),
83  htopoToken_(esConsumes<edm::Transition::BeginRun>()),
84  buggedCondToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", "bugged"))),
85  calogeomTokenRun_(esConsumes<edm::Transition::BeginRun>()),
86  calogeomTokenEvent_(esConsumes()),
87  shortFibreThr_(iConfig.getParameter<double>("shortFibreThr")),
88  longFibreThr_(iConfig.getParameter<double>("longFibreThr")) {
89  produces<reco::PFCandidateCollection>();
90  produces<reco::PFCandidateCollection>("discarded");
91  produces<edm::ValueMap<reco::PFCandidateRef>>();
92 }
93 
95  auto newDb = hcalDbWatcher_.check(iSetup);
96  auto newRC = hcalRCWatcher_.check(iSetup);
97  if (newDb || newRC) {
98  //Get Calib Constants from current GT
99  HcalDbService const& gtCond = iSetup.getData(gtCondToken_);
100 
101  //Get Calib Constants from bugged tag
102  const HcalTopology& theHBHETopology = iSetup.getData(htopoToken_);
103 
104  HcalRespCorrs buggedRespCorrs = iSetup.getData(buggedCondToken_);
105  buggedRespCorrs.setTopo(&theHBHETopology);
106 
107  //access calogeometry
108  const CaloGeometry& cgeo = iSetup.getData(calogeomTokenRun_);
109  const HcalGeometry* hgeom = static_cast<const HcalGeometry*>(cgeo.getSubdetectorGeometry(DetId::Hcal, HcalForward));
110 
111  //reset the bad channel containers
112  badChHE_.clear();
113  badChHF_.clear();
114 
115  //fill bad cells HE (use eta, phi)
116  const std::vector<DetId>& cellsHE = hgeom->getValidDetIds(DetId::Detector::Hcal, HcalEndcap);
117  for (auto id : cellsHE) {
118  float currentRespCorr = gtCond.getHcalRespCorr(id)->getValue();
119  float buggedRespCorr = buggedRespCorrs.getValues(id)->getValue();
120  if (buggedRespCorr == 0.)
121  continue;
122 
123  float ratio = currentRespCorr / buggedRespCorr;
124  if (std::abs(ratio - 1.f) > 0.001) {
125  GlobalPoint pos = hgeom->getPosition(id);
126  badChHE_.push_back(HEChannel(pos.eta(), pos.phi(), ratio));
127  }
128  }
129 
130  //fill bad cells HF (use ieta, iphi)
131  auto const& cellsHF = hgeom->getValidDetIds(DetId::Detector::Hcal, HcalForward);
132  for (auto id : cellsHF) {
133  float currentRespCorr = gtCond.getHcalRespCorr(id)->getValue();
134  float buggedRespCorr = buggedRespCorrs.getValues(id)->getValue();
135  if (buggedRespCorr == 0.)
136  continue;
137 
138  float ratio = currentRespCorr / buggedRespCorr;
139  if (std::abs(ratio - 1.f) > 0.001) {
140  HcalDetId dummyId(id);
141  badChHF_.push_back(HFChannel(dummyId.ieta(), dummyId.iphi(), dummyId.depth(), ratio));
142  }
143  }
144  }
145 }
146 
147 void PFCandidateRecalibrator::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {}
148 
150  //access calogeometry
151  const CaloGeometry& cgeo = iSetup.getData(calogeomTokenEvent_);
152  const HcalGeometry* hgeom = static_cast<const HcalGeometry*>(cgeo.getSubdetectorGeometry(DetId::Hcal, HcalForward));
153 
154  //access PFCandidates
156  iEvent.getByToken(pfcandidates_, pfcandidates);
157 
158  int nPfCand = pfcandidates->size();
159  std::unique_ptr<reco::PFCandidateCollection> copy(new reco::PFCandidateCollection());
160  std::unique_ptr<reco::PFCandidateCollection> discarded(new reco::PFCandidateCollection());
161  copy->reserve(nPfCand);
162  std::vector<int> oldToNew(nPfCand), newToOld, badToOld;
163  newToOld.reserve(nPfCand);
164 
165  LogDebug("PFCandidateRecalibrator") << "NEW EV:";
166 
167  //loop over PFCandidates
168  int i = -1;
169  for (const reco::PFCandidate& pf : *pfcandidates) {
170  ++i;
171  float absEta = std::abs(pf.eta());
172 
173  //deal with HE
174  if (pf.particleId() == reco::PFCandidate::ParticleType::h0 &&
175  !badChHE_.empty() && //don't touch if no miscalibration is found
176  absEta > 1.4 && absEta < 3.) {
177  bool toKill = false;
178  for (auto const& badIt : badChHE_)
179  if (reco::deltaR2(pf.eta(), pf.phi(), badIt.eta, badIt.phi) < 0.07)
180  toKill = true;
181 
182  if (toKill) {
183  discarded->push_back(pf);
184  oldToNew[i] = (-discarded->size());
185  badToOld.push_back(i);
186  continue;
187  } else {
188  copy->push_back(pf);
189  oldToNew[i] = (copy->size());
190  newToOld.push_back(i);
191  }
192  }
193  //deal with HF
194  else if ((pf.particleId() == reco::PFCandidate::ParticleType::h_HF ||
195  pf.particleId() == reco::PFCandidate::ParticleType::egamma_HF) &&
196  !badChHF_.empty() && //don't touch if no miscalibration is found
197  absEta >= 3.) {
198  const math::XYZPointF& ecalPoint = pf.positionAtECALEntrance();
199  GlobalPoint ecalGPoint(ecalPoint.X(), ecalPoint.Y(), ecalPoint.Z());
200  HcalDetId closestDetId(hgeom->getClosestCell(ecalGPoint));
201 
202  if (closestDetId.subdet() == HcalForward) {
203  HcalDetId hDetId(closestDetId.subdet(), closestDetId.ieta(), closestDetId.iphi(), 1); //depth1
204 
205  //raw*calEnergy() is the same as *calEnergy() - no corrections are done for HF
206  float longE = pf.rawEcalEnergy() + pf.rawHcalEnergy() / 2.; //depth1
207  float shortE = pf.rawHcalEnergy() / 2.; //depth2
208 
209  float ecalEnergy = pf.rawEcalEnergy();
210  float hcalEnergy = pf.rawHcalEnergy();
211  float totEnergy = ecalEnergy + hcalEnergy;
212  float totEnergyOrig = totEnergy;
213 
214  bool toKill = false;
215 
216  for (auto const& badIt : badChHF_) {
217  if (hDetId.ieta() == badIt.ieta && hDetId.iphi() == badIt.iphi) {
218  LogDebug("PFCandidateRecalibrator")
219  << "==> orig en (tot,H,E): " << pf.energy() << " " << pf.rawHcalEnergy() << " " << pf.rawEcalEnergy();
220  if (badIt.depth == 1) //depth1
221  {
222  longE *= badIt.ratio;
223  ecalEnergy = ((longE - shortE) > 0.) ? (longE - shortE) : 0.;
224  totEnergy = ecalEnergy + hcalEnergy;
225  } else //depth2
226  {
227  shortE *= badIt.ratio;
228  hcalEnergy = 2 * shortE;
229  ecalEnergy = ((longE - shortE) > 0.) ? (longE - shortE) : 0.;
230  totEnergy = ecalEnergy + hcalEnergy;
231  }
232  //kill candidate if goes below thr
233  if ((pf.particleId() == reco::PFCandidate::ParticleType::h_HF && shortE < shortFibreThr_) ||
234  (pf.particleId() == reco::PFCandidate::ParticleType::egamma_HF && longE < longFibreThr_))
235  toKill = true;
236 
237  LogDebug("PFCandidateRecalibrator") << "====> ieta,iphi,depth: " << badIt.ieta << " " << badIt.iphi << " "
238  << badIt.depth << " corr: " << badIt.ratio;
239  LogDebug("PFCandidateRecalibrator")
240  << "====> recal en (tot,H,E): " << totEnergy << " " << hcalEnergy << " " << ecalEnergy;
241  }
242  }
243 
244  if (toKill) {
245  discarded->push_back(pf);
246  oldToNew[i] = (-discarded->size());
247  badToOld.push_back(i);
248 
249  LogDebug("PFCandidateRecalibrator") << "==> KILLED ";
250  } else {
251  copy->push_back(pf);
252  oldToNew[i] = (copy->size());
253  newToOld.push_back(i);
254 
255  copy->back().setHcalEnergy(hcalEnergy, hcalEnergy);
256  copy->back().setEcalEnergy(ecalEnergy, ecalEnergy);
257 
258  float scalingFactor = totEnergy / totEnergyOrig;
259  math::XYZTLorentzVector recalibP4 = pf.p4() * scalingFactor;
260  copy->back().setP4(recalibP4);
261 
262  LogDebug("PFCandidateRecalibrator") << "====> stored en (tot,H,E): " << copy->back().energy() << " "
263  << copy->back().hcalEnergy() << " " << copy->back().ecalEnergy();
264  }
265  } else {
266  copy->push_back(pf);
267  oldToNew[i] = (copy->size());
268  newToOld.push_back(i);
269  }
270  } else {
271  copy->push_back(pf);
272  oldToNew[i] = (copy->size());
273  newToOld.push_back(i);
274  }
275  }
276 
277  // Now we put things in the event
279  edm::OrphanHandle<reco::PFCandidateCollection> badpf = iEvent.put(std::move(discarded), "discarded");
280 
281  std::unique_ptr<edm::ValueMap<reco::PFCandidateRef>> pf2pf(new edm::ValueMap<reco::PFCandidateRef>());
283  std::vector<reco::PFCandidateRef> refs;
284  refs.reserve(nPfCand);
285 
286  // old to new
287  for (auto iOldToNew : oldToNew) {
288  if (iOldToNew > 0) {
289  refs.push_back(reco::PFCandidateRef(newpf, iOldToNew - 1));
290  } else {
291  refs.push_back(reco::PFCandidateRef(badpf, -iOldToNew - 1));
292  }
293  }
294  filler.insert(pfcandidates, refs.begin(), refs.end());
295  // new good to old
296  refs.clear();
297  for (int i : newToOld) {
298  refs.push_back(reco::PFCandidateRef(pfcandidates, i));
299  }
300  filler.insert(newpf, refs.begin(), refs.end());
301  // new bad to old
302  refs.clear();
303  for (int i : badToOld) {
304  refs.push_back(reco::PFCandidateRef(pfcandidates, i));
305  }
306  filler.insert(badpf, refs.begin(), refs.end());
307  // done
308  filler.fill();
309  iEvent.put(std::move(pf2pf));
310 }
311 
314 
315  desc.add<edm::InputTag>("pfcandidates", edm::InputTag("particleFlow"));
316  desc.add<double>("shortFibreThr", 1.4);
317  desc.add<double>("longFibreThr", 1.4);
318 
319  descriptions.add("pfCandidateRecalibrator", desc);
320 }
321 
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > calogeomTokenRun_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
DetId getClosestCell(const GlobalPoint &r) const override
edm::ESWatcher< HcalRecNumberingRecord > hcalDbWatcher_
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:82
void endRun(const edm::Run &iRun, edm::EventSetup const &iSetup) override
std::vector< HEChannel > badChHE_
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
PFCandidateRecalibrator(const edm::ParameterSet &)
const Item * getValues(DetId fId, bool throwOnFail=true) const
HEChannel(float eta, float phi, float ratio)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const edm::ESGetToken< HcalRespCorrs, HcalRespCorrsRcd > buggedCondToken_
int iEvent
Definition: GenABIO.cc:224
const HcalRespCorr * getHcalRespCorr(const HcalGenericDetId &fId) const
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > calogeomTokenEvent_
edm::ESWatcher< HcalRespCorrsRcd > hcalRCWatcher_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
float getValue() const
Definition: HcalRespCorr.h:19
std::vector< HFChannel > badChHF_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HFChannel(int ieta, int iphi, int depth, float ratio)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
edm::EDGetTokenT< reco::PFCandidateCollection > pfcandidates_
const edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
void beginRun(const edm::Run &iRun, edm::EventSetup const &iSetup) override
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
HLT enums.
GlobalPoint getPosition(const DetId &id) const
const edm::ESGetToken< HcalDbService, HcalDbRecord > gtCondToken_
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
def move(src, dest)
Definition: eostools.py:511
void setTopo(const HcalTopology *topo)
Definition: Run.h:45
#define LogDebug(id)
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164