CMS 3D CMS Logo

HBHEIsolatedNoiseReflagger.cc
Go to the documentation of this file.
1 /*
2 Description: "Reflags" HB/HE hits based on their ECAL, HCAL, and tracking isolation.
3 
4 Original Author: John Paul Chou (Brown University)
5  Thursday, September 2, 2010
6 */
7 
9 
22 
24  :
25 
26  LooseHcalIsol_(iConfig.getParameter<double>("LooseHcalIsol")),
27  LooseEcalIsol_(iConfig.getParameter<double>("LooseEcalIsol")),
28  LooseTrackIsol_(iConfig.getParameter<double>("LooseTrackIsol")),
29  TightHcalIsol_(iConfig.getParameter<double>("TightHcalIsol")),
30  TightEcalIsol_(iConfig.getParameter<double>("TightEcalIsol")),
31  TightTrackIsol_(iConfig.getParameter<double>("TightTrackIsol")),
32 
33  LooseRBXEne1_(iConfig.getParameter<double>("LooseRBXEne1")),
34  LooseRBXEne2_(iConfig.getParameter<double>("LooseRBXEne2")),
35  LooseRBXHits1_(iConfig.getParameter<int>("LooseRBXHits1")),
36  LooseRBXHits2_(iConfig.getParameter<int>("LooseRBXHits2")),
37  TightRBXEne1_(iConfig.getParameter<double>("TightRBXEne1")),
38  TightRBXEne2_(iConfig.getParameter<double>("TightRBXEne2")),
39  TightRBXHits1_(iConfig.getParameter<int>("TightRBXHits1")),
40  TightRBXHits2_(iConfig.getParameter<int>("TightRBXHits2")),
41 
42  LooseHPDEne1_(iConfig.getParameter<double>("LooseHPDEne1")),
43  LooseHPDEne2_(iConfig.getParameter<double>("LooseHPDEne2")),
44  LooseHPDHits1_(iConfig.getParameter<int>("LooseHPDHits1")),
45  LooseHPDHits2_(iConfig.getParameter<int>("LooseHPDHits2")),
46  TightHPDEne1_(iConfig.getParameter<double>("TightHPDEne1")),
47  TightHPDEne2_(iConfig.getParameter<double>("TightHPDEne2")),
48  TightHPDHits1_(iConfig.getParameter<int>("TightHPDHits1")),
49  TightHPDHits2_(iConfig.getParameter<int>("TightHPDHits2")),
50 
51  LooseDiHitEne_(iConfig.getParameter<double>("LooseDiHitEne")),
52  TightDiHitEne_(iConfig.getParameter<double>("TightDiHitEne")),
53  LooseMonoHitEne_(iConfig.getParameter<double>("LooseMonoHitEne")),
54  TightMonoHitEne_(iConfig.getParameter<double>("TightMonoHitEne")),
55 
56  RBXEneThreshold_(iConfig.getParameter<double>("RBXEneThreshold")),
57 
58  debug_(iConfig.getUntrackedParameter<bool>("debug", true)),
59  objvalidator_(iConfig) {
60  tok_hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
61  tok_EB_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ebInput"));
62  tok_EE_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("eeInput"));
64  consumes<std::vector<reco::TrackExtrapolation> >(iConfig.getParameter<edm::InputTag>("trackExtrapolationInput"));
65 
66  // ES tokens
67  ecalChStatusToken_ = esConsumes<EcalChannelStatus, EcalChannelStatusRcd>();
68  hcalChStatusToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
69  hcalSevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
70  ecalSevToken_ = esConsumes<EcalSeverityLevelAlgo, EcalSeverityLevelAlgoRcd>();
71  ctcmToken_ = esConsumes<CaloTowerConstituentsMap, CaloGeometryRecord>();
72  hfemapToken_ = esConsumes<HcalFrontEndMap, HcalFrontEndMapRcd>();
74 
75  produces<HBHERecHitCollection>();
76 }
77 
79 
81  // get the ECAL channel status map
82  const EcalChannelStatus* dbEcalChStatus = &evSetup.getData(ecalChStatusToken_);
83 
84  // get the HCAL channel status map
85  const HcalChannelQuality* dbHcalChStatus = &evSetup.getData(hcalChStatusToken_);
86 
87  // get the severity level computers
88  const HcalSeverityLevelComputer* hcalSevLvlComputer = &evSetup.getData(hcalSevToken_);
89  const EcalSeverityLevelAlgo* ecalSevLvlAlgo = &evSetup.getData(ecalSevToken_);
90 
91  // get the calotower mappings
92  const CaloTowerConstituentsMap& ctcm = evSetup.getData(ctcmToken_);
93 
94  // get hcal frontend map
95  const HcalFrontEndMap* hfemap = &evSetup.getData(hfemapToken_);
96 
97  // get the HB/HE hits
99  iEvent.getByToken(tok_hbhe_, hbhehits_h);
100 
101  // get the ECAL hits
103  iEvent.getByToken(tok_EB_, ebhits_h);
105  iEvent.getByToken(tok_EE_, eehits_h);
106 
107  // get the tracks
109  iEvent.getByToken(tok_trackExt_, trackextraps_h);
110 
111  // set the status maps and severity level computers for the hit validator
112  objvalidator_.setHcalChannelQuality(dbHcalChStatus);
113  objvalidator_.setEcalChannelStatus(dbEcalChStatus);
114  objvalidator_.setHcalSeverityLevelComputer(hcalSevLvlComputer);
115  objvalidator_.setEcalSeverityLevelAlgo(ecalSevLvlAlgo);
116  objvalidator_.setEBRecHitCollection(&(*ebhits_h));
117  objvalidator_.setEERecHitCollection(&(*eehits_h));
118 
119  // organizer the hits
121  hbhehits_h, ebhits_h, eehits_h, trackextraps_h, objvalidator_, ctcm, evSetup.getData(geoToken_));
122  HBHEHitMapOrganizer organizer(hbhehits_h, objvalidator_, pto, hfemap);
123 
124  // get the rbxs, hpds, dihits, and monohits
125  std::vector<HBHEHitMap> rbxs;
126  std::vector<HBHEHitMap> hpds;
127  std::vector<HBHEHitMap> dihits;
128  std::vector<HBHEHitMap> monohits;
129  organizer.getRBXs(rbxs, LooseRBXEne1_ < TightRBXEne1_ ? LooseRBXEne1_ : TightRBXEne1_);
130  organizer.getHPDs(hpds, LooseHPDEne1_ < TightHPDEne1_ ? LooseHPDEne1_ : TightHPDEne1_);
131  organizer.getDiHits(dihits, LooseDiHitEne_ < TightDiHitEne_ ? LooseDiHitEne_ : TightDiHitEne_);
132  organizer.getMonoHits(monohits, LooseMonoHitEne_ < TightMonoHitEne_ ? LooseMonoHitEne_ : TightMonoHitEne_);
133 
134  if (debug_ && (!rbxs.empty() || !hpds.empty() || !dihits.empty() || !monohits.empty())) {
135  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "RBXs:" << std::endl;
136  DumpHBHEHitMap(rbxs);
137  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nHPDs:" << std::endl;
138  DumpHBHEHitMap(hpds);
139  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nDiHits:" << std::endl;
140  DumpHBHEHitMap(dihits);
141  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nMonoHits:" << std::endl;
142  DumpHBHEHitMap(monohits);
143  }
144 
145  // bool result=true;
146 
147  // determine which hits are noisy
148  std::set<const HBHERecHit*> noisehits;
149  for (int i = 0; i < static_cast<int>(rbxs.size()); i++) {
150  int nhits = rbxs[i].nHits();
151  double ene = rbxs[i].hitEnergy();
152  double trkfide = rbxs[i].hitEnergyTrackFiducial();
153  double isolhcale = rbxs[i].hcalEnergySameTowers() + rbxs[i].hcalEnergyNeighborTowers();
154  double isolecale = rbxs[i].ecalEnergySameTowers();
155  double isoltrke = rbxs[i].trackEnergySameTowers() + rbxs[i].trackEnergyNeighborTowers();
156  //
157  // RBX mistag reduction
158  bool isLooseIso = false;
159  bool isTightIso = false;
160  if (ene > RBXEneThreshold_ && ene > 0) { // New absolute iso-cut for high energy RBX clusters
161  if (isolhcale < LooseHcalIsol_ * RBXEneThreshold_ && isolecale < LooseEcalIsol_ * RBXEneThreshold_ &&
162  isoltrke < LooseTrackIsol_ * RBXEneThreshold_)
163  isLooseIso = true;
164  if (isolhcale < TightHcalIsol_ * RBXEneThreshold_ && isolecale < TightEcalIsol_ * RBXEneThreshold_ &&
165  isoltrke < TightTrackIsol_ * RBXEneThreshold_)
166  isTightIso = true;
167  }
168  if (ene <= RBXEneThreshold_ && ene > 0) { // Old relative iso-cut for low energy RBX clusters
169  if (isolhcale / ene < LooseHcalIsol_ && isolecale / ene < LooseEcalIsol_ && isoltrke / ene < LooseTrackIsol_)
170  isLooseIso = true;
171  if (isolhcale / ene < TightHcalIsol_ && isolecale / ene < TightEcalIsol_ && isoltrke / ene < TightTrackIsol_)
172  isTightIso = true;
173  }
174  //
175  if ((isLooseIso && ((trkfide > LooseRBXEne1_ && nhits >= LooseRBXHits1_) ||
176  (trkfide > LooseRBXEne2_ && nhits >= LooseRBXHits2_))) ||
177  (isTightIso && ((trkfide > TightRBXEne1_ && nhits >= TightRBXHits1_) ||
178  (trkfide > TightRBXEne2_ && nhits >= TightRBXHits2_)))) {
179  for (HBHEHitMap::hitmap_const_iterator it = rbxs[i].beginHits(); it != rbxs[i].endHits(); ++it)
180  noisehits.insert(it->first);
181  // result=false;
182  }
183  }
184 
185  for (int i = 0; i < static_cast<int>(hpds.size()); i++) {
186  int nhits = hpds[i].nHits();
187  double ene = hpds[i].hitEnergy();
188  double trkfide = hpds[i].hitEnergyTrackFiducial();
189  double isolhcale = hpds[i].hcalEnergySameTowers() + hpds[i].hcalEnergyNeighborTowers();
190  double isolecale = hpds[i].ecalEnergySameTowers();
191  double isoltrke = hpds[i].trackEnergySameTowers() + hpds[i].trackEnergyNeighborTowers();
192  if ((ene > 0 && isolhcale / ene < LooseHcalIsol_ && isolecale / ene < LooseEcalIsol_ &&
193  isoltrke / ene < LooseTrackIsol_ &&
194  ((trkfide > LooseHPDEne1_ && nhits >= LooseHPDHits1_) ||
195  (trkfide > LooseHPDEne2_ && nhits >= LooseHPDHits2_))) ||
196  (ene > 0 && isolhcale / ene < TightHcalIsol_ && isolecale / ene < TightEcalIsol_ &&
197  isoltrke / ene < TightTrackIsol_ &&
198  ((trkfide > TightHPDEne1_ && nhits >= TightHPDHits1_) ||
199  (trkfide > TightHPDEne2_ && nhits >= TightHPDHits2_)))) {
200  for (HBHEHitMap::hitmap_const_iterator it = hpds[i].beginHits(); it != hpds[i].endHits(); ++it)
201  noisehits.insert(it->first);
202  // result=false;
203  }
204  }
205 
206  for (int i = 0; i < static_cast<int>(dihits.size()); i++) {
207  double ene = dihits[i].hitEnergy();
208  double trkfide = dihits[i].hitEnergyTrackFiducial();
209  double isolhcale = dihits[i].hcalEnergySameTowers() + dihits[i].hcalEnergyNeighborTowers();
210  double isolecale = dihits[i].ecalEnergySameTowers();
211  double isoltrke = dihits[i].trackEnergySameTowers() + dihits[i].trackEnergyNeighborTowers();
212  if ((ene > 0 && isolhcale / ene < LooseHcalIsol_ && isolecale / ene < LooseEcalIsol_ &&
213  isoltrke / ene < LooseTrackIsol_ && trkfide > 0.99 * ene && trkfide > LooseDiHitEne_) ||
214  (ene > 0 && isolhcale / ene < TightHcalIsol_ && isolecale / ene < TightEcalIsol_ &&
215  isoltrke / ene < TightTrackIsol_ && ene > TightDiHitEne_)) {
216  for (HBHEHitMap::hitmap_const_iterator it = dihits[i].beginHits(); it != dihits[i].endHits(); ++it)
217  noisehits.insert(it->first);
218  // result=false;
219  }
220  }
221 
222  for (int i = 0; i < static_cast<int>(monohits.size()); i++) {
223  double ene = monohits[i].hitEnergy();
224  double trkfide = monohits[i].hitEnergyTrackFiducial();
225  double isolhcale = monohits[i].hcalEnergySameTowers() + monohits[i].hcalEnergyNeighborTowers();
226  double isolecale = monohits[i].ecalEnergySameTowers();
227  double isoltrke = monohits[i].trackEnergySameTowers() + monohits[i].trackEnergyNeighborTowers();
228  if ((ene > 0 && isolhcale / ene < LooseHcalIsol_ && isolecale / ene < LooseEcalIsol_ &&
229  isoltrke / ene < LooseTrackIsol_ && trkfide > 0.99 * ene && trkfide > LooseMonoHitEne_) ||
230  (ene > 0 && isolhcale / ene < TightHcalIsol_ && isolecale / ene < TightEcalIsol_ &&
231  isoltrke / ene < TightTrackIsol_ && ene > TightMonoHitEne_)) {
232  for (HBHEHitMap::hitmap_const_iterator it = monohits[i].beginHits(); it != monohits[i].endHits(); ++it)
233  noisehits.insert(it->first);
234  // result=false;
235  }
236  }
237 
238  // prepare the output HBHE RecHit collection
239  auto pOut = std::make_unique<HBHERecHitCollection>();
240  // loop over rechits, and set the new bit you wish to use
241  for (HBHERecHitCollection::const_iterator it = hbhehits_h->begin(); it != hbhehits_h->end(); ++it) {
242  const HBHERecHit* hit = &(*it);
243  HBHERecHit newhit(*hit);
244  if (noisehits.end() != noisehits.find(hit)) {
246  }
247  pOut->push_back(newhit);
248  }
249 
250  iEvent.put(std::move(pOut));
251 
252  return;
253 }
254 
255 void HBHEIsolatedNoiseReflagger::DumpHBHEHitMap(std::vector<HBHEHitMap>& i) const {
256  for (std::vector<HBHEHitMap>::const_iterator it = i.begin(); it != i.end(); ++it) {
257  edm::LogInfo("HBHEIsolatedNoiseReflagger")
258  << "hit energy=" << it->hitEnergy() << "; # hits=" << it->nHits()
259  << "; hcal energy same=" << it->hcalEnergySameTowers() << "; ecal energy same=" << it->ecalEnergySameTowers()
260  << "; track energy same=" << it->trackEnergySameTowers()
261  << "; neighbor hcal energy=" << it->hcalEnergyNeighborTowers() << std::endl;
262  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "hits:" << std::endl;
263  for (HBHEHitMap::hitmap_const_iterator it2 = it->beginHits(); it2 != it->endHits(); ++it2) {
264  const HBHERecHit* hit = it2->first;
265  edm::LogInfo("HBHEIsolatedNoiseReflagger")
266  << "RBX #=" << hfemap->lookupRBX(hit->id()) << "; HPD #=" << hfemap->lookupRMIndex(hit->id()) << "; "
267  << (*hit) << std::endl;
268  }
269  }
270  return;
271 }
272 
273 //define this as a plug-in
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::map< const HBHERecHit *, const PhysicsTower * >::const_iterator hitmap_const_iterator
edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > ecalChStatusToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setHcalChannelQuality(const HcalChannelQuality *q)
edm::ESGetToken< CaloTowerConstituentsMap, CaloGeometryRecord > ctcmToken_
void DumpHBHEHitMap(std::vector< HBHEHitMap > &i) const
std::vector< T >::const_iterator const_iterator
edm::ESGetToken< EcalSeverityLevelAlgo, EcalSeverityLevelAlgoRcd > ecalSevToken_
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
HBHEIsolatedNoiseReflagger(const edm::ParameterSet &)
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > hcalSevToken_
edm::ESGetToken< HcalFrontEndMap, HcalFrontEndMapRcd > hfemapToken_
constexpr void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.h:36
void setHcalSeverityLevelComputer(const HcalSeverityLevelComputer *q)
edm::EDGetTokenT< EcalRecHitCollection > tok_EE_
const int lookupRMIndex(DetId fId) const
int iEvent
Definition: GenABIO.cc:224
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > hcalChStatusToken_
edm::EDGetTokenT< std::vector< reco::TrackExtrapolation > > tok_trackExt_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geoToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const_iterator begin() const
unsigned int id
const_iterator end() const
void setEcalSeverityLevelAlgo(const EcalSeverityLevelAlgo *q)
void setEcalChannelStatus(const EcalChannelStatus *q)
Log< level::Info, false > LogInfo
void setEERecHitCollection(const EcalRecHitCollection *q)
void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< EcalRecHitCollection > tok_EB_
void setEBRecHitCollection(const EcalRecHitCollection *q)
const std::string lookupRBX(DetId fId) const
brief lookup the RBX associated with the given logical id
def move(src, dest)
Definition: eostools.py:511