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
00026 if ( sev == EcalSeverityLevelAlgo::kProblematic
00027 || sev == EcalSeverityLevelAlgo::kRecovered || sev == EcalSeverityLevelAlgo::kBad )
00028 {
00029
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
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
00055 for ( it = hitsAndFracs.begin(); it != hitsAndFracs.end(); ++it )
00056 {
00057 DetId id = (*it).first;
00058 if ( id != (*itn) )
00059 continue;
00060
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
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
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
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
00098 uint32_t sev = sevlv.severityLevel( *it, recHits);
00099 if (sev == EcalSeverityLevelAlgo::kGood)
00100 continue;
00101
00102
00103 int deta=EBDetId::distanceEta(EBDetId(seed),EBDetId(*it));
00104 int dphi=EBDetId::distancePhi(EBDetId(seed),EBDetId(*it));
00105 if (sqrt(deta*deta + dphi*dphi) < minDist)
00106 closestProb = *it;
00107 }
00108
00109 return closestProb;
00110 }
00111
00112 std::pair<int,int> EcalClusterSeverityLevelAlgo::etaphiDistanceClosestProblematic( const reco::CaloCluster & cluster,
00113 const EcalRecHitCollection & recHits,
00114 const CaloTopology* topology, const EcalSeverityLevelAlgo& sevlv )
00115 {
00116 DetId seed=EcalClusterTools::getMaximum(cluster,&recHits).first;
00117 if ( (seed.det() != DetId::Ecal) ||
00118 (EcalSubdetector(seed.subdetId()) != EcalBarrel) )
00119 {
00120 edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster seed is not in the BARREL";
00121
00122 return std::pair<int,int>(-1,-1);
00123 }
00124
00125 DetId closestProb = closestProblematic(cluster , recHits, topology, sevlv);
00126
00127 if (! closestProb.null())
00128 return std::pair<int,int>(EBDetId::distanceEta(EBDetId(seed),EBDetId(closestProb)),EBDetId::distancePhi(EBDetId(seed),EBDetId(closestProb)));
00129 else
00130 return std::pair<int,int>(-1,-1);
00131 }
00132
00133
00134
00135
00136
00137
00138