CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 
22  hbheLabel_(iConfig.getParameter<edm::InputTag>("hbheInput")),
23  ebLabel_(iConfig.getParameter<edm::InputTag>("ebInput")),
24  eeLabel_(iConfig.getParameter<edm::InputTag>("eeInput")),
25  trackExtrapolationLabel_(iConfig.getParameter<edm::InputTag>("trackExtrapolationInput")),
26 
27  LooseHcalIsol_(iConfig.getParameter<double>("LooseHcalIsol")),
28  LooseEcalIsol_(iConfig.getParameter<double>("LooseEcalIsol")),
29  LooseTrackIsol_(iConfig.getParameter<double>("LooseTrackIsol")),
30  TightHcalIsol_(iConfig.getParameter<double>("TightHcalIsol")),
31  TightEcalIsol_(iConfig.getParameter<double>("TightEcalIsol")),
32  TightTrackIsol_(iConfig.getParameter<double>("TightTrackIsol")),
33 
34  LooseRBXEne1_(iConfig.getParameter<double>("LooseRBXEne1")),
35  LooseRBXEne2_(iConfig.getParameter<double>("LooseRBXEne2")),
36  LooseRBXHits1_(iConfig.getParameter<int>("LooseRBXHits1")),
37  LooseRBXHits2_(iConfig.getParameter<int>("LooseRBXHits2")),
38  TightRBXEne1_(iConfig.getParameter<double>("TightRBXEne1")),
39  TightRBXEne2_(iConfig.getParameter<double>("TightRBXEne2")),
40  TightRBXHits1_(iConfig.getParameter<int>("TightRBXHits1")),
41  TightRBXHits2_(iConfig.getParameter<int>("TightRBXHits2")),
42 
43  LooseHPDEne1_(iConfig.getParameter<double>("LooseHPDEne1")),
44  LooseHPDEne2_(iConfig.getParameter<double>("LooseHPDEne2")),
45  LooseHPDHits1_(iConfig.getParameter<int>("LooseHPDHits1")),
46  LooseHPDHits2_(iConfig.getParameter<int>("LooseHPDHits2")),
47  TightHPDEne1_(iConfig.getParameter<double>("TightHPDEne1")),
48  TightHPDEne2_(iConfig.getParameter<double>("TightHPDEne2")),
49  TightHPDHits1_(iConfig.getParameter<int>("TightHPDHits1")),
50  TightHPDHits2_(iConfig.getParameter<int>("TightHPDHits2")),
51 
52  LooseDiHitEne_(iConfig.getParameter<double>("LooseDiHitEne")),
53  TightDiHitEne_(iConfig.getParameter<double>("TightDiHitEne")),
54  LooseMonoHitEne_(iConfig.getParameter<double>("LooseMonoHitEne")),
55  TightMonoHitEne_(iConfig.getParameter<double>("TightMonoHitEne")),
56 
57  debug_(iConfig.getUntrackedParameter<bool>("debug",true)),
58  objvalidator_(iConfig)
59 {
60  produces<HBHERecHitCollection>();
61 }
62 
64 {
65 }
66 
67 
68 void
70 {
71  // get the ECAL channel status map
73  evSetup.get<EcalChannelStatusRcd>().get( ecalChStatus );
74  const EcalChannelStatus* dbEcalChStatus = ecalChStatus.product();
75 
76  // get the HCAL channel status map
78  evSetup.get<HcalChannelQualityRcd>().get( hcalChStatus );
79  const HcalChannelQuality* dbHcalChStatus = hcalChStatus.product();
80 
81  // get the severity level computers
82  edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
83  evSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
84  const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
85  const EcalSeverityLevelAlgo* ecalSevLvlAlgo = new EcalSeverityLevelAlgo();
86 
87  // get the calotower mappings
89  evSetup.get<IdealGeometryRecord>().get(ctcm);
90 
91  // get the HB/HE hits
93  iEvent.getByLabel(hbheLabel_, hbhehits_h);
94 
95  // get the ECAL hits
97  iEvent.getByLabel(ebLabel_, ebhits_h);
99  iEvent.getByLabel(eeLabel_, eehits_h);
100 
101  // get the tracks
103  iEvent.getByLabel(trackExtrapolationLabel_, trackextraps_h);
104 
105  // set the status maps and severity level computers for the hit validator
106  objvalidator_.setHcalChannelQuality(dbHcalChStatus);
107  objvalidator_.setEcalChannelStatus(dbEcalChStatus);
108  objvalidator_.setHcalSeverityLevelComputer(hcalSevLvlComputer);
109  objvalidator_.setEcalSeverityLevelAlgo(ecalSevLvlAlgo);
110  objvalidator_.setEBRecHitCollection(&(*ebhits_h));
111  objvalidator_.setEERecHitCollection(&(*eehits_h));
112 
113  // organizer the hits
114  PhysicsTowerOrganizer pto(iEvent, evSetup, hbhehits_h, ebhits_h, eehits_h, trackextraps_h, objvalidator_, *(ctcm.product()));
115  HBHEHitMapOrganizer organizer(hbhehits_h, objvalidator_, pto);
116 
117  // get the rbxs, hpds, dihits, and monohits
118  std::vector<HBHEHitMap> rbxs;
119  std::vector<HBHEHitMap> hpds;
120  std::vector<HBHEHitMap> dihits;
121  std::vector<HBHEHitMap> monohits;
126 
127  if(debug_ && (rbxs.size()>0 || hpds.size()>0 || dihits.size()>0 || monohits.size()>0)) {
128  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "RBXs:" << std::endl;
129  DumpHBHEHitMap(rbxs);
130  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nHPDs:" << std::endl;
131  DumpHBHEHitMap(hpds);
132  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nDiHits:" << std::endl;
133  DumpHBHEHitMap(dihits);
134  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "\nMonoHits:" << std::endl;
135  DumpHBHEHitMap(monohits);
136  }
137 
138  bool result=true;
139 
140  // determine which hits are noisy
141  std::set<const HBHERecHit*> noisehits;
142  for(int i=0; i<static_cast<int>(rbxs.size()); i++) {
143  int nhits=rbxs[i].nHits();
144  double ene=rbxs[i].hitEnergy();
145  double trkfide=rbxs[i].hitEnergyTrackFiducial();
146  double isolhcale=rbxs[i].hcalEnergySameTowers()+rbxs[i].hcalEnergyNeighborTowers();
147  double isolecale=rbxs[i].ecalEnergySameTowers();
148  double isoltrke=rbxs[i].trackEnergySameTowers()+rbxs[i].trackEnergyNeighborTowers();
149  if((isolhcale/ene<LooseHcalIsol_ && isolecale/ene<LooseEcalIsol_ && isoltrke/ene<LooseTrackIsol_ && ((trkfide>LooseRBXEne1_ && nhits>=LooseRBXHits1_) || (trkfide>LooseRBXEne2_ && nhits>=LooseRBXHits2_))) ||
150  (isolhcale/ene<TightHcalIsol_ && isolecale/ene<TightEcalIsol_ && isoltrke/ene<TightTrackIsol_ && ((trkfide>TightRBXEne1_ && nhits>=TightRBXHits1_) || (trkfide>TightRBXEne2_ && nhits>=TightRBXHits2_)))) {
151  for(HBHEHitMap::hitmap_const_iterator it=rbxs[i].beginHits(); it!=rbxs[i].endHits(); ++it)
152  noisehits.insert(it->first);
153  result=false;
154  }
155  }
156 
157  for(int i=0; i<static_cast<int>(hpds.size()); i++) {
158  int nhits=hpds[i].nHits();
159  double ene=hpds[i].hitEnergy();
160  double trkfide=hpds[i].hitEnergyTrackFiducial();
161  double isolhcale=hpds[i].hcalEnergySameTowers()+hpds[i].hcalEnergyNeighborTowers();
162  double isolecale=hpds[i].ecalEnergySameTowers();
163  double isoltrke=hpds[i].trackEnergySameTowers()+hpds[i].trackEnergyNeighborTowers();
164  if((isolhcale/ene<LooseHcalIsol_ && isolecale/ene<LooseEcalIsol_ && isoltrke/ene<LooseTrackIsol_ && ((trkfide>LooseHPDEne1_ && nhits>=LooseHPDHits1_) || (trkfide>LooseHPDEne2_ && nhits>=LooseHPDHits2_))) ||
165  (isolhcale/ene<TightHcalIsol_ && isolecale/ene<TightEcalIsol_ && isoltrke/ene<TightTrackIsol_ && ((trkfide>TightHPDEne1_ && nhits>=TightHPDHits1_) || (trkfide>TightHPDEne2_ && nhits>=TightHPDHits2_)))) {
166  for(HBHEHitMap::hitmap_const_iterator it=hpds[i].beginHits(); it!=hpds[i].endHits(); ++it)
167  noisehits.insert(it->first);
168  result=false;
169  }
170  }
171 
172  for(int i=0; i<static_cast<int>(dihits.size()); i++) {
173  double ene=dihits[i].hitEnergy();
174  double trkfide=dihits[i].hitEnergyTrackFiducial();
175  double isolhcale=dihits[i].hcalEnergySameTowers()+dihits[i].hcalEnergyNeighborTowers();
176  double isolecale=dihits[i].ecalEnergySameTowers();
177  double isoltrke=dihits[i].trackEnergySameTowers()+dihits[i].trackEnergyNeighborTowers();
178  if((isolhcale/ene<LooseHcalIsol_ && isolecale/ene<LooseEcalIsol_ && isoltrke/ene<LooseTrackIsol_ && trkfide>0.99*ene && trkfide>LooseDiHitEne_) ||
179  (isolhcale/ene<TightHcalIsol_ && isolecale/ene<TightEcalIsol_ && isoltrke/ene<TightTrackIsol_ && ene>TightDiHitEne_)) {
180  for(HBHEHitMap::hitmap_const_iterator it=dihits[i].beginHits(); it!=dihits[i].endHits(); ++it)
181  noisehits.insert(it->first);
182  result=false;
183  }
184  }
185 
186  for(int i=0; i<static_cast<int>(monohits.size()); i++) {
187  double ene=monohits[i].hitEnergy();
188  double trkfide=monohits[i].hitEnergyTrackFiducial();
189  double isolhcale=monohits[i].hcalEnergySameTowers()+monohits[i].hcalEnergyNeighborTowers();
190  double isolecale=monohits[i].ecalEnergySameTowers();
191  double isoltrke=monohits[i].trackEnergySameTowers()+monohits[i].trackEnergyNeighborTowers();
192  if((isolhcale/ene<LooseHcalIsol_ && isolecale/ene<LooseEcalIsol_ && isoltrke/ene<LooseTrackIsol_ && trkfide>0.99*ene && trkfide>LooseMonoHitEne_) ||
193  (isolhcale/ene<TightHcalIsol_ && isolecale/ene<TightEcalIsol_ && isoltrke/ene<TightTrackIsol_ && ene>TightMonoHitEne_)) {
194  for(HBHEHitMap::hitmap_const_iterator it=monohits[i].beginHits(); it!=monohits[i].endHits(); ++it)
195  noisehits.insert(it->first);
196  result=false;
197  }
198  }
199 
200  // prepare the output HBHE RecHit collection
201  std::auto_ptr<HBHERecHitCollection> pOut(new HBHERecHitCollection());
202  // loop over rechits, and set the new bit you wish to use
203  for(HBHERecHitCollection::const_iterator it=hbhehits_h->begin(); it!=hbhehits_h->end(); ++it) {
204  const HBHERecHit* hit=&(*it);
205  HBHERecHit newhit(*hit);
206  if(noisehits.end()!=noisehits.find(hit)) {
208  }
209  pOut->push_back(newhit);
210  }
211 
212  iEvent.put(pOut);
213 
214  delete ecalSevLvlAlgo;
215 
216  return;
217 }
218 
219 
220 void HBHEIsolatedNoiseReflagger::DumpHBHEHitMap(std::vector<HBHEHitMap>& i) const
221 {
222  for(std::vector<HBHEHitMap>::const_iterator it=i.begin(); it!=i.end(); ++it) {
223  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "hit energy=" << it->hitEnergy()
224  << "; # hits=" << it->nHits()
225  << "; hcal energy same=" << it->hcalEnergySameTowers()
226  << "; ecal energy same=" << it->ecalEnergySameTowers()
227  << "; track energy same=" << it->trackEnergySameTowers()
228  << "; neighbor hcal energy=" << it->hcalEnergyNeighborTowers() << std::endl;
229  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "hits:" << std::endl;
230  for(HBHEHitMap::hitmap_const_iterator it2=it->beginHits(); it2!=it->endHits(); ++it2) {
231  const HBHERecHit *hit=it2->first;
232  edm::LogInfo("HBHEIsolatedNoiseReflagger") << "RBX #=" << HcalHPDRBXMap::indexRBX(hit->id())
233  << "; HPD #=" << HcalHPDRBXMap::indexHPD(hit->id())
234  << "; " << (*hit) << std::endl;
235  }
236  }
237  return;
238 }
239 
240 
241 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
void DumpHBHEHitMap(std::vector< HBHEHitMap > &i) const
std::map< const HBHERecHit *, const PhysicsTower * >::const_iterator hitmap_const_iterator
virtual void produce(edm::Event &, const edm::EventSetup &)
static int indexRBX(const HcalDetId &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setHcalChannelQuality(const HcalChannelQuality *q)
void getHPDs(std::vector< HBHEHitMap > &v, double energy) const
HcalDetId id() const
get the id
Definition: HBHERecHit.h:21
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:22
std::vector< T >::const_iterator const_iterator
HBHEIsolatedNoiseReflagger(const edm::ParameterSet &)
void getRBXs(std::vector< HBHEHitMap > &v, double energy) const
void getMonoHits(std::vector< HBHEHitMap > &v, double energy) const
void setHcalSeverityLevelComputer(const HcalSeverityLevelComputer *q)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
tuple result
Definition: query.py:137
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
edm::SortedCollection< HBHERecHit > HBHERecHitCollection
void setEcalSeverityLevelAlgo(const EcalSeverityLevelAlgo *q)
void setEcalChannelStatus(const EcalChannelStatus *q)
void setEERecHitCollection(const EcalRecHitCollection *q)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
void getDiHits(std::vector< HBHEHitMap > &v, double energy) const
static int indexHPD(const HcalDetId &)
void setEBRecHitCollection(const EcalRecHitCollection *q)