CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoEcal/EgammaCoreTools/src/EcalClusterSeverityLevelAlgo.cc

Go to the documentation of this file.
00001 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
00004 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterSeverityLevelAlgo.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00007 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00008 #include "Geometry/CaloTopology/interface/CaloTopology.h"
00009 
00010 float EcalClusterSeverityLevelAlgo::goodFraction( const reco::CaloCluster & cluster, 
00011                                                   const EcalRecHitCollection & recHits, const EcalSeverityLevelAlgo& sevlv) 
00012 {
00013         float fraction = 0.;
00014         std::vector< std::pair<DetId, float> > hitsAndFracs = cluster.hitsAndFractions();
00015         std::vector< std::pair<DetId, float> >::const_iterator it;
00016         for ( it = hitsAndFracs.begin(); it != hitsAndFracs.end(); ++it ) {
00017                 DetId id = (*it).first;
00018                 EcalRecHitCollection::const_iterator jrh = recHits.find( id );
00019                 if ( jrh == recHits.end() ) {
00020                         edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster DetId " << id.rawId() << " is not in the recHit collection!!";
00021                         return -1;
00022                 }
00023                       
00024                 uint32_t sev = sevlv.severityLevel( id, recHits);
00025                 //                if ( sev == EcalSeverityLevelAlgo::kBad ) ++recoveryFailed;
00026                 if ( sev == EcalSeverityLevel::kProblematic 
00027                      || sev == EcalSeverityLevel::kRecovered || sev == EcalSeverityLevel::kBad ) 
00028                   {
00029 //                  std::cout << "[goodFraction] Found a problematic channel " << EBDetId(id) << " " << flag << " energy: " <<  (*jrh).energy() << std::endl;
00030                     fraction += (*jrh).energy() * (*it).second / cluster.energy();
00031                   }
00032         }
00033         return 1. - fraction;
00034 }
00035 
00036 float EcalClusterSeverityLevelAlgo::fractionAroundClosestProblematic( const reco::CaloCluster & cluster, 
00037                                                                       const EcalRecHitCollection & recHits, const CaloTopology* topology, const EcalSeverityLevelAlgo& sevlv )
00038 {
00039   DetId closestProb = closestProblematic(cluster , recHits, topology, sevlv);
00040   //  std::cout << "%%%%%%%%%%% Closest prob is " << EBDetId(closestProb) << std::endl;
00041   if (closestProb.null())
00042     return 0.;
00043   
00044   std::vector<DetId> neighbours = topology->getWindow(closestProb,3,3);
00045   std::vector<DetId>::const_iterator itn;
00046 
00047   std::vector< std::pair<DetId, float> > hitsAndFracs = cluster.hitsAndFractions();
00048   std::vector< std::pair<DetId, float> >::const_iterator it;
00049 
00050   float fraction = 0.;  
00051 
00052   for ( itn = neighbours.begin(); itn != neighbours.end(); ++itn )
00053     {
00054       //      std::cout << "Checking detId " << EBDetId((*itn)) << std::endl;
00055       for ( it = hitsAndFracs.begin(); it != hitsAndFracs.end(); ++it ) 
00056         {
00057           DetId id = (*it).first;
00058           if ( id != (*itn) )
00059             continue;
00060           //      std::cout << "Is in cluster detId " << EBDetId(id) << std::endl;
00061           EcalRecHitCollection::const_iterator jrh = recHits.find( id );
00062           if ( jrh == recHits.end() ) 
00063             {
00064               edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster DetId " << id.rawId() << " is not in the recHit collection!!";
00065               return -1;
00066             }
00067 
00068           fraction += (*jrh).energy() * (*it).second  / cluster.energy();
00069         }
00070     }
00071   //  std::cout << "%%%%%%%%%%% Fraction is " << fraction << std::endl;
00072   return fraction;
00073 }
00074 
00075 DetId EcalClusterSeverityLevelAlgo::closestProblematic(const reco::CaloCluster & cluster, 
00076                                                      const EcalRecHitCollection & recHits, 
00077                                                        const CaloTopology* topology ,  const EcalSeverityLevelAlgo& sevlv)
00078 {
00079   DetId seed=EcalClusterTools::getMaximum(cluster,&recHits).first;
00080   if ( (seed.det() != DetId::Ecal) || 
00081        (EcalSubdetector(seed.subdetId()) != EcalBarrel) )
00082     {
00083       //method not supported if not in Barrel
00084       edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster seed is not in the BARREL";
00085       return DetId(0);
00086     }
00087 
00088   int minDist=9999; DetId closestProb(0);   
00089   //Get a window of DetId around the seed crystal
00090   std::vector<DetId> neighbours = topology->getWindow(seed,51,11);
00091 
00092   for ( std::vector<DetId>::const_iterator it = neighbours.begin(); it != neighbours.end(); ++it ) 
00093     {
00094       EcalRecHitCollection::const_iterator jrh = recHits.find(*it);
00095       if ( jrh == recHits.end() ) 
00096         continue;
00097       //Now checking rh flag   
00098       uint32_t sev = sevlv.severityLevel( *it, recHits);
00099       if (sev == EcalSeverityLevel::kGood)
00100         continue;
00101       //      std::cout << "[closestProblematic] Found a problematic channel " << EBDetId(*it) << " " << flag << std::endl;
00102       //Find the closest DetId in eta,phi space (distance defined by deta^2 + dphi^2)
00103       int deta=EBDetId::distanceEta(EBDetId(seed),EBDetId(*it));
00104       int dphi=EBDetId::distancePhi(EBDetId(seed),EBDetId(*it));
00105       double r = sqrt(deta*deta + dphi*dphi);
00106       if (r < minDist){
00107         closestProb = *it;
00108         minDist = r;
00109       }
00110     }
00111       
00112   return closestProb;
00113 }
00114 
00115 std::pair<int,int> EcalClusterSeverityLevelAlgo::etaphiDistanceClosestProblematic( const reco::CaloCluster & cluster, 
00116                                                      const EcalRecHitCollection & recHits, 
00117                                                                                    const CaloTopology* topology,  const EcalSeverityLevelAlgo& sevlv )
00118 {
00119   DetId seed=EcalClusterTools::getMaximum(cluster,&recHits).first;
00120   if ( (seed.det() != DetId::Ecal) || 
00121        (EcalSubdetector(seed.subdetId()) != EcalBarrel) )
00122     {
00123       edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster seed is not in the BARREL";
00124       //method not supported if not in Barrel
00125       return std::pair<int,int>(-1,-1);
00126     }
00127 
00128   DetId closestProb = closestProblematic(cluster , recHits, topology, sevlv);
00129 
00130   if (! closestProb.null())
00131     return std::pair<int,int>(EBDetId::distanceEta(EBDetId(seed),EBDetId(closestProb)),EBDetId::distancePhi(EBDetId(seed),EBDetId(closestProb)));
00132   else
00133     return std::pair<int,int>(-1,-1);
00134 } 
00135 
00136 
00137 
00138 
00139 
00140 
00141