CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoLocalCalo/HcalRecAlgos/interface/HBHEIsolatedNoiseAlgos.h

Go to the documentation of this file.
00001 #ifndef __HBHE_ISOLATED_NOISE_ALGOS_H__
00002 #define __HBHE_ISOLATED_NOISE_ALGOS_H__
00003 
00004 /*
00005 Description: Isolation algorithms used to identify anomalous noise in the HB/HE.
00006              These algorithms will be used to reflag HB/HE rechits as noise.
00007 
00008              There are 6 objects defined here:
00009              1) ObjectValidatorAbs
00010              2) ObjectValidator
00011              3) PhysicsTower
00012              4) PhysicsTowerOrganizer
00013              5) HBHEHitMap
00014              6) HBHEHitMapOrganizer
00015              See comments below for details.
00016 
00017 Original Author: John Paul Chou (Brown University)
00018                  Thursday, September 2, 2010
00019 */
00020 
00021 
00022 // system include files
00023 #include <memory>
00024 #include <string>
00025 #include <vector>
00026 
00027 // user include files
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00030 
00031 #include "CondFormats/EcalObjects/interface/EcalChannelStatus.h"
00032 #include "Geometry/CaloTopology/interface/CaloTowerConstituentsMap.h"
00033 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00034 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00035 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00036 #include "DataFormats/JetReco/interface/TrackExtrapolation.h"
00037 
00038 #include <vector>
00039 #include <set>
00040 #include <map>
00041 
00042 // forward declarations
00043 class HBHERecHit;
00044 class EcalRecHit;
00045 class HcalChannelQuality;
00046 class HcalSeverityLevelComputer;
00047 class EcalSeverityLevelAlgo;
00048 
00050 //
00051 // ObjectValidator
00052 //
00053 // determines if an ECAL hit, HCAL hit, or track is valid or not.
00054 // Note that various objects need to be passed to the validator before use.
00055 // ObjectValidatorAbs provides a base class in case you want to change the
00056 // hit/track validation algorithms.
00058 
00059 class ObjectValidatorAbs
00060 {
00061  public:
00062   ObjectValidatorAbs() {}
00063   virtual ~ObjectValidatorAbs() {}
00064 
00065   virtual bool validHit(const HBHERecHit&) const = 0;
00066   virtual bool validHit(const EcalRecHit&) const = 0;
00067   virtual bool validTrack(const reco::Track&) const = 0;
00068 };
00069 
00070 class ObjectValidator : public ObjectValidatorAbs
00071 {
00072  public:
00073   explicit ObjectValidator(const edm::ParameterSet&);
00074   ObjectValidator(double HBThreshold, double HESThreshold, double HEDThreshold, double EBThreshold, double EEThreshold,
00075                   uint32_t HcalAcceptSeverityLevel, uint32_t EcalAcceptSeverityLevel, bool UseHcalRecoveredHits, bool UseEcalRecoveredHits,
00076                   double MinValidTrackPt, double MinValidTrackPtBarrel, int MinValidTrackNHits) :
00077     HBThreshold_(HBThreshold),
00078     HESThreshold_(HESThreshold),
00079     HEDThreshold_(HEDThreshold),
00080     EBThreshold_(EBThreshold),
00081     EEThreshold_(EEThreshold),
00082     HcalAcceptSeverityLevel_(HcalAcceptSeverityLevel),
00083     EcalAcceptSeverityLevel_(EcalAcceptSeverityLevel),
00084     UseHcalRecoveredHits_(UseHcalRecoveredHits),
00085     UseEcalRecoveredHits_(UseEcalRecoveredHits),
00086     MinValidTrackPt_(MinValidTrackPt),
00087     MinValidTrackPtBarrel_(MinValidTrackPtBarrel),
00088     MinValidTrackNHits_(MinValidTrackNHits),
00089     theHcalChStatus_(0),
00090     theEcalChStatus_(0),
00091     theHcalSevLvlComputer_(0),
00092     theEcalSevLvlAlgo_(0),
00093     theEBRecHitCollection_(0),
00094     theEERecHitCollection_(0) {}
00095   virtual ~ObjectValidator();
00096   
00097   inline void setHcalChannelQuality(const HcalChannelQuality* q) { theHcalChStatus_=q; }
00098   inline void setEcalChannelStatus(const EcalChannelStatus* q) { theEcalChStatus_=q; }
00099   inline void setHcalSeverityLevelComputer(const HcalSeverityLevelComputer* q) { theHcalSevLvlComputer_=q; }
00100   inline void setEcalSeverityLevelAlgo(const EcalSeverityLevelAlgo* q) { theEcalSevLvlAlgo_=q; }
00101   inline void setEBRecHitCollection(const EcalRecHitCollection* q) { theEBRecHitCollection_=q; }
00102   inline void setEERecHitCollection(const EcalRecHitCollection* q) { theEERecHitCollection_=q; }
00103 
00104 
00105   bool validHit(const HBHERecHit&) const;
00106   bool validHit(const EcalRecHit&) const;
00107   bool validTrack(const reco::Track&) const;
00108 
00109  private:
00110 
00111   double HBThreshold_;  // energy threshold for HB hits
00112   double HESThreshold_; // energy threshold for 5-degree (phi) HE hits
00113   double HEDThreshold_; // energy threshold for 10-degree (phi) HE hits
00114   double EBThreshold_;  // energy threshold for EB hits
00115   double EEThreshold_;  // energy threshold for EE hits
00116 
00117   uint32_t HcalAcceptSeverityLevel_; // severity level to accept HCAL hits
00118   uint32_t EcalAcceptSeverityLevel_; // severity level to accept ECAL hits
00119   bool UseHcalRecoveredHits_; // whether or not to use recovered HCAL hits
00120   bool UseEcalRecoveredHits_; // whether or not to use recovered HCAL hits
00121 
00122   double MinValidTrackPt_; // minimum valid track pT
00123   double MinValidTrackPtBarrel_; // minimum valid track pT in the Barrel
00124   int MinValidTrackNHits_; // minimum number of hits needed for a valid track
00125   
00126   // for checking the status of ECAL and HCAL channels stored in the DB 
00127   const HcalChannelQuality* theHcalChStatus_;
00128   const EcalChannelStatus* theEcalChStatus_;
00129   
00130   // calculator of severety level for ECAL and HCAL
00131   const HcalSeverityLevelComputer* theHcalSevLvlComputer_;
00132   const EcalSeverityLevelAlgo* theEcalSevLvlAlgo_;
00133 
00134   // needed to determine an ECAL channel's validity
00135   const EcalRecHitCollection* theEBRecHitCollection_;
00136   const EcalRecHitCollection* theEERecHitCollection_;
00137 };
00138 
00139 
00141 //
00142 // PhysicsTower
00143 //
00144 // basic structure to hold information relevant at the tower level
00145 // consists of hcal hits, ecal hits, and tracks
00147 
00148 struct PhysicsTower {
00149   CaloTowerDetId id;
00150   std::set<const HBHERecHit*> hcalhits;
00151   std::set<const EcalRecHit*> ecalhits;
00152   std::set<const reco::Track*> tracks;
00153 };
00154 
00156 //
00157 // PhysicsTowerOrganizer
00158 //
00159 // Organizers the PhysicsTowers into CaloTower (ieta,iphi) space.
00160 // Also provides methods to find a towers nearest neighbors.
00161 // For simplicity, merges the ieta=28 and ieta=29 towers together (calls them
00162 // tower "28", collectively).
00164 
00165 class PhysicsTowerOrganizer {
00166  public:
00167   struct towercmp {
00168     bool operator() (const PhysicsTower& lhs, const PhysicsTower& rhs) const {
00169       return (lhs.id < rhs.id);
00170     }
00171   };
00172 
00173   PhysicsTowerOrganizer(const edm::Event& iEvent,
00174                         const edm::EventSetup& evSetup,
00175                         const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
00176                         const edm::Handle<EcalRecHitCollection>& ebhitcoll_h,
00177                         const edm::Handle<EcalRecHitCollection>& eehitcoll_h,
00178                         const edm::Handle<std::vector<reco::TrackExtrapolation> >& trackextrapcoll_h,
00179                         const ObjectValidatorAbs& objectvalidator,
00180                         const CaloTowerConstituentsMap& ctcm);
00181 
00182   virtual ~PhysicsTowerOrganizer() {}
00183 
00184   // find a PhysicsTower by some coordinate
00185   inline const PhysicsTower* findTower(const CaloTowerDetId& id) const;
00186   inline const PhysicsTower* findTower(int ieta, int iphi) const;
00187   
00188   // get the neighbors +/- 1 in eta-space or +/- 1 in phi-space
00189   // (accounts for change in phi-segmentation starting with eta=21)
00190   void findNeighbors(const CaloTowerDetId& id, std::set<const PhysicsTower*>& neighbors) const;
00191   void findNeighbors(const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors) const;
00192   void findNeighbors(int ieta, int iphi, std::set<const PhysicsTower*>& neighbors) const;  
00193 
00194  private:
00195   // the non-const, private version of findTower()
00196   PhysicsTower* findTower(const CaloTowerDetId& id);
00197   PhysicsTower* findTower(int ieta, int iphi);
00198 
00199   void insert_(CaloTowerDetId& id, const HBHERecHit* hit);
00200   void insert_(CaloTowerDetId& id, const EcalRecHit* hit);
00201   void insert_(CaloTowerDetId& id, const reco::Track* hit);
00202 
00203   std::set<PhysicsTower, towercmp> towers_;
00204   
00205 };
00206 
00208 //
00209 // HBHEHitMap
00210 //
00211 // Collection of HBHERecHits and their associated PhysicsTowers.
00212 // Typically, these maps are organized into RBXs, HPDs, dihits, or monohits.
00213 // From this one may calculate the hcal, ecal, and track energy associated
00214 // with these hits.
00215 //
00216 // N.B. Many, of the operations can be computationally expensive and should be
00217 // used with care.
00219 
00220 
00221 class HBHEHitMap {
00222 
00223  public:
00224 
00225   typedef std::map<const HBHERecHit*, const PhysicsTower*>::const_iterator hitmap_const_iterator;
00226   typedef std::set<const PhysicsTower*>::const_iterator neighbor_const_iterator;
00227 
00228   struct twrinfo {
00229     double hcalInMap;
00230     double hcalOutOfMap;
00231     double ecal;
00232     double track;
00233   };
00234 
00235   HBHEHitMap();
00236   virtual ~HBHEHitMap() {}
00237 
00238   // energy of the hits in this collection
00239   double hitEnergy(void) const;
00240 
00241   // energy of the hits in a region fiducial to tracks
00242   double hitEnergyTrackFiducial(void) const;
00243 
00244   // number of hits in this collection
00245   int nHits(void) const;
00246 
00247   // same as above, except for the HCAL hits, ECAL hits, and tracks in the same towers as the hits
00248   // note that HCAL hits may be present in the same tower, but not be a part of the collection
00249   double hcalEnergySameTowers(void) const;
00250   double ecalEnergySameTowers(void) const;
00251   double trackEnergySameTowers(void) const;
00252   int nHcalHitsSameTowers(void) const;
00253   int nEcalHitsSameTowers(void) const;
00254   int nTracksSameTowers(void) const;
00255   void hcalHitsSameTowers(std::set<const HBHERecHit*>& v) const;
00256   void ecalHitsSameTowers(std::set<const EcalRecHit*>& v) const;
00257   void tracksSameTowers(std::set<const reco::Track*>& v) const;
00258 
00259   // same as above except for the hits and tracks in the neighboring towers to the hits in the collection
00260   double hcalEnergyNeighborTowers(void) const;
00261   double ecalEnergyNeighborTowers(void) const;
00262   double trackEnergyNeighborTowers(void) const;
00263   int nHcalHitsNeighborTowers(void) const;
00264   int nEcalHitsNeighborTowers(void) const;
00265   int nTracksNeighborTowers(void) const;
00266   void hcalHitsNeighborTowers(std::set<const HBHERecHit*>& v) const;
00267   void ecalHitsNeighborTowers(std::set<const EcalRecHit*>& v) const;
00268   void tracksNeighborTowers(std::set<const reco::Track*>& v) const;
00269 
00270   // gives the total energy due to hcal, ecal, and tracks tower-by-tower
00271   // separates the hcal energy into that which is a hit in the collection, and that which is not
00272   void byTowers(std::vector<twrinfo>& v) const;
00273 
00274   // find a hit in the hitmap
00275   hitmap_const_iterator findHit(const HBHERecHit* hit) const { return hits_.find(hit); }
00276 
00277   // find a neighbor
00278   neighbor_const_iterator findNeighbor(const PhysicsTower* twr) const { return neighbors_.find(twr); }
00279 
00280   // add a hit to the hitmap
00281   void insert(const HBHERecHit* hit, const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors);
00282 
00283   // access to the private maps and sets
00284   inline hitmap_const_iterator beginHits(void) const { return hits_.begin(); }
00285   inline hitmap_const_iterator endHits(void) const { return hits_.end(); }
00286 
00287   inline neighbor_const_iterator beginNeighbors(void) const { return neighbors_.begin(); }
00288   inline neighbor_const_iterator endNeighbors(void) const { return neighbors_.end(); }
00289 
00290  private:
00291   std::map<const HBHERecHit*, const PhysicsTower*> hits_;
00292   std::set<const PhysicsTower*> neighbors_;
00293 
00294   void calcHits_(void) const;
00295   mutable double hitEnergy_;
00296   mutable double hitEnergyTrkFid_;
00297   mutable int nHits_;
00298 
00299   void calcHcalSameTowers_(void) const;
00300   mutable double hcalEnergySameTowers_;
00301   mutable int nHcalHitsSameTowers_;
00302 
00303   void calcEcalSameTowers_(void) const;
00304   mutable double ecalEnergySameTowers_;
00305   mutable int nEcalHitsSameTowers_;
00306 
00307   void calcTracksSameTowers_(void) const;
00308   mutable double trackEnergySameTowers_;
00309   mutable int nTracksSameTowers_;
00310 
00311   void calcHcalNeighborTowers_(void) const;
00312   mutable double hcalEnergyNeighborTowers_;
00313   mutable int nHcalHitsNeighborTowers_;
00314 
00315   void calcEcalNeighborTowers_(void) const;
00316   mutable double ecalEnergyNeighborTowers_;
00317   mutable int nEcalHitsNeighborTowers_;
00318 
00319   void calcTracksNeighborTowers_(void) const;
00320   mutable double trackEnergyNeighborTowers_;
00321   mutable int nTracksNeighborTowers_;
00322 
00323 };
00324 
00326 //
00327 // HBHEHitMapOrganizer
00328 //
00329 // Organizers the HBHEHitMaps into RBXs, HPDs, dihits, and monohits
00331 
00332 class HBHEHitMapOrganizer
00333 {
00334  public:
00335   HBHEHitMapOrganizer(const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
00336                       const ObjectValidatorAbs& objvalidator,
00337                       const PhysicsTowerOrganizer& pto);
00338 
00339   virtual ~HBHEHitMapOrganizer() {}
00340 
00341   void getRBXs(std::vector<HBHEHitMap>& v, double energy) const;
00342   void getHPDs(std::vector<HBHEHitMap>& v, double energy) const;
00343   void getDiHits(std::vector<HBHEHitMap>& v, double energy) const;
00344   void getMonoHits(std::vector<HBHEHitMap>& v, double energy) const;
00345 
00346  private:
00347   
00348   std::map<int, HBHEHitMap> rbxs_, hpds_;
00349   std::vector<HBHEHitMap> dihits_, monohits_;
00350 
00351   // helper functions
00352   // finds all of the hits which are neighbors and in the same HPD as the reference hit
00353   void getHPDNeighbors(const HBHERecHit* hit, std::vector<const HBHERecHit*>& neighbors, const PhysicsTowerOrganizer& pto);
00354 
00355 };
00356 
00357 
00358 #endif