CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalCalo/HcalRecAlgos/src/HBHEIsolatedNoiseAlgos.cc

Go to the documentation of this file.
00001 /*
00002 Description: Isolation algorithms used to identify anomalous noise in the HB/HE.
00003              These algorithms will be used to reflag HB/HE rechits as noise.
00004 
00005              There are 4 objects implemented here:
00006              1) ObjectValidator
00007              2) PhysicsTowerOrganizer
00008              3) HBHEHitMap
00009              4) HBHEHitMapOrganizer
00010              See comments below for details.
00011 
00012 Original Author: John Paul Chou (Brown University)
00013                  Thursday, September 2, 2010
00014 */
00015 
00016 #include "RecoLocalCalo/HcalRecAlgos/interface/HBHEIsolatedNoiseAlgos.h"
00017 
00018 #include "FWCore/Utilities/interface/EDMException.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 
00021 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00022 #include "DataFormats/METReco/interface/CaloMET.h"
00023 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00024 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00025 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00026 #include "DataFormats/TrackReco/interface/Track.h"
00027 
00028 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00029 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00030 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00031 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00032 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
00033 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00034 #include "RecoMET/METAlgorithms/interface/HcalHPDRBXMap.h"
00035 
00037 //
00038 // ObjectValidator
00039 //
00041 
00042 ObjectValidator::ObjectValidator(const edm::ParameterSet& iConfig)
00043 {
00044   HBThreshold_  = iConfig.getParameter<double>("HBThreshold");
00045   HESThreshold_ = iConfig.getParameter<double>("HESThreshold");
00046   HEDThreshold_ = iConfig.getParameter<double>("HEDThreshold");
00047   EBThreshold_  = iConfig.getParameter<double>("EBThreshold");
00048   EEThreshold_  = iConfig.getParameter<double>("EEThreshold");
00049 
00050   HcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("HcalAcceptSeverityLevel");
00051   EcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("EcalAcceptSeverityLevel");
00052   UseHcalRecoveredHits_ = iConfig.getParameter<bool>("UseHcalRecoveredHits");
00053   UseEcalRecoveredHits_ = iConfig.getParameter<bool>("UseEcalRecoveredHits");
00054 
00055   MinValidTrackPt_ = iConfig.getParameter<double>("MinValidTrackPt");
00056   MinValidTrackPtBarrel_ = iConfig.getParameter<double>("MinValidTrackPtBarrel");
00057   MinValidTrackNHits_ = iConfig.getParameter<int>("MinValidTrackNHits");
00058 
00059   theHcalChStatus_=0;
00060   theEcalChStatus_=0;
00061   theHcalSevLvlComputer_=0;
00062   theEcalSevLvlAlgo_=0;
00063   theEBRecHitCollection_=0;
00064   theEERecHitCollection_=0;
00065 
00066   return;
00067 }
00068 
00069 ObjectValidator::~ObjectValidator()
00070 {
00071 }
00072 
00073 bool ObjectValidator::validHit(const HBHERecHit& hit) const
00074 {
00075   assert(theHcalSevLvlComputer_!=0 && theHcalChStatus_!=0);
00076 
00077   // require the hit to pass a certain energy threshold
00078   if(hit.id().subdet()==HcalBarrel && hit.energy()<HBThreshold_) return false;
00079   else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()<=20 && hit.energy()<HESThreshold_) return false;
00080   else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()>20 && hit.energy()<HEDThreshold_) return false;
00081 
00082   // determine if the hit is good, bad, or recovered
00083   const DetId id = hit.detid();
00084   const uint32_t recHitFlag = hit.flags();
00085   const uint32_t dbStatusFlag = theHcalChStatus_->getValues(id)->getValue();
00086   int severityLevel = theHcalSevLvlComputer_->getSeverityLevel(id, recHitFlag, dbStatusFlag);
00087   bool isRecovered  = theHcalSevLvlComputer_->recoveredRecHit(id, recHitFlag);
00088 
00089   if(severityLevel==0) return true;
00090   if(isRecovered) return UseHcalRecoveredHits_;
00091   if(severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) return false;
00092   else return true;
00093 }
00094 
00095 bool ObjectValidator::validHit(const EcalRecHit& hit) const
00096 {
00097   assert(theEcalSevLvlAlgo_!=0 && theEcalChStatus_!=0);
00098 
00099   // require the hit to pass a certain energy threshold
00100   const DetId id = hit.detid();
00101   if(id.subdetId()==EcalBarrel && hit.energy()<EBThreshold_) return false;
00102   else if(id.subdetId()==EcalEndcap && hit.energy()<EEThreshold_) return false;
00103 
00104   // determine if the hit is good, bad, or recovered
00105   int severityLevel = 999;
00106   if     (id.subdetId() == EcalBarrel && theEBRecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(hit);//id, *theEBRecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999);
00107   else if(id.subdetId() == EcalEndcap && theEERecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(hit);//id, *theEERecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999);
00108   else return false;
00109   
00110   if(severityLevel == EcalSeverityLevelAlgo::kGood) return true;
00111   if(severityLevel == EcalSeverityLevelAlgo::kRecovered) return UseEcalRecoveredHits_;
00112   if(severityLevel > static_cast<int>(EcalAcceptSeverityLevel_)) return false;
00113   else return true;
00114 }
00115 
00116 bool ObjectValidator::validTrack(const reco::Track& trk) const
00117 {
00118   if(trk.pt()<MinValidTrackPt_) return false;
00119   if(trk.pt()<MinValidTrackPtBarrel_  && std::fabs(trk.momentum().eta())<1.479) return false;
00120   if(trk.numberOfValidHits()<MinValidTrackNHits_) return false;
00121   return true;
00122 }
00123 
00125 //
00126 // PhysicsTowerOrganizer
00127 //
00129 
00130 PhysicsTowerOrganizer::PhysicsTowerOrganizer(const edm::Event& iEvent,
00131                                              const edm::EventSetup& evSetup,
00132                                              const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
00133                                              const edm::Handle<EcalRecHitCollection>& ebhitcoll_h,
00134                                              const edm::Handle<EcalRecHitCollection>& eehitcoll_h,
00135                                              const edm::Handle<std::vector<reco::TrackExtrapolation> >& trackextrapcoll_h,
00136                                              const ObjectValidatorAbs& objectvalidator,
00137                                              const CaloTowerConstituentsMap& ctcm)
00138 {
00139   // get some geometries
00140   edm::ESHandle<CaloGeometry> pG;
00141   evSetup.get<CaloGeometryRecord>().get(pG);
00142   const CaloGeometry* geo = pG.product();
00143   const CaloSubdetectorGeometry* gEB = geo->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);
00144   const CaloSubdetectorGeometry* gEE = geo->getSubdetectorGeometry(DetId::Ecal,EcalEndcap);
00145 
00146   // do the HCAL hits
00147   for(HBHERecHitCollection::const_iterator it=hbhehitcoll_h->begin(); it!=hbhehitcoll_h->end(); ++it) {
00148     const HBHERecHit* hit=&(*it);
00149       
00150     // check that the hit is valid
00151     if(!objectvalidator.validHit(*hit)) continue;
00152       
00153     // add the hit to the organizer
00154     CaloTowerDetId tid = ctcm.towerOf(hit->id());
00155     insert_(tid, hit);
00156   }
00157 
00158   // do the EB hits
00159   for(EcalRecHitCollection::const_iterator it=ebhitcoll_h->begin(); it!=ebhitcoll_h->end(); ++it) {
00160     const EcalRecHit* hit=&(*it);
00161       
00162     if(!objectvalidator.validHit(*hit)) continue;
00163     CaloTowerDetId tid = ctcm.towerOf(hit->id());
00164     insert_(tid, hit);
00165   }
00166 
00167   // do the EE hits
00168   for(EcalRecHitCollection::const_iterator it=eehitcoll_h->begin(); it!=eehitcoll_h->end(); ++it) {
00169     const EcalRecHit* hit=&(*it);
00170     
00171     if(!objectvalidator.validHit(*hit)) continue;
00172     CaloTowerDetId tid = ctcm.towerOf(hit->id());
00173     insert_(tid, hit);
00174   }
00175   
00176   // do the tracks
00177   for(std::vector<reco::TrackExtrapolation>::const_iterator it=trackextrapcoll_h->begin(); it!=trackextrapcoll_h->end(); ++it) {
00178     const reco::TrackExtrapolation* extrap=&(*it);
00179     const reco::Track* track = &(*(extrap->track()));
00180 
00181     // validate track
00182     if(!objectvalidator.validTrack(*track)) continue;
00183     
00184     // get the point
00185     if ( extrap->positions().size()==0 ) continue; 
00186     const GlobalPoint point(extrap->positions().front().x(),
00187                             extrap->positions().front().y(),
00188                             extrap->positions().front().z());
00189 
00190     
00191     if(std::fabs(point.eta())<1.479) {
00192       EBDetId cell = gEB->getClosestCell(point);
00193       CaloTowerDetId tid = ctcm.towerOf(cell);
00194       insert_(tid, track);
00195     } else {
00196       EEDetId cell = gEE->getClosestCell(point);
00197       CaloTowerDetId tid = ctcm.towerOf(cell);
00198       insert_(tid, track);
00199     }
00200   }
00201 
00202   return;
00203 }
00204 
00205 PhysicsTower* PhysicsTowerOrganizer::findTower(const CaloTowerDetId& id)
00206 {
00207   // create dummy PhysicsTower
00208   PhysicsTower dummy;
00209 
00210   // correct for the merging of the |ieta|=28-29 towers
00211   if(id.ietaAbs()==29) dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
00212   else dummy.id=id;
00213 
00214   // search on the dummy
00215   std::set<PhysicsTower, towercmp>::iterator it=towers_.find(dummy);
00216   
00217   if(it==towers_.end()) return 0;
00218 
00219   // for whatever reason, I can't get a non-const out of the find method
00220   PhysicsTower &twr = const_cast<PhysicsTower&>(*it);
00221   return &twr;
00222 }
00223 
00224 const PhysicsTower* PhysicsTowerOrganizer::findTower(const CaloTowerDetId& id) const
00225 {
00226   // create dummy PhysicsTower
00227   PhysicsTower dummy;
00228 
00229   // correct for the merging of the |ieta|=28-29 towers
00230   if(id.ietaAbs()==29) dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
00231   else dummy.id=id;
00232 
00233   // search on the dummy
00234   std::set<PhysicsTower, towercmp>::iterator it=towers_.find(dummy);
00235   
00236   if(it==towers_.end()) return 0;
00237   return &(*it);
00238 }
00239 
00240 const PhysicsTower* PhysicsTowerOrganizer::findTower(int ieta, int iphi) const
00241 {
00242   CaloTowerDetId tid(ieta, iphi);
00243   return findTower(tid);
00244 }
00245 
00246 PhysicsTower* PhysicsTowerOrganizer::findTower(int ieta, int iphi)
00247 {
00248   CaloTowerDetId tid(ieta, iphi);
00249   return findTower(tid);
00250 }
00251 
00252 void PhysicsTowerOrganizer::findNeighbors(const CaloTowerDetId& tempid, std::set<const PhysicsTower*>& neighbors) const
00253 {
00254   // correct for the merging of the |ieta|=28-29 towers
00255   CaloTowerDetId id(tempid);
00256   if(tempid.ietaAbs()==29) id = CaloTowerDetId((tempid.ietaAbs()-1)*tempid.zside(), tempid.iphi());
00257 
00258   std::vector<CaloTowerDetId> ids;
00259   // get the neighbor with higher iphi
00260   if(id.ietaAbs()<=20) {
00261     if(id.iphi()==72) ids.push_back(CaloTowerDetId(id.ieta(), 1));
00262     else              ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()+1));
00263   } else {
00264     if(id.iphi()==71) ids.push_back(CaloTowerDetId(id.ieta(), 1));
00265     else              ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()+2));
00266   }
00267 
00268   // get the neighbor with the lower iphi
00269   if(id.ietaAbs()<=20) {
00270     if(id.iphi()==1) ids.push_back(CaloTowerDetId(id.ieta(), 72));
00271     else             ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()-1));
00272   } else {
00273     if(id.iphi()==1) ids.push_back(CaloTowerDetId(id.ieta(), 71));
00274     else             ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()-2));
00275   }
00276 
00277   // get the neighbor with the higher ietaAbs
00278   if(id.ietaAbs()==20 && (id.iphi()%2)==0)
00279     ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-1));
00280   else
00281     ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()));
00282 
00283   // get the neighbor(s) with the lower ietaAbs
00284   if(id.ietaAbs()==21) {
00285     ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()));
00286     ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+1));
00287   } else if(id.ietaAbs()==1) {
00288     ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()));
00289   } else {
00290     ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()));
00291   }
00292   
00293   // get the neighbor with higher ieta and higher iphi
00294   if(id.ietaAbs()<=19 || (id.ietaAbs()==20 && (id.iphi()%2)==0)) {
00295     if(id.iphi()==72) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 1));
00296     else              ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()+1));
00297   } else if(id.ietaAbs()>=21) {
00298     if(id.iphi()==71) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 1));
00299     else              ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()+2));
00300   }
00301 
00302   // get the neighbor with higher ieta and lower iphi
00303   if(id.ietaAbs()<=19) {
00304     if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 72));
00305     else             ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-1));
00306   } else if(id.ietaAbs()>=21 || (id.ietaAbs()==20 && (id.iphi()%2)==1)) {
00307     if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 71));
00308     else             ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-2));
00309   }
00310 
00311   // get the neighbor with lower ieta and higher iphi
00312   if(id.ietaAbs()==1) {
00313     if(id.iphi()==72) ids.push_back(CaloTowerDetId(-id.ieta(), 1));
00314     else              ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()+1));
00315   } else if(id.ietaAbs()<=20) {
00316     if(id.iphi()==72) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 1));
00317     else              ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+1));
00318   } else if(id.ietaAbs()>=21) {
00319     if(id.iphi()==71) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 1));
00320     else              ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+2));
00321   }
00322 
00323   // get the neighbor with lower ieta and lower iphi
00324   if(id.ietaAbs()==1) {
00325     if(id.iphi()==1) ids.push_back(CaloTowerDetId(-id.ieta(), 72));
00326     else             ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()-1));
00327   } else if(id.ietaAbs()<=20) {
00328     if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 72));
00329     else             ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-1));
00330   } else if(id.ietaAbs()>=22) {
00331     if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 71));
00332     else             ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-2));
00333   } else if(id.ietaAbs()==21) {
00334     if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 72));
00335     else             ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-1));
00336   }
00337   
00338   // clear neighbors
00339   neighbors.clear();
00340 
00341   // find the neighbors and add them to the eponymous set
00342   for(std::vector<CaloTowerDetId>::const_iterator it=ids.begin(); it!=ids.end(); ++it) {
00343     const PhysicsTower* twr=findTower(*it);
00344     if(twr) neighbors.insert(twr);
00345   }
00346 
00347   return;
00348 }
00349 
00350 void PhysicsTowerOrganizer::findNeighbors(const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors) const
00351 {
00352   findNeighbors(twr->id, neighbors);
00353   return;
00354 }
00355 
00356 void PhysicsTowerOrganizer::findNeighbors(int ieta, int iphi, std::set<const PhysicsTower*>& neighbors) const
00357 {
00358   findNeighbors(CaloTowerDetId(ieta, iphi), neighbors);
00359   return;
00360 }
00361 
00362 void PhysicsTowerOrganizer::insert_(CaloTowerDetId& id, const HBHERecHit* hit)
00363 {
00364   PhysicsTower* twr=findTower(id);
00365   if(twr==0) {
00366     PhysicsTower dummy;
00367     if(id.ietaAbs()==29)
00368       dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
00369     else
00370       dummy.id = id;
00371     dummy.hcalhits.insert(hit);
00372     towers_.insert(dummy);
00373   } else {
00374     twr->hcalhits.insert(hit);
00375   }
00376   return;
00377 }
00378 
00379 void PhysicsTowerOrganizer::insert_(CaloTowerDetId& id, const EcalRecHit* hit)
00380 {
00381   PhysicsTower* twr=findTower(id);
00382   if(twr==0) {
00383     PhysicsTower dummy;
00384     if(id.ietaAbs()==29)
00385       dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
00386     else
00387       dummy.id = id;
00388     dummy.ecalhits.insert(hit);
00389     towers_.insert(dummy);
00390   } else {
00391     twr->ecalhits.insert(hit);
00392   }
00393   return;
00394 }
00395 
00396 void PhysicsTowerOrganizer::insert_(CaloTowerDetId& id, const reco::Track* track)
00397 {
00398   PhysicsTower* twr=findTower(id);
00399   if(twr==0) {
00400     PhysicsTower dummy;
00401     if(id.ietaAbs()==29)
00402       dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
00403     else
00404       dummy.id = id;
00405     dummy.tracks.insert(track);
00406     towers_.insert(dummy);
00407   } else {
00408     twr->tracks.insert(track);
00409   }
00410   return;
00411 }
00412 
00413 
00415 //
00416 // HBHEHitMap
00417 //
00419 
00420 HBHEHitMap::HBHEHitMap()
00421 {
00422   hitEnergy_=hitEnergyTrkFid_=-999.;
00423   nHits_=-999;
00424   hcalEnergySameTowers_=ecalEnergySameTowers_=trackEnergySameTowers_=-999.;
00425   nHcalHitsSameTowers_=nEcalHitsSameTowers_=nTracksSameTowers_=-999;
00426   hcalEnergyNeighborTowers_=ecalEnergyNeighborTowers_=trackEnergyNeighborTowers_=-999.;
00427   nHcalHitsNeighborTowers_=nEcalHitsNeighborTowers_=nTracksNeighborTowers_=-999;
00428 
00429 }
00430 
00431 double HBHEHitMap::hitEnergy(void) const
00432 {
00433   if(hitEnergy_<-900) calcHits_();
00434   return hitEnergy_;
00435 }
00436 
00437 int HBHEHitMap::nHits(void) const
00438 {
00439   if(nHits_<-900) calcHits_();
00440   return nHits_;
00441 }
00442 
00443 double HBHEHitMap::hitEnergyTrackFiducial(void) const
00444 {
00445   if(hitEnergyTrkFid_<-900) calcHits_();
00446   return hitEnergyTrkFid_;
00447 }
00448 
00449 
00450 double HBHEHitMap::hcalEnergySameTowers(void) const
00451 {
00452   if(hcalEnergySameTowers_<-900) calcHcalSameTowers_();
00453   return hcalEnergySameTowers_;
00454 }
00455 
00456 int HBHEHitMap::nHcalHitsSameTowers(void) const
00457 {
00458   if(nHcalHitsSameTowers_<-900) calcHcalSameTowers_();
00459   return nHcalHitsSameTowers_;
00460 }
00461 
00462 double HBHEHitMap::ecalEnergySameTowers(void) const
00463 {
00464   if(ecalEnergySameTowers_<-900) calcEcalSameTowers_();
00465   return ecalEnergySameTowers_;
00466 }
00467 
00468 int HBHEHitMap::nEcalHitsSameTowers(void) const
00469 {
00470   if(nEcalHitsSameTowers_<-900) calcEcalSameTowers_();
00471   return nEcalHitsSameTowers_;
00472 }
00473 
00474 double HBHEHitMap::trackEnergySameTowers(void) const
00475 {
00476   if(trackEnergySameTowers_<-900) calcTracksSameTowers_();
00477   return trackEnergySameTowers_;
00478 }
00479 
00480 int HBHEHitMap::nTracksSameTowers(void) const
00481 {
00482   if(nTracksSameTowers_<-900) calcTracksSameTowers_();
00483   return nTracksSameTowers_;
00484 }
00485 
00486 void HBHEHitMap::hcalHitsSameTowers(std::set<const HBHERecHit*>& v) const
00487 {
00488   v.clear();
00489   for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
00490     for(std::set<const HBHERecHit*>::const_iterator it2=it1->second->hcalhits.begin(); it2!=it1->second->hcalhits.end(); ++it2) {
00491       const HBHERecHit* hit=(*it2);
00492       // if the hit in the tower is already in the hitmap, don't include it
00493       if(findHit(hit)==endHits()) v.insert(hit);
00494     }
00495   }
00496   return;
00497 }
00498 
00499 void HBHEHitMap::ecalHitsSameTowers(std::set<const EcalRecHit*>& v) const
00500 {
00501   v.clear();
00502   for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
00503     v.insert(it1->second->ecalhits.begin(), it1->second->ecalhits.end());
00504   }
00505   return;
00506 }
00507 
00508 void HBHEHitMap::tracksSameTowers(std::set<const reco::Track*>& v) const
00509 {
00510   v.clear();
00511   for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
00512     v.insert(it1->second->tracks.begin(), it1->second->tracks.end());
00513   }
00514   return;
00515 }
00516 
00517 double HBHEHitMap::hcalEnergyNeighborTowers(void) const
00518 {
00519   if(hcalEnergyNeighborTowers_<-900) calcHcalNeighborTowers_();
00520   return hcalEnergyNeighborTowers_;
00521 }
00522 
00523 int HBHEHitMap::nHcalHitsNeighborTowers(void) const
00524 {
00525   if(nHcalHitsNeighborTowers_<-900) calcHcalNeighborTowers_();
00526   return nHcalHitsNeighborTowers_;
00527 }
00528 
00529 double HBHEHitMap::ecalEnergyNeighborTowers(void) const
00530 {
00531   if(ecalEnergyNeighborTowers_<-900) calcEcalNeighborTowers_();
00532   return ecalEnergyNeighborTowers_;
00533 }
00534 
00535 int HBHEHitMap::nEcalHitsNeighborTowers(void) const
00536 {
00537   if(nEcalHitsNeighborTowers_<-900) calcEcalNeighborTowers_();
00538   return nEcalHitsNeighborTowers_;
00539 }
00540 
00541 double HBHEHitMap::trackEnergyNeighborTowers(void) const
00542 {
00543   if(trackEnergyNeighborTowers_<-900) calcTracksNeighborTowers_();
00544   return trackEnergyNeighborTowers_;
00545 }
00546 
00547 int HBHEHitMap::nTracksNeighborTowers(void) const
00548 {
00549   if(nTracksNeighborTowers_<-900) calcTracksNeighborTowers_();
00550   return nTracksNeighborTowers_;
00551 }
00552 
00553 void HBHEHitMap::hcalHitsNeighborTowers(std::set<const HBHERecHit*>& v) const
00554 {
00555   v.clear();
00556   for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
00557     const PhysicsTower* twr=(*it1);
00558     v.insert(twr->hcalhits.begin(), twr->hcalhits.end());
00559   }
00560   return;
00561 }
00562 
00563 void HBHEHitMap::ecalHitsNeighborTowers(std::set<const EcalRecHit*>& v) const
00564 {
00565   v.clear();
00566   for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
00567     const PhysicsTower* twr=(*it1);
00568     v.insert(twr->ecalhits.begin(), twr->ecalhits.end());
00569   }
00570 
00571   return;
00572 }
00573 
00574 void HBHEHitMap::tracksNeighborTowers(std::set<const reco::Track*>& v) const
00575 {
00576   v.clear();
00577   for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
00578     const PhysicsTower* twr=(*it1);
00579     v.insert(twr->tracks.begin(), twr->tracks.end());
00580   }
00581   return;
00582 }
00583 
00584 void HBHEHitMap::byTowers(std::vector<twrinfo>& v) const
00585 {
00586   assert(0);
00587 }
00588 
00589 void HBHEHitMap::insert(const HBHERecHit* hit, const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors)
00590 {
00591   hits_[hit]=twr;
00592   neighbors_.insert(neighbors.begin(), neighbors.end());
00593 
00594   // make sure none of the neighbors are also are part of the hitmap
00595   for(hitmap_const_iterator it=beginHits(); it!=endHits(); ++it) {
00596     const PhysicsTower* t=it->second;
00597     neighbor_const_iterator find=findNeighbor(t);
00598 
00599     // if a hit is also a neighbor, remove the neighbor
00600     if(find!=endNeighbors()) neighbors_.erase(find);
00601   }
00602   return;
00603 }
00604 
00605 
00606 void HBHEHitMap::calcHits_(void) const
00607 {
00608   hitEnergy_=0;
00609   nHits_=0;
00610   hitEnergyTrkFid_=0;
00611   for(hitmap_const_iterator it=hits_.begin(); it!=hits_.end(); ++it) {
00612     const HBHERecHit* hit=it->first;
00613     if(hit->id().ietaAbs()<=26) hitEnergyTrkFid_+=hit->energy();
00614     hitEnergy_+=hit->energy();
00615     ++nHits_;
00616   }
00617   return;
00618 }
00619 
00620 void HBHEHitMap::calcHcalSameTowers_(void) const
00621 {
00622   hcalEnergySameTowers_=0;
00623   nHcalHitsSameTowers_=0;
00624   std::set<const HBHERecHit*> v;
00625   hcalHitsSameTowers(v);
00626   for(std::set<const HBHERecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00627     const HBHERecHit* hit=(*it);
00628     hcalEnergySameTowers_+=hit->energy();
00629     ++nHcalHitsSameTowers_;
00630   }
00631   return;
00632 }
00633 
00634 void HBHEHitMap::calcEcalSameTowers_(void) const
00635 {
00636   ecalEnergySameTowers_=0;
00637   nEcalHitsSameTowers_=0;
00638   std::set<const EcalRecHit*> v;
00639   ecalHitsSameTowers(v);
00640   for(std::set<const EcalRecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00641     const EcalRecHit* hit=(*it);
00642     ecalEnergySameTowers_+=hit->energy();
00643     ++nEcalHitsSameTowers_;
00644   }
00645   return;
00646 }
00647 
00648 void HBHEHitMap::calcTracksSameTowers_(void) const
00649 {
00650   trackEnergySameTowers_=0;
00651   nTracksSameTowers_=0;
00652   std::set<const reco::Track*> v;
00653   tracksSameTowers(v);
00654   for(std::set<const reco::Track*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00655     const reco::Track* trk=(*it);
00656     trackEnergySameTowers_+=trk->p();
00657     ++nTracksSameTowers_;
00658   }
00659   return;
00660 }
00661 
00662 void HBHEHitMap::calcHcalNeighborTowers_(void) const
00663 {
00664   hcalEnergyNeighborTowers_=0;
00665   nHcalHitsNeighborTowers_=0;
00666   std::set<const HBHERecHit*> v;
00667   hcalHitsNeighborTowers(v);
00668   for(std::set<const HBHERecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00669     const HBHERecHit* hit=(*it);
00670     hcalEnergyNeighborTowers_+=hit->energy();
00671     ++nHcalHitsNeighborTowers_;
00672   }
00673   return;
00674 }
00675 
00676 void HBHEHitMap::calcEcalNeighborTowers_(void) const
00677 {
00678   ecalEnergyNeighborTowers_=0;
00679   nEcalHitsNeighborTowers_=0;
00680   std::set<const EcalRecHit*> v;
00681   ecalHitsNeighborTowers(v);
00682   for(std::set<const EcalRecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00683     const EcalRecHit* hit=(*it);
00684     ecalEnergyNeighborTowers_+=hit->energy();
00685     ++nEcalHitsNeighborTowers_;
00686   }
00687   return;
00688 }
00689 
00690 void HBHEHitMap::calcTracksNeighborTowers_(void) const
00691 {
00692   trackEnergyNeighborTowers_=0;
00693   nTracksNeighborTowers_=0;
00694   std::set<const reco::Track*> v;
00695   tracksNeighborTowers(v);
00696   for(std::set<const reco::Track*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
00697     const reco::Track* trk=(*it);
00698     trackEnergyNeighborTowers_+=trk->p();
00699     ++nTracksNeighborTowers_;
00700   }
00701   return;
00702 }
00703 
00705 //
00706 // HBHEHitMapOrganizer
00707 //
00709 
00710 
00711 HBHEHitMapOrganizer::HBHEHitMapOrganizer(const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
00712                                          const ObjectValidatorAbs& objvalidator,
00713                                          const PhysicsTowerOrganizer& pto)
00714 {
00715   // loop over the hits
00716   for(HBHERecHitCollection::const_iterator it=hbhehitcoll_h->begin(); it!=hbhehitcoll_h->end(); ++it) {
00717     const HBHERecHit *hit=&(*it);
00718     if(!objvalidator.validHit(*hit)) continue;
00719     
00720     // get the Physics Tower and the neighbors
00721     const PhysicsTower* tower=pto.findTower(hit->id().ieta(), hit->id().iphi());
00722     
00723     std::set<const PhysicsTower*> neighbors;
00724     pto.findNeighbors(hit->id().ieta(), hit->id().iphi(), neighbors);
00725     
00726     // organize the RBXs
00727     int rbxidnum = HcalHPDRBXMap::indexRBX(hit->id());
00728     rbxs_[rbxidnum].insert(hit, tower, neighbors);
00729     
00730     // organize the HPDs
00731     int hpdidnum = HcalHPDRBXMap::indexHPD(hit->id());
00732     hpds_[hpdidnum].insert(hit, tower, neighbors);
00733     
00734     
00735     // organize the dihits
00736     std::vector<const HBHERecHit*> hpdneighbors;
00737     getHPDNeighbors(hit, hpdneighbors, pto);
00738     
00739     if(hpdneighbors.size()==1) {
00740       std::vector<const HBHERecHit*> hpdneighborsneighbors;
00741       getHPDNeighbors(hpdneighbors[0], hpdneighborsneighbors, pto);
00742       
00743       if(hpdneighborsneighbors.size()==1 && hpdneighborsneighbors[0]==hit && hit->energy()>hpdneighbors[0]->energy()) {
00744         // we've found two hits who are neighbors in the same HPD, but who have no other
00745         // neighbors (in the same HPD) in common.  In order not to double-count, we
00746         // require that the first hit has more energy
00747         
00748         const PhysicsTower* tower2=pto.findTower(hpdneighbors[0]->id().ieta(), hpdneighbors[0]->id().iphi());
00749         std::set<const PhysicsTower*> neighbors2;
00750         pto.findNeighbors(hpdneighbors[0]->id().ieta(), hpdneighbors[0]->id().iphi(), neighbors2);
00751         
00752         HBHEHitMap dihit;
00753         dihit.insert(hit, tower, neighbors);
00754         dihit.insert(hpdneighbors[0], tower2, neighbors2);
00755         dihits_.push_back(dihit);
00756       }
00757     } else if(hpdneighbors.size()==0) {
00758       
00759       // organize the monohits
00760       HBHEHitMap monohit;
00761       monohit.insert(hit, tower, neighbors);
00762       monohits_.push_back(monohit);
00763     }
00764     
00765   } // finished looping over HBHERecHits
00766   return;
00767 }
00768 
00769 void HBHEHitMapOrganizer::getRBXs(std::vector<HBHEHitMap>& v, double energy) const
00770 {
00771   for(std::map<int, HBHEHitMap>::const_iterator it=rbxs_.begin(); it!=rbxs_.end(); ++it) {
00772     const HBHEHitMap &map=it->second;
00773     if(map.hitEnergy()>energy) v.push_back(map);
00774   }
00775   return;
00776 }
00777 
00778 void HBHEHitMapOrganizer::getHPDs(std::vector<HBHEHitMap>& v, double energy) const
00779 {
00780   for(std::map<int, HBHEHitMap>::const_iterator it=hpds_.begin(); it!=hpds_.end(); ++it) {
00781     const HBHEHitMap &map=it->second;
00782     if(map.hitEnergy()>energy) v.push_back(map);
00783   }
00784   return;
00785 }
00786 
00787 void HBHEHitMapOrganizer::getDiHits(std::vector<HBHEHitMap>& v, double energy) const
00788 {
00789   for(std::vector<HBHEHitMap>::const_iterator it=dihits_.begin(); it!=dihits_.end(); ++it) {
00790     if(it->hitEnergy()>energy) v.push_back(*it);
00791   }
00792   return;
00793 }
00794 
00795 void HBHEHitMapOrganizer::getMonoHits(std::vector<HBHEHitMap>& v, double energy) const
00796 {
00797   for(std::vector<HBHEHitMap>::const_iterator it=monohits_.begin(); it!=monohits_.end(); ++it) {
00798     if(it->hitEnergy()>energy) v.push_back(*it);
00799   }
00800   return;
00801 }
00802 
00803 
00804 
00805 void HBHEHitMapOrganizer::getHPDNeighbors(const HBHERecHit* hit, std::vector<const HBHERecHit*>& neighbors, const PhysicsTowerOrganizer& pto)
00806 {
00807   std::set<const PhysicsTower*> temp;
00808   pto.findNeighbors(hit->id().ieta(), hit->id().iphi(), temp);
00809 
00810   // make sure to include the same tower that the hit is in
00811   temp.insert(pto.findTower(hit->id().ieta(), hit->id().iphi()));
00812 
00813   // loop over the rechits in the temp neighbors
00814   for(std::set<const PhysicsTower*>::const_iterator it1=temp.begin(); it1!=temp.end(); ++it1) {
00815     for(std::set<const HBHERecHit*>::const_iterator it2=(*it1)->hcalhits.begin(); it2!=(*it1)->hcalhits.end(); ++it2) {
00816       const HBHERecHit* hit2(*it2);
00817       if(hit!=hit2 && HcalHPDRBXMap::indexHPD(hit->id())==HcalHPDRBXMap::indexHPD(hit2->id())) {
00818         neighbors.push_back(hit2);
00819       }
00820     }
00821   }
00822   return;
00823 }
00824 
00825 
00826