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