CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/RecoEcal/EgammaClusterAlgos/src/Multi5x5ClusterAlgo.cc

Go to the documentation of this file.
00001 
00002 #include "RecoEcal/EgammaClusterAlgos/interface/Multi5x5ClusterAlgo.h"
00003 
00004 // Geometry
00005 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00006 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00007 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00008 #include "Geometry/CaloTopology/interface/EcalEndcapTopology.h"
00009 #include "Geometry/CaloTopology/interface/EcalBarrelTopology.h"
00010 #include "RecoEcal/EgammaCoreTools/interface/ClusterEtLess.h"
00011 
00012 #include "DataFormats/CaloRecHit/interface/CaloID.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 //temporary sorter, I'm sure this must exist already in CMSSW
00016 template <class T1,class T2,typename Comp=std::less<T1> > struct PairSortByFirst { 
00017   Comp comp;
00018   bool operator()(const std::pair<T1,T2>& lhs,const std::pair<T1,T2>&rhs){return comp(lhs.first,rhs.first);}
00019   bool operator()(const T1& lhs,const std::pair<T1,T2>&rhs){return comp(lhs,rhs.first);}
00020   bool operator()(const std::pair<T1,T2>& lhs,const T1 &rhs){return comp(lhs.first,rhs);}
00021   bool operator()(const T1& lhs,const T1 &rhs){return comp(lhs,rhs);}
00022 };
00023 
00024 struct EBDetIdSorter{
00025   bool operator()(const EBDetId& lhs,const EBDetId& rhs){
00026     if(lhs.ieta()<rhs.ieta()) return true;
00027     else if(lhs.ieta()>rhs.ieta()) return false;
00028     else return lhs.iphi()<rhs.iphi(); 
00029   }
00030 };
00031 
00032 struct EEDetIdSorter{
00033   bool operator()(const EEDetId& lhs,const EEDetId& rhs){
00034     if(lhs.zside()<rhs.zside()) return true;
00035     else if(lhs.zside()>rhs.zside()) return false;
00036     else { //z is equal, onto ix
00037       if(lhs.ix()<rhs.ix()) return true;
00038       else if(lhs.ix()>rhs.ix()) return false;
00039       else return lhs.iy()<rhs.iy();
00040     }
00041   }
00042 };
00043 
00044 // Return a vector of clusters from a collection of EcalRecHits:
00045 //
00046 std::vector<reco::BasicCluster> Multi5x5ClusterAlgo::makeClusters(
00047         const EcalRecHitCollection* hits,
00048         const CaloSubdetectorGeometry *geometry_p,
00049         const CaloSubdetectorTopology *topology_p,
00050         const CaloSubdetectorGeometry *geometryES_p,
00051         reco::CaloID::Detectors detector,
00052         bool regional,
00053         const std::vector<EcalEtaPhiRegion>& regions
00054         )
00055 {
00056     seeds.clear();
00057     used_s.clear();
00058     canSeed_s.clear();
00059     clusters_v.clear();
00060 
00061     recHits_ = hits;
00062 
00063     double threshold = 0;
00064     std::string ecalPart_string;
00065     detector_ = reco::CaloID::DET_NONE;
00066     if (detector == reco::CaloID::DET_ECAL_ENDCAP) 
00067     {
00068         detector_ = reco::CaloID::DET_ECAL_ENDCAP;
00069         threshold = ecalEndcapSeedThreshold;
00070         ecalPart_string = "EndCap";
00071     }
00072     if (detector == reco::CaloID::DET_ECAL_BARREL) 
00073     {
00074         detector_ = reco::CaloID::DET_ECAL_BARREL;
00075         threshold = ecalBarrelSeedThreshold;
00076         ecalPart_string = "Barrel";
00077     }
00078 
00079     LogTrace("EcalClusters") << "-------------------------------------------------------------";
00080     LogTrace("EcalClusters") << "Island algorithm invoked for ECAL" << ecalPart_string ;
00081     LogTrace("EcalClusters") << "Looking for seeds, energy threshold used = " << threshold << " GeV";
00082 
00083 
00084     int nregions=0;
00085     if(regional) nregions=regions.size();
00086 
00087     if(!regional || nregions) {
00088 
00089         EcalRecHitCollection::const_iterator it;
00090         for(it = hits->begin(); it != hits->end(); it++)
00091         {
00092             double energy = it->energy();
00093             if (energy < threshold) continue; // need to check to see if this line is useful!
00094 
00095             const CaloCellGeometry *thisCell = geometry_p->getGeometry(it->id());
00096             GlobalPoint position = thisCell->getPosition();
00097 
00098             // Require that RecHit is within clustering region in case
00099             // of regional reconstruction
00100             bool withinRegion = false;
00101             if (regional) {
00102                 std::vector<EcalEtaPhiRegion>::const_iterator region;
00103                 for (region=regions.begin(); region!=regions.end(); region++) {
00104                     if (region->inRegion(position)) {
00105                         withinRegion =  true;
00106                         break;
00107                     }
00108                 }
00109             }
00110 
00111             if (!regional || withinRegion) {
00112                 float ET = it->energy() * sin(position.theta());
00113                 if (ET > threshold) seeds.push_back(*it);
00114             }
00115         }
00116 
00117     }
00118 
00119     sort(seeds.begin(), seeds.end(), EcalRecHitLess());
00120 
00121     LogTrace("EcalClusters") << "Total number of seeds found in event = " << seeds.size();
00122 
00123 
00124     mainSearch(hits, geometry_p, topology_p, geometryES_p);
00125     sort(clusters_v.rbegin(), clusters_v.rend(), ClusterEtLess());
00126 
00127     LogTrace("EcalClusters") << "---------- end of main search. clusters have been sorted ----";
00128 
00129 
00130     return clusters_v;
00131 
00132 }
00133 
00134 // Search for clusters
00135 //
00136 void Multi5x5ClusterAlgo::mainSearch(const EcalRecHitCollection* hits,
00137         const CaloSubdetectorGeometry *geometry_p,
00138         const CaloSubdetectorTopology *topology_p,
00139         const CaloSubdetectorGeometry *geometryES_p
00140         )
00141 {
00142 
00143     LogTrace("EcalClusters") << "Building clusters............";
00144 
00145     // Loop over seeds:
00146     std::vector<EcalRecHit>::iterator it;
00147     for (it = seeds.begin(); it != seeds.end(); it++)
00148     {
00149 
00150         // check if this crystal is able to seed
00151         // (event though it is already used)
00152         bool usedButCanSeed = false;
00153         if (canSeed_s.find(it->id()) != canSeed_s.end()) usedButCanSeed = true;
00154 
00155         // avoid seeding for anomalous channels (recoFlag based)
00156         uint32_t rhFlag = (*it).recoFlag();
00157         std::vector<int>::const_iterator vit = std::find( v_chstatus_.begin(), v_chstatus_.end(), rhFlag );
00158         if ( vit != v_chstatus_.end() ) continue; // the recHit has to be excluded from seeding
00159 
00160         // make sure the current seed does not belong to a cluster already.
00161         if ((used_s.find(it->id()) != used_s.end()) && (usedButCanSeed == false))
00162         {
00163             if (it == seeds.begin())
00164             {
00165                 LogTrace("EcalClusters") << "##############################################################" ;
00166                 LogTrace("EcalClusters") << "DEBUG ALERT: Highest energy seed already belongs to a cluster!";
00167                 LogTrace("EcalClusters") << "##############################################################";
00168 
00169             }
00170 
00171             // seed crystal is used or is used and cannot seed a cluster
00172             // so continue to the next seed crystal...
00173             continue;
00174         }
00175 
00176         // clear the vector of hits in current cluster
00177         current_v.clear();
00178 
00179         // Create a navigator at the seed and get seed
00180         // energy
00181         CaloNavigator<DetId> navigator(it->id(), topology_p);
00182         DetId seedId = navigator.pos();
00183         EcalRecHitCollection::const_iterator seedIt = hits->find(seedId);
00184         navigator.setHome(seedId);
00185 
00186         // Is the seed a local maximum?
00187         bool localMaxima = checkMaxima(navigator, hits);
00188 
00189         if (localMaxima)
00190         {
00191             // build the 5x5 taking care over which crystals
00192             // can seed new clusters and which can't
00193             prepareCluster(navigator, hits, geometry_p);
00194         }
00195 
00196         // If some crystals in the current vector then 
00197         // make them into a cluster 
00198         if (current_v.size() > 0) 
00199         {
00200             makeCluster(hits, geometry_p, geometryES_p, seedIt, usedButCanSeed);
00201         }
00202 
00203     }  // End loop on seed crystals
00204 
00205     
00206     if(reassignSeedCrysToClusterItSeeds_){
00207       std::sort(whichClusCrysBelongsTo_.begin(),whichClusCrysBelongsTo_.end(),PairSortByFirst<DetId,int>());
00208       
00209   
00210       for(size_t clusNr=0;clusNr<protoClusters_.size();clusNr++){
00211         if(!protoClusters_[clusNr].containsSeed()){
00212           const EcalRecHit& seedHit =protoClusters_[clusNr].seed();
00213           typedef std::vector<std::pair<DetId,int> >::iterator It;
00214           std::pair<It,It> result = std::equal_range(whichClusCrysBelongsTo_.begin(),whichClusCrysBelongsTo_.end(),seedHit.id(),PairSortByFirst<DetId,int>());
00215           
00216           if(result.first!=result.second) protoClusters_[result.first->second].removeHit(seedHit);
00217           protoClusters_[clusNr].addSeed();
00218           
00219         }
00220       }
00221     }
00222     
00223     
00224       
00225     for(size_t clusNr=0;clusNr<protoClusters_.size();clusNr++){
00226       const ProtoBasicCluster& protoCluster= protoClusters_[clusNr];
00227       Point position;
00228       position = posCalculator_.Calculate_Location(protoCluster.hits(), hits,geometry_p, geometryES_p);
00229       clusters_v.push_back(reco::BasicCluster(protoCluster.energy(), position, reco::CaloID(detector_), protoCluster.hits(),
00230                                               reco::CaloCluster::multi5x5, protoCluster.seed().id()));
00231     }
00232       
00233     protoClusters_.clear();
00234     whichClusCrysBelongsTo_.clear();
00235 }
00236 
00237 void Multi5x5ClusterAlgo::makeCluster(const EcalRecHitCollection* hits,
00238         const CaloSubdetectorGeometry *geometry,
00239         const CaloSubdetectorGeometry *geometryES,
00240         const EcalRecHitCollection::const_iterator &seedIt,
00241         bool seedOutside)
00242 {
00243 
00244     double energy = 0;
00245     //double chi2   = 0;
00246     reco::CaloID caloID;
00247     Point position;
00248     position = posCalculator_.Calculate_Location(current_v, hits,geometry, geometryES);
00249 
00250     std::vector<std::pair<DetId, float> >::iterator it;
00251     for (it = current_v.begin(); it != current_v.end(); it++)
00252     {
00253         EcalRecHitCollection::const_iterator itt = hits->find( (*it).first );
00254         EcalRecHit hit_p = *itt;
00255         energy += hit_p.energy();
00256         //chi2 += 0;
00257         if ( (*it).first.subdetId() == EcalBarrel ) {
00258             caloID = reco::CaloID::DET_ECAL_BARREL;
00259         } else {
00260             caloID = reco::CaloID::DET_ECAL_ENDCAP;
00261         }
00262 
00263     }
00264     //chi2 /= energy;
00265 
00266     LogTrace("EcalClusters") << "******** NEW CLUSTER ********";
00267     LogTrace("EcalClusters") << "No. of crystals = " << current_v.size();
00268     LogTrace("EcalClusters") << "     Energy     = " << energy ;
00269     LogTrace("EcalClusters") << "     Phi        = " << position.phi();
00270     LogTrace("EcalClusters") << "     Eta " << position.eta();
00271     LogTrace("EcalClusters") << "*****************************";  
00272 
00273 
00274     // to be a valid cluster the cluster energy
00275     // must be at least the seed energy
00276     double seedEnergy = seedIt->energy();
00277     if ((seedOutside && energy>=0) || (!seedOutside && energy >= seedEnergy)) 
00278     {
00279       if(reassignSeedCrysToClusterItSeeds_){ //if we're not doing this, we dont need this info so lets not bother filling it
00280         for(size_t hitNr=0;hitNr<current_v.size();hitNr++) whichClusCrysBelongsTo_.push_back(std::pair<DetId,int>(current_v[hitNr].first,protoClusters_.size()));
00281       }
00282         protoClusters_.push_back(ProtoBasicCluster(energy,*seedIt,current_v));
00283       
00284       // clusters_v.push_back(reco::BasicCluster(energy, position, reco::CaloID(detector_), current_v, reco::CaloCluster::multi5x5, seedIt->id()));
00285         
00286     // if no valid cluster was built,
00287     // then free up these crystals to be used in the next...
00288     } else {
00289         std::vector<std::pair<DetId, float> >::iterator iter;
00290         for (iter = current_v.begin(); iter != current_v.end(); iter++)
00291         {
00292             used_s.erase(iter->first);
00293         } //for(iter)
00294     } //else
00295 
00296 }
00297 
00298 bool Multi5x5ClusterAlgo::checkMaxima(CaloNavigator<DetId> &navigator,
00299         const EcalRecHitCollection *hits)
00300 {
00301 
00302     bool maxima = true;
00303     EcalRecHitCollection::const_iterator thisHit;
00304     EcalRecHitCollection::const_iterator seedHit = hits->find(navigator.pos());
00305     double seedEnergy = seedHit->energy();
00306 
00307     std::vector<DetId> swissCrossVec;
00308     swissCrossVec.clear();
00309 
00310     swissCrossVec.push_back(navigator.west());
00311     navigator.home();
00312     swissCrossVec.push_back(navigator.east());
00313     navigator.home();
00314     swissCrossVec.push_back(navigator.north());
00315     navigator.home();
00316     swissCrossVec.push_back(navigator.south());
00317     navigator.home();
00318 
00319     std::vector<DetId>::const_iterator detItr;
00320     for (unsigned int i = 0; i < swissCrossVec.size(); ++i)
00321     {
00322 
00323         // look for this hit
00324         thisHit = recHits_->find(swissCrossVec[i]);
00325 
00326         // continue if this hit was not found
00327         if  ((swissCrossVec[i] == DetId(0)) || thisHit == recHits_->end()) continue; 
00328 
00329         // the recHit has to be skipped in the local maximum search if it was found
00330         // in the map of channels to be excluded 
00331         uint32_t rhFlag = thisHit->recoFlag();
00332         std::vector<int>::const_iterator vit = std::find(v_chstatus_.begin(), v_chstatus_.end(), rhFlag);
00333         if (vit != v_chstatus_.end()) continue;
00334 
00335         // if this crystal has more energy than the seed then we do 
00336         // not have a local maxima
00337         if (thisHit->energy() > seedEnergy)
00338         {
00339             maxima = false;
00340             break;
00341         }
00342     }
00343 
00344     return maxima;
00345 
00346 }
00347 
00348 void Multi5x5ClusterAlgo::prepareCluster(CaloNavigator<DetId> &navigator, 
00349         const EcalRecHitCollection *hits, 
00350         const CaloSubdetectorGeometry *geometry)
00351 {
00352 
00353     DetId thisDet;
00354     std::set<DetId>::iterator setItr;
00355 
00356     // now add the 5x5 taking care to mark the edges
00357     // as able to seed and where overlapping in the central
00358     // region with crystals that were previously able to seed
00359     // change their status so they are not able to seed
00360     //std::cout << std::endl;
00361     for (int dx = -2; dx < 3; ++dx)
00362     {
00363         for (int dy = -2; dy < 3; ++ dy)
00364         {
00365 
00366             // navigate in free steps forming
00367             // a full 5x5
00368             thisDet = navigator.offsetBy(dx, dy);
00369             navigator.home();
00370 
00371             // add the current crystal
00372             //std::cout << "adding " << dx << ", " << dy << std::endl;
00373             addCrystal(thisDet);
00374 
00375             // now consider if we are in an edge (outer 16)
00376             // or central (inner 9) region
00377             if ((abs(dx) > 1) || (abs(dy) > 1))
00378             {    
00379                 // this is an "edge" so should be allowed to seed
00380                 // provided it is not already used
00381                 //std::cout << "   setting can seed" << std::endl;
00382                 canSeed_s.insert(thisDet);
00383             }  // end if "edge"
00384             else 
00385             {
00386                 // or else we are in the central 3x3
00387                 // and must remove any of these crystals from the canSeed set
00388                 setItr = canSeed_s.find(thisDet);
00389                 if (setItr != canSeed_s.end())
00390                 {
00391                     //std::cout << "   unsetting can seed" << std::endl;
00392                     canSeed_s.erase(setItr);
00393                 }
00394             }  // end if "centre"
00395 
00396 
00397         } // end loop on dy
00398 
00399     } // end loop on dx
00400 
00401     //std::cout << "*** " << std::endl;
00402     //std::cout << " current_v contains " << current_v.size() << std::endl;
00403     //std::cout << "*** " << std::endl;
00404 }
00405 
00406 
00407 void Multi5x5ClusterAlgo::addCrystal(const DetId &det)
00408 {   
00409 
00410     EcalRecHitCollection::const_iterator thisIt =  recHits_->find(det);
00411     if ((thisIt != recHits_->end()) && (thisIt->id() != DetId(0)))
00412     { 
00413         if ((used_s.find(thisIt->id()) == used_s.end())) 
00414         {
00415             //std::cout << "   ... this is a good crystal and will be added" << std::endl;
00416             current_v.push_back( std::pair<DetId, float>(det, 1.) ); // by default hit energy fractions are set at 1.
00417             used_s.insert(det);
00418         }
00419     } 
00420 
00421 }
00422 
00423 bool Multi5x5ClusterAlgo::ProtoBasicCluster::removeHit(const EcalRecHit& hitToRM)
00424 {
00425   std::vector<std::pair<DetId,float> >::iterator hitPos;
00426   for(hitPos=hits_.begin();hitPos<hits_.end();hitPos++){
00427     if(hitToRM.id()==hitPos->first) break;
00428   }
00429   if(hitPos!=hits_.end()){
00430     hits_.erase(hitPos);
00431     energy_-=hitToRM.energy();
00432     return true;
00433   }return false;
00434 }
00435 
00436 //now the tricky thing is to make sure we insert it in the right place in the hit vector
00437 //order is from -2,-2, -2,-1 to 2,2 with seed being 0,0 (eta,phi)
00438 bool Multi5x5ClusterAlgo::ProtoBasicCluster::addSeed()
00439 {
00440   typedef std::vector<std::pair<DetId,float> >::iterator It;
00441   std::pair<It,It> hitPos;
00442  
00443   if(seed_.id().subdetId()==EcalBarrel){
00444     hitPos = std::equal_range(hits_.begin(),hits_.end(),seed_.id(),PairSortByFirst<DetId,float,EBDetIdSorter>());
00445   }else{
00446     hitPos = std::equal_range(hits_.begin(),hits_.end(),seed_.id(),PairSortByFirst<DetId,float,EEDetIdSorter>());
00447   }
00448 
00449   if(hitPos.first==hitPos.second){//it doesnt already exist in the vec, add it
00450     hits_.insert(hitPos.first,std::pair<DetId,float>(seed_.id(),1.));
00451     energy_+=seed_.energy(); 
00452     containsSeed_=true;
00453 
00454     return true;
00455   }else return false;
00456  
00457 }
00458 
00459 bool Multi5x5ClusterAlgo::ProtoBasicCluster::isSeedCrysInHits_()const
00460 {
00461   for(size_t hitNr=0;hitNr<hits_.size();hitNr++){
00462     if(seed_.id()==hits_[hitNr].first) return true;
00463   }
00464   return false;
00465 }