CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HBHEIsolatedNoiseAlgos.h
Go to the documentation of this file.
1 #ifndef __HBHE_ISOLATED_NOISE_ALGOS_H__
2 #define __HBHE_ISOLATED_NOISE_ALGOS_H__
3 
4 /*
5 Description: Isolation algorithms used to identify anomalous noise in the HB/HE.
6  These algorithms will be used to reflag HB/HE rechits as noise.
7 
8  There are 6 objects defined here:
9  1) ObjectValidatorAbs
10  2) ObjectValidator
11  3) PhysicsTower
12  4) PhysicsTowerOrganizer
13  5) HBHEHitMap
14  6) HBHEHitMapOrganizer
15  See comments below for details.
16 
17 Original Author: John Paul Chou (Brown University)
18  Thursday, September 2, 2010
19 */
20 
21 
22 // system include files
23 #include <memory>
24 #include <string>
25 #include <vector>
26 
27 // user include files
30 
38 
39 #include <vector>
40 #include <set>
41 #include <map>
42 
43 // forward declarations
44 class HBHERecHit;
45 class EcalRecHit;
46 class HcalChannelQuality;
49 
51 //
52 // ObjectValidator
53 //
54 // determines if an ECAL hit, HCAL hit, or track is valid or not.
55 // Note that various objects need to be passed to the validator before use.
56 // ObjectValidatorAbs provides a base class in case you want to change the
57 // hit/track validation algorithms.
59 
61 {
62  public:
64  virtual ~ObjectValidatorAbs() {}
65 
66  virtual bool validHit(const HBHERecHit&) const = 0;
67  virtual bool validHit(const EcalRecHit&) const = 0;
68  virtual bool validTrack(const reco::Track&) const = 0;
69 };
70 
72 {
73  public:
74  explicit ObjectValidator(const edm::ParameterSet&);
76  uint32_t HcalAcceptSeverityLevel, uint32_t EcalAcceptSeverityLevel, bool UseHcalRecoveredHits, bool UseEcalRecoveredHits,
77  double MinValidTrackPt, double MinValidTrackPtBarrel, int MinValidTrackNHits) :
78  HBThreshold_(HBThreshold),
79  HESThreshold_(HESThreshold),
80  HEDThreshold_(HEDThreshold),
81  EBThreshold_(EBThreshold),
82  EEThreshold_(EEThreshold),
83  HcalAcceptSeverityLevel_(HcalAcceptSeverityLevel),
84  EcalAcceptSeverityLevel_(EcalAcceptSeverityLevel),
85  UseHcalRecoveredHits_(UseHcalRecoveredHits),
86  UseEcalRecoveredHits_(UseEcalRecoveredHits),
87  MinValidTrackPt_(MinValidTrackPt),
88  MinValidTrackPtBarrel_(MinValidTrackPtBarrel),
89  MinValidTrackNHits_(MinValidTrackNHits),
96  virtual ~ObjectValidator();
97 
104 
105 
106  bool validHit(const HBHERecHit&) const;
107  bool validHit(const EcalRecHit&) const;
108  bool validTrack(const reco::Track&) const;
109 
110  private:
111 
112  double HBThreshold_; // energy threshold for HB hits
113  double HESThreshold_; // energy threshold for 5-degree (phi) HE hits
114  double HEDThreshold_; // energy threshold for 10-degree (phi) HE hits
115  double EBThreshold_; // energy threshold for EB hits
116  double EEThreshold_; // energy threshold for EE hits
117 
118  uint32_t HcalAcceptSeverityLevel_; // severity level to accept HCAL hits
119  uint32_t EcalAcceptSeverityLevel_; // severity level to accept ECAL hits
120  bool UseHcalRecoveredHits_; // whether or not to use recovered HCAL hits
121  bool UseEcalRecoveredHits_; // whether or not to use recovered HCAL hits
122 
123  double MinValidTrackPt_; // minimum valid track pT
124  double MinValidTrackPtBarrel_; // minimum valid track pT in the Barrel
125  int MinValidTrackNHits_; // minimum number of hits needed for a valid track
126 
127  // for checking the status of ECAL and HCAL channels stored in the DB
130 
131  // calculator of severety level for ECAL and HCAL
134 
135  // needed to determine an ECAL channel's validity
138 };
139 
140 
142 //
143 // PhysicsTower
144 //
145 // basic structure to hold information relevant at the tower level
146 // consists of hcal hits, ecal hits, and tracks
148 
149 struct PhysicsTower {
151  std::set<const HBHERecHit*> hcalhits;
152  std::set<const EcalRecHit*> ecalhits;
153  std::set<const reco::Track*> tracks;
154 };
155 
157 //
158 // PhysicsTowerOrganizer
159 //
160 // Organizers the PhysicsTowers into CaloTower (ieta,iphi) space.
161 // Also provides methods to find a towers nearest neighbors.
162 // For simplicity, merges the ieta=28 and ieta=29 towers together (calls them
163 // tower "28", collectively).
165 
167  public:
168  struct towercmp {
169  bool operator() (const PhysicsTower& lhs, const PhysicsTower& rhs) const {
170  return (lhs.id < rhs.id);
171  }
172  };
173 
175  const edm::EventSetup& evSetup,
176  const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
177  const edm::Handle<EcalRecHitCollection>& ebhitcoll_h,
178  const edm::Handle<EcalRecHitCollection>& eehitcoll_h,
179  const edm::Handle<std::vector<reco::TrackExtrapolation> >& trackextrapcoll_h,
180  const ObjectValidatorAbs& objectvalidator,
181  const CaloTowerConstituentsMap& ctcm);
182 
184 
185  // find a PhysicsTower by some coordinate
186  inline const PhysicsTower* findTower(const CaloTowerDetId& id) const;
187  inline const PhysicsTower* findTower(int ieta, int iphi) const;
188 
189  // get the neighbors +/- 1 in eta-space or +/- 1 in phi-space
190  // (accounts for change in phi-segmentation starting with eta=21)
191  void findNeighbors(const CaloTowerDetId& id, std::set<const PhysicsTower*>& neighbors) const;
192  void findNeighbors(const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors) const;
193  void findNeighbors(int ieta, int iphi, std::set<const PhysicsTower*>& neighbors) const;
194 
195  private:
196  // the non-const, private version of findTower()
198  PhysicsTower* findTower(int ieta, int iphi);
199 
200  void insert_(CaloTowerDetId& id, const HBHERecHit* hit);
201  void insert_(CaloTowerDetId& id, const EcalRecHit* hit);
202  void insert_(CaloTowerDetId& id, const reco::Track* hit);
203 
204  std::set<PhysicsTower, towercmp> towers_;
205 
206 };
207 
209 //
210 // HBHEHitMap
211 //
212 // Collection of HBHERecHits and their associated PhysicsTowers.
213 // Typically, these maps are organized into RBXs, HPDs, dihits, or monohits.
214 // From this one may calculate the hcal, ecal, and track energy associated
215 // with these hits.
216 //
217 // N.B. Many, of the operations can be computationally expensive and should be
218 // used with care.
220 
221 
222 class HBHEHitMap {
223 
224  public:
225 
226  typedef std::map<const HBHERecHit*, const PhysicsTower*>::const_iterator hitmap_const_iterator;
227  typedef std::set<const PhysicsTower*>::const_iterator neighbor_const_iterator;
228 
229  struct twrinfo {
230  double hcalInMap;
231  double hcalOutOfMap;
232  double ecal;
233  double track;
234  };
235 
236  HBHEHitMap();
237  virtual ~HBHEHitMap() {}
238 
239  // energy of the hits in this collection
240  double hitEnergy(void) const;
241 
242  // energy of the hits in a region fiducial to tracks
243  double hitEnergyTrackFiducial(void) const;
244 
245  // number of hits in this collection
246  int nHits(void) const;
247 
248  // same as above, except for the HCAL hits, ECAL hits, and tracks in the same towers as the hits
249  // note that HCAL hits may be present in the same tower, but not be a part of the collection
250  double hcalEnergySameTowers(void) const;
251  double ecalEnergySameTowers(void) const;
252  double trackEnergySameTowers(void) const;
253  int nHcalHitsSameTowers(void) const;
254  int nEcalHitsSameTowers(void) const;
255  int nTracksSameTowers(void) const;
256  void hcalHitsSameTowers(std::set<const HBHERecHit*>& v) const;
257  void ecalHitsSameTowers(std::set<const EcalRecHit*>& v) const;
258  void tracksSameTowers(std::set<const reco::Track*>& v) const;
259 
260  // same as above except for the hits and tracks in the neighboring towers to the hits in the collection
261  double hcalEnergyNeighborTowers(void) const;
262  double ecalEnergyNeighborTowers(void) const;
263  double trackEnergyNeighborTowers(void) const;
264  int nHcalHitsNeighborTowers(void) const;
265  int nEcalHitsNeighborTowers(void) const;
266  int nTracksNeighborTowers(void) const;
267  void hcalHitsNeighborTowers(std::set<const HBHERecHit*>& v) const;
268  void ecalHitsNeighborTowers(std::set<const EcalRecHit*>& v) const;
269  void tracksNeighborTowers(std::set<const reco::Track*>& v) const;
270 
271  // gives the total energy due to hcal, ecal, and tracks tower-by-tower
272  // separates the hcal energy into that which is a hit in the collection, and that which is not
273  void byTowers(std::vector<twrinfo>& v) const;
274 
275  // find a hit in the hitmap
276  hitmap_const_iterator findHit(const HBHERecHit* hit) const { return hits_.find(hit); }
277 
278  // find a neighbor
279  neighbor_const_iterator findNeighbor(const PhysicsTower* twr) const { return neighbors_.find(twr); }
280 
281  // add a hit to the hitmap
282  void insert(const HBHERecHit* hit, const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors);
283 
284  // access to the private maps and sets
285  inline hitmap_const_iterator beginHits(void) const { return hits_.begin(); }
286  inline hitmap_const_iterator endHits(void) const { return hits_.end(); }
287 
288  inline neighbor_const_iterator beginNeighbors(void) const { return neighbors_.begin(); }
289  inline neighbor_const_iterator endNeighbors(void) const { return neighbors_.end(); }
290 
291  private:
292  std::map<const HBHERecHit*, const PhysicsTower*> hits_;
293  std::set<const PhysicsTower*> neighbors_;
294 
295  void calcHits_(void) const;
296  mutable double hitEnergy_;
297  mutable double hitEnergyTrkFid_;
298  mutable int nHits_;
299 
300  void calcHcalSameTowers_(void) const;
301  mutable double hcalEnergySameTowers_;
302  mutable int nHcalHitsSameTowers_;
303 
304  void calcEcalSameTowers_(void) const;
305  mutable double ecalEnergySameTowers_;
306  mutable int nEcalHitsSameTowers_;
307 
308  void calcTracksSameTowers_(void) const;
309  mutable double trackEnergySameTowers_;
310  mutable int nTracksSameTowers_;
311 
312  void calcHcalNeighborTowers_(void) const;
315 
316  void calcEcalNeighborTowers_(void) const;
319 
320  void calcTracksNeighborTowers_(void) const;
323 
324 };
325 
327 //
328 // HBHEHitMapOrganizer
329 //
330 // Organizers the HBHEHitMaps into RBXs, HPDs, dihits, and monohits
332 
334 {
335  public:
337  const ObjectValidatorAbs& objvalidator,
338  const PhysicsTowerOrganizer& pto,
339  const HcalFrontEndMap* hfemap);
340 
341  virtual ~HBHEHitMapOrganizer() {}
342 
343  void getRBXs(std::vector<HBHEHitMap>& v, double energy) const;
344  void getHPDs(std::vector<HBHEHitMap>& v, double energy) const;
345  void getDiHits(std::vector<HBHEHitMap>& v, double energy) const;
346  void getMonoHits(std::vector<HBHEHitMap>& v, double energy) const;
347 
348  private:
349 
351  std::map<int, HBHEHitMap> rbxs_, hpds_;
352  std::vector<HBHEHitMap> dihits_, monohits_;
353 
354  // helper functions
355  // finds all of the hits which are neighbors and in the same HPD as the reference hit
356  void getHPDNeighbors(const HBHERecHit* hit, std::vector<const HBHERecHit*>& neighbors, const PhysicsTowerOrganizer& pto);
357 
358 };
359 
360 
361 #endif
void getHPDNeighbors(const HBHERecHit *hit, std::vector< const HBHERecHit * > &neighbors, const PhysicsTowerOrganizer &pto)
double hcalEnergySameTowers(void) const
std::set< const HBHERecHit * > hcalhits
void tracksNeighborTowers(std::set< const reco::Track * > &v) const
int nHcalHitsSameTowers(void) const
std::map< const HBHERecHit *, const PhysicsTower * >::const_iterator hitmap_const_iterator
std::vector< HBHEHitMap > dihits_
hitmap_const_iterator findHit(const HBHERecHit *hit) const
tuple HcalAcceptSeverityLevel
void setHcalChannelQuality(const HcalChannelQuality *q)
void getHPDs(std::vector< HBHEHitMap > &v, double energy) const
HBHEHitMapOrganizer(const edm::Handle< HBHERecHitCollection > &hbhehitcoll_h, const ObjectValidatorAbs &objvalidator, const PhysicsTowerOrganizer &pto, const HcalFrontEndMap *hfemap)
neighbor_const_iterator findNeighbor(const PhysicsTower *twr) const
virtual bool validTrack(const reco::Track &) const =0
int nEcalHitsNeighborTowers(void) const
bool validTrack(const reco::Track &) const
neighbor_const_iterator beginNeighbors(void) const
double hcalEnergyNeighborTowers_
void getRBXs(std::vector< HBHEHitMap > &v, double energy) const
ObjectValidator(const edm::ParameterSet &)
void ecalHitsNeighborTowers(std::set< const EcalRecHit * > &v) const
hitmap_const_iterator endHits(void) const
tuple HESThreshold
std::map< const HBHERecHit *, const PhysicsTower * > hits_
void getMonoHits(std::vector< HBHEHitMap > &v, double energy) const
void ecalHitsSameTowers(std::set< const EcalRecHit * > &v) const
void setHcalSeverityLevelComputer(const HcalSeverityLevelComputer *q)
double trackEnergyNeighborTowers_
PhysicsTowerOrganizer(const edm::Event &iEvent, const edm::EventSetup &evSetup, const edm::Handle< HBHERecHitCollection > &hbhehitcoll_h, const edm::Handle< EcalRecHitCollection > &ebhitcoll_h, const edm::Handle< EcalRecHitCollection > &eehitcoll_h, const edm::Handle< std::vector< reco::TrackExtrapolation > > &trackextrapcoll_h, const ObjectValidatorAbs &objectvalidator, const CaloTowerConstituentsMap &ctcm)
void calcTracksSameTowers_(void) const
double ecalEnergyNeighborTowers(void) const
std::set< const PhysicsTower * > neighbors_
void insert(const HBHERecHit *hit, const PhysicsTower *twr, std::set< const PhysicsTower * > &neighbors)
const PhysicsTower * findTower(const CaloTowerDetId &id) const
void tracksSameTowers(std::set< const reco::Track * > &v) const
int iEvent
Definition: GenABIO.cc:230
tuple UseHcalRecoveredHits
int nHcalHitsNeighborTowers(void) const
void calcHcalSameTowers_(void) const
std::map< int, HBHEHitMap > rbxs_
const EcalChannelStatus * theEcalChStatus_
const EcalRecHitCollection * theEERecHitCollection_
double hcalEnergyNeighborTowers(void) const
double trackEnergySameTowers(void) const
void insert_(CaloTowerDetId &id, const HBHERecHit *hit)
tuple HEDThreshold
const EcalRecHitCollection * theEBRecHitCollection_
ObjectValidator(double HBThreshold, double HESThreshold, double HEDThreshold, double EBThreshold, double EEThreshold, uint32_t HcalAcceptSeverityLevel, uint32_t EcalAcceptSeverityLevel, bool UseHcalRecoveredHits, bool UseEcalRecoveredHits, double MinValidTrackPt, double MinValidTrackPtBarrel, int MinValidTrackNHits)
std::set< const PhysicsTower * >::const_iterator neighbor_const_iterator
std::map< int, HBHEHitMap > hpds_
void setEcalSeverityLevelAlgo(const EcalSeverityLevelAlgo *q)
std::set< PhysicsTower, towercmp > towers_
void setEcalChannelStatus(const EcalChannelStatus *q)
double hitEnergy(void) const
void calcEcalSameTowers_(void) const
virtual bool validHit(const HBHERecHit &) const =0
double ecalEnergyNeighborTowers_
void setEERecHitCollection(const EcalRecHitCollection *q)
void calcTracksNeighborTowers_(void) const
void calcHcalNeighborTowers_(void) const
const EcalSeverityLevelAlgo * theEcalSevLvlAlgo_
std::set< const EcalRecHit * > ecalhits
int nHits(void) const
double trackEnergyNeighborTowers(void) const
bool operator()(const PhysicsTower &lhs, const PhysicsTower &rhs) const
double hitEnergyTrackFiducial(void) const
void calcHits_(void) const
std::set< const reco::Track * > tracks
void calcEcalNeighborTowers_(void) const
std::vector< HBHEHitMap > monohits_
int nTracksSameTowers(void) const
void getDiHits(std::vector< HBHEHitMap > &v, double energy) const
double ecalEnergySameTowers(void) const
const HcalChannelQuality * theHcalChStatus_
int nEcalHitsSameTowers(void) const
int nTracksNeighborTowers(void) const
const HcalFrontEndMap * hfemap_
void hcalHitsSameTowers(std::set< const HBHERecHit * > &v) const
bool validHit(const HBHERecHit &) const
void setEBRecHitCollection(const EcalRecHitCollection *q)
const HcalSeverityLevelComputer * theHcalSevLvlComputer_
void byTowers(std::vector< twrinfo > &v) const
neighbor_const_iterator endNeighbors(void) const
void findNeighbors(const CaloTowerDetId &id, std::set< const PhysicsTower * > &neighbors) const
hitmap_const_iterator beginHits(void) const
void hcalHitsNeighborTowers(std::set< const HBHERecHit * > &v) const
tuple UseEcalRecoveredHits