test
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 
37 
38 #include <vector>
39 #include <set>
40 #include <map>
41 
42 // forward declarations
43 class HBHERecHit;
44 class EcalRecHit;
45 class HcalChannelQuality;
48 
50 //
51 // ObjectValidator
52 //
53 // determines if an ECAL hit, HCAL hit, or track is valid or not.
54 // Note that various objects need to be passed to the validator before use.
55 // ObjectValidatorAbs provides a base class in case you want to change the
56 // hit/track validation algorithms.
58 
60 {
61  public:
63  virtual ~ObjectValidatorAbs() {}
64 
65  virtual bool validHit(const HBHERecHit&) const = 0;
66  virtual bool validHit(const EcalRecHit&) const = 0;
67  virtual bool validTrack(const reco::Track&) const = 0;
68 };
69 
71 {
72  public:
73  explicit ObjectValidator(const edm::ParameterSet&);
75  uint32_t HcalAcceptSeverityLevel, uint32_t EcalAcceptSeverityLevel, bool UseHcalRecoveredHits, bool UseEcalRecoveredHits,
76  double MinValidTrackPt, double MinValidTrackPtBarrel, int MinValidTrackNHits) :
77  HBThreshold_(HBThreshold),
78  HESThreshold_(HESThreshold),
79  HEDThreshold_(HEDThreshold),
80  EBThreshold_(EBThreshold),
81  EEThreshold_(EEThreshold),
82  HcalAcceptSeverityLevel_(HcalAcceptSeverityLevel),
83  EcalAcceptSeverityLevel_(EcalAcceptSeverityLevel),
84  UseHcalRecoveredHits_(UseHcalRecoveredHits),
85  UseEcalRecoveredHits_(UseEcalRecoveredHits),
86  MinValidTrackPt_(MinValidTrackPt),
87  MinValidTrackPtBarrel_(MinValidTrackPtBarrel),
88  MinValidTrackNHits_(MinValidTrackNHits),
95  virtual ~ObjectValidator();
96 
103 
104 
105  bool validHit(const HBHERecHit&) const;
106  bool validHit(const EcalRecHit&) const;
107  bool validTrack(const reco::Track&) const;
108 
109  private:
110 
111  double HBThreshold_; // energy threshold for HB hits
112  double HESThreshold_; // energy threshold for 5-degree (phi) HE hits
113  double HEDThreshold_; // energy threshold for 10-degree (phi) HE hits
114  double EBThreshold_; // energy threshold for EB hits
115  double EEThreshold_; // energy threshold for EE hits
116 
117  uint32_t HcalAcceptSeverityLevel_; // severity level to accept HCAL hits
118  uint32_t EcalAcceptSeverityLevel_; // severity level to accept ECAL hits
119  bool UseHcalRecoveredHits_; // whether or not to use recovered HCAL hits
120  bool UseEcalRecoveredHits_; // whether or not to use recovered HCAL hits
121 
122  double MinValidTrackPt_; // minimum valid track pT
123  double MinValidTrackPtBarrel_; // minimum valid track pT in the Barrel
124  int MinValidTrackNHits_; // minimum number of hits needed for a valid track
125 
126  // for checking the status of ECAL and HCAL channels stored in the DB
129 
130  // calculator of severety level for ECAL and HCAL
133 
134  // needed to determine an ECAL channel's validity
137 };
138 
139 
141 //
142 // PhysicsTower
143 //
144 // basic structure to hold information relevant at the tower level
145 // consists of hcal hits, ecal hits, and tracks
147 
148 struct PhysicsTower {
150  std::set<const HBHERecHit*> hcalhits;
151  std::set<const EcalRecHit*> ecalhits;
152  std::set<const reco::Track*> tracks;
153 };
154 
156 //
157 // PhysicsTowerOrganizer
158 //
159 // Organizers the PhysicsTowers into CaloTower (ieta,iphi) space.
160 // Also provides methods to find a towers nearest neighbors.
161 // For simplicity, merges the ieta=28 and ieta=29 towers together (calls them
162 // tower "28", collectively).
164 
166  public:
167  struct towercmp {
168  bool operator() (const PhysicsTower& lhs, const PhysicsTower& rhs) const {
169  return (lhs.id < rhs.id);
170  }
171  };
172 
174  const edm::EventSetup& evSetup,
175  const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
176  const edm::Handle<EcalRecHitCollection>& ebhitcoll_h,
177  const edm::Handle<EcalRecHitCollection>& eehitcoll_h,
178  const edm::Handle<std::vector<reco::TrackExtrapolation> >& trackextrapcoll_h,
179  const ObjectValidatorAbs& objectvalidator,
180  const CaloTowerConstituentsMap& ctcm);
181 
183 
184  // find a PhysicsTower by some coordinate
185  inline const PhysicsTower* findTower(const CaloTowerDetId& id) const;
186  inline const PhysicsTower* findTower(int ieta, int iphi) const;
187 
188  // get the neighbors +/- 1 in eta-space or +/- 1 in phi-space
189  // (accounts for change in phi-segmentation starting with eta=21)
190  void findNeighbors(const CaloTowerDetId& id, std::set<const PhysicsTower*>& neighbors) const;
191  void findNeighbors(const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors) const;
192  void findNeighbors(int ieta, int iphi, std::set<const PhysicsTower*>& neighbors) const;
193 
194  private:
195  // the non-const, private version of findTower()
197  PhysicsTower* findTower(int ieta, int iphi);
198 
199  void insert_(CaloTowerDetId& id, const HBHERecHit* hit);
200  void insert_(CaloTowerDetId& id, const EcalRecHit* hit);
201  void insert_(CaloTowerDetId& id, const reco::Track* hit);
202 
203  std::set<PhysicsTower, towercmp> towers_;
204 
205 };
206 
208 //
209 // HBHEHitMap
210 //
211 // Collection of HBHERecHits and their associated PhysicsTowers.
212 // Typically, these maps are organized into RBXs, HPDs, dihits, or monohits.
213 // From this one may calculate the hcal, ecal, and track energy associated
214 // with these hits.
215 //
216 // N.B. Many, of the operations can be computationally expensive and should be
217 // used with care.
219 
220 
221 class HBHEHitMap {
222 
223  public:
224 
225  typedef std::map<const HBHERecHit*, const PhysicsTower*>::const_iterator hitmap_const_iterator;
226  typedef std::set<const PhysicsTower*>::const_iterator neighbor_const_iterator;
227 
228  struct twrinfo {
229  double hcalInMap;
230  double hcalOutOfMap;
231  double ecal;
232  double track;
233  };
234 
235  HBHEHitMap();
236  virtual ~HBHEHitMap() {}
237 
238  // energy of the hits in this collection
239  double hitEnergy(void) const;
240 
241  // energy of the hits in a region fiducial to tracks
242  double hitEnergyTrackFiducial(void) const;
243 
244  // number of hits in this collection
245  int nHits(void) const;
246 
247  // same as above, except for the HCAL hits, ECAL hits, and tracks in the same towers as the hits
248  // note that HCAL hits may be present in the same tower, but not be a part of the collection
249  double hcalEnergySameTowers(void) const;
250  double ecalEnergySameTowers(void) const;
251  double trackEnergySameTowers(void) const;
252  int nHcalHitsSameTowers(void) const;
253  int nEcalHitsSameTowers(void) const;
254  int nTracksSameTowers(void) const;
255  void hcalHitsSameTowers(std::set<const HBHERecHit*>& v) const;
256  void ecalHitsSameTowers(std::set<const EcalRecHit*>& v) const;
257  void tracksSameTowers(std::set<const reco::Track*>& v) const;
258 
259  // same as above except for the hits and tracks in the neighboring towers to the hits in the collection
260  double hcalEnergyNeighborTowers(void) const;
261  double ecalEnergyNeighborTowers(void) const;
262  double trackEnergyNeighborTowers(void) const;
263  int nHcalHitsNeighborTowers(void) const;
264  int nEcalHitsNeighborTowers(void) const;
265  int nTracksNeighborTowers(void) const;
266  void hcalHitsNeighborTowers(std::set<const HBHERecHit*>& v) const;
267  void ecalHitsNeighborTowers(std::set<const EcalRecHit*>& v) const;
268  void tracksNeighborTowers(std::set<const reco::Track*>& v) const;
269 
270  // gives the total energy due to hcal, ecal, and tracks tower-by-tower
271  // separates the hcal energy into that which is a hit in the collection, and that which is not
272  void byTowers(std::vector<twrinfo>& v) const;
273 
274  // find a hit in the hitmap
275  hitmap_const_iterator findHit(const HBHERecHit* hit) const { return hits_.find(hit); }
276 
277  // find a neighbor
278  neighbor_const_iterator findNeighbor(const PhysicsTower* twr) const { return neighbors_.find(twr); }
279 
280  // add a hit to the hitmap
281  void insert(const HBHERecHit* hit, const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors);
282 
283  // access to the private maps and sets
284  inline hitmap_const_iterator beginHits(void) const { return hits_.begin(); }
285  inline hitmap_const_iterator endHits(void) const { return hits_.end(); }
286 
287  inline neighbor_const_iterator beginNeighbors(void) const { return neighbors_.begin(); }
288  inline neighbor_const_iterator endNeighbors(void) const { return neighbors_.end(); }
289 
290  private:
291  std::map<const HBHERecHit*, const PhysicsTower*> hits_;
292  std::set<const PhysicsTower*> neighbors_;
293 
294  void calcHits_(void) const;
295  mutable double hitEnergy_;
296  mutable double hitEnergyTrkFid_;
297  mutable int nHits_;
298 
299  void calcHcalSameTowers_(void) const;
300  mutable double hcalEnergySameTowers_;
301  mutable int nHcalHitsSameTowers_;
302 
303  void calcEcalSameTowers_(void) const;
304  mutable double ecalEnergySameTowers_;
305  mutable int nEcalHitsSameTowers_;
306 
307  void calcTracksSameTowers_(void) const;
308  mutable double trackEnergySameTowers_;
309  mutable int nTracksSameTowers_;
310 
311  void calcHcalNeighborTowers_(void) const;
314 
315  void calcEcalNeighborTowers_(void) const;
318 
319  void calcTracksNeighborTowers_(void) const;
322 
323 };
324 
326 //
327 // HBHEHitMapOrganizer
328 //
329 // Organizers the HBHEHitMaps into RBXs, HPDs, dihits, and monohits
331 
333 {
334  public:
336  const ObjectValidatorAbs& objvalidator,
337  const PhysicsTowerOrganizer& pto);
338 
339  virtual ~HBHEHitMapOrganizer() {}
340 
341  void getRBXs(std::vector<HBHEHitMap>& v, double energy) const;
342  void getHPDs(std::vector<HBHEHitMap>& v, double energy) const;
343  void getDiHits(std::vector<HBHEHitMap>& v, double energy) const;
344  void getMonoHits(std::vector<HBHEHitMap>& v, double energy) const;
345 
346  private:
347 
348  std::map<int, HBHEHitMap> rbxs_, hpds_;
349  std::vector<HBHEHitMap> dihits_, monohits_;
350 
351  // helper functions
352  // finds all of the hits which are neighbors and in the same HPD as the reference hit
353  void getHPDNeighbors(const HBHERecHit* hit, std::vector<const HBHERecHit*>& neighbors, const PhysicsTowerOrganizer& pto);
354 
355 };
356 
357 
358 #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
void setHcalChannelQuality(const HcalChannelQuality *q)
void getHPDs(std::vector< HBHEHitMap > &v, double energy) const
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
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
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)
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
HBHEHitMapOrganizer(const edm::Handle< HBHERecHitCollection > &hbhehitcoll_h, const ObjectValidatorAbs &objvalidator, const PhysicsTowerOrganizer &pto)
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
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