#include <EcalClusterSeverityLevelAlgo.h>
Static Public Member Functions | |
static DetId | closestProblematic (const reco::CaloCluster &, const EcalRecHitCollection &, const CaloTopology *topology, const EcalSeverityLevelAlgo &) |
static std::pair< int, int > | etaphiDistanceClosestProblematic (const reco::CaloCluster &, const EcalRecHitCollection &, const CaloTopology *topology, const EcalSeverityLevelAlgo &) |
static float | fractionAroundClosestProblematic (const reco::CaloCluster &, const EcalRecHitCollection &, const CaloTopology *topology, const EcalSeverityLevelAlgo &) |
static float | goodFraction (const reco::CaloCluster &, const EcalRecHitCollection &, const EcalSeverityLevelAlgo &) |
Definition at line 14 of file EcalClusterSeverityLevelAlgo.h.
DetId EcalClusterSeverityLevelAlgo::closestProblematic | ( | const reco::CaloCluster & | cluster, |
const EcalRecHitCollection & | recHits, | ||
const CaloTopology * | topology, | ||
const EcalSeverityLevelAlgo & | sevlv | ||
) | [static] |
Definition at line 75 of file EcalClusterSeverityLevelAlgo.cc.
References DetId::det(), EBDetId::distanceEta(), EBDetId::distancePhi(), DetId::Ecal, EcalBarrel, edm::SortedCollection< T, SORT >::end(), edm::SortedCollection< T, SORT >::find(), EcalClusterTools::getMaximum(), CaloTopology::getWindow(), EcalSeverityLevelAlgo::kGood, EcalSeverityLevelAlgo::severityLevel(), mathSSE::sqrt(), and DetId::subdetId().
Referenced by etaphiDistanceClosestProblematic(), and fractionAroundClosestProblematic().
{ DetId seed=EcalClusterTools::getMaximum(cluster,&recHits).first; if ( (seed.det() != DetId::Ecal) || (EcalSubdetector(seed.subdetId()) != EcalBarrel) ) { //method not supported if not in Barrel edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster seed is not in the BARREL"; return DetId(0); } int minDist=9999; DetId closestProb(0); //Get a window of DetId around the seed crystal std::vector<DetId> neighbours = topology->getWindow(seed,51,11); for ( std::vector<DetId>::const_iterator it = neighbours.begin(); it != neighbours.end(); ++it ) { EcalRecHitCollection::const_iterator jrh = recHits.find(*it); if ( jrh == recHits.end() ) continue; //Now checking rh flag uint32_t sev = sevlv.severityLevel( *it, recHits); if (sev == EcalSeverityLevelAlgo::kGood) continue; // std::cout << "[closestProblematic] Found a problematic channel " << EBDetId(*it) << " " << flag << std::endl; //Find the closest DetId in eta,phi space (distance defined by deta^2 + dphi^2) int deta=EBDetId::distanceEta(EBDetId(seed),EBDetId(*it)); int dphi=EBDetId::distancePhi(EBDetId(seed),EBDetId(*it)); if (sqrt(deta*deta + dphi*dphi) < minDist) closestProb = *it; } return closestProb; }
std::pair< int, int > EcalClusterSeverityLevelAlgo::etaphiDistanceClosestProblematic | ( | const reco::CaloCluster & | cluster, |
const EcalRecHitCollection & | recHits, | ||
const CaloTopology * | topology, | ||
const EcalSeverityLevelAlgo & | sevlv | ||
) | [static] |
Definition at line 112 of file EcalClusterSeverityLevelAlgo.cc.
References closestProblematic(), DetId::det(), EBDetId::distanceEta(), EBDetId::distancePhi(), DetId::Ecal, EcalBarrel, EcalClusterTools::getMaximum(), DetId::null(), and DetId::subdetId().
{ DetId seed=EcalClusterTools::getMaximum(cluster,&recHits).first; if ( (seed.det() != DetId::Ecal) || (EcalSubdetector(seed.subdetId()) != EcalBarrel) ) { edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster seed is not in the BARREL"; //method not supported if not in Barrel return std::pair<int,int>(-1,-1); } DetId closestProb = closestProblematic(cluster , recHits, topology, sevlv); if (! closestProb.null()) return std::pair<int,int>(EBDetId::distanceEta(EBDetId(seed),EBDetId(closestProb)),EBDetId::distancePhi(EBDetId(seed),EBDetId(closestProb))); else return std::pair<int,int>(-1,-1); }
float EcalClusterSeverityLevelAlgo::fractionAroundClosestProblematic | ( | const reco::CaloCluster & | cluster, |
const EcalRecHitCollection & | recHits, | ||
const CaloTopology * | topology, | ||
const EcalSeverityLevelAlgo & | sevlv | ||
) | [static] |
Definition at line 36 of file EcalClusterSeverityLevelAlgo.cc.
References closestProblematic(), edm::SortedCollection< T, SORT >::end(), reco::CaloCluster::energy(), edm::SortedCollection< T, SORT >::find(), CaloTopology::getWindow(), reco::CaloCluster::hitsAndFractions(), and DetId::null().
{ DetId closestProb = closestProblematic(cluster , recHits, topology, sevlv); // std::cout << "%%%%%%%%%%% Closest prob is " << EBDetId(closestProb) << std::endl; if (closestProb.null()) return 0.; std::vector<DetId> neighbours = topology->getWindow(closestProb,3,3); std::vector<DetId>::const_iterator itn; std::vector< std::pair<DetId, float> > hitsAndFracs = cluster.hitsAndFractions(); std::vector< std::pair<DetId, float> >::const_iterator it; float fraction = 0.; for ( itn = neighbours.begin(); itn != neighbours.end(); ++itn ) { // std::cout << "Checking detId " << EBDetId((*itn)) << std::endl; for ( it = hitsAndFracs.begin(); it != hitsAndFracs.end(); ++it ) { DetId id = (*it).first; if ( id != (*itn) ) continue; // std::cout << "Is in cluster detId " << EBDetId(id) << std::endl; EcalRecHitCollection::const_iterator jrh = recHits.find( id ); if ( jrh == recHits.end() ) { edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster DetId " << id.rawId() << " is not in the recHit collection!!"; return -1; } fraction += (*jrh).energy() * (*it).second / cluster.energy(); } } // std::cout << "%%%%%%%%%%% Fraction is " << fraction << std::endl; return fraction; }
float EcalClusterSeverityLevelAlgo::goodFraction | ( | const reco::CaloCluster & | cluster, |
const EcalRecHitCollection & | recHits, | ||
const EcalSeverityLevelAlgo & | sevlv | ||
) | [static] |
Definition at line 10 of file EcalClusterSeverityLevelAlgo.cc.
References edm::SortedCollection< T, SORT >::end(), reco::CaloCluster::energy(), edm::SortedCollection< T, SORT >::find(), reco::CaloCluster::hitsAndFractions(), EcalSeverityLevelAlgo::kBad, EcalSeverityLevelAlgo::kProblematic, EcalSeverityLevelAlgo::kRecovered, and EcalSeverityLevelAlgo::severityLevel().
{ float fraction = 0.; std::vector< std::pair<DetId, float> > hitsAndFracs = cluster.hitsAndFractions(); std::vector< std::pair<DetId, float> >::const_iterator it; for ( it = hitsAndFracs.begin(); it != hitsAndFracs.end(); ++it ) { DetId id = (*it).first; EcalRecHitCollection::const_iterator jrh = recHits.find( id ); if ( jrh == recHits.end() ) { edm::LogError("EcalClusterSeverityLevelAlgo") << "The cluster DetId " << id.rawId() << " is not in the recHit collection!!"; return -1; } uint32_t sev = sevlv.severityLevel( id, recHits); // if ( sev == EcalSeverityLevelAlgo::kBad ) ++recoveryFailed; if ( sev == EcalSeverityLevelAlgo::kProblematic || sev == EcalSeverityLevelAlgo::kRecovered || sev == EcalSeverityLevelAlgo::kBad ) { // std::cout << "[goodFraction] Found a problematic channel " << EBDetId(id) << " " << flag << " energy: " << (*jrh).energy() << std::endl; fraction += (*jrh).energy() * (*it).second / cluster.energy(); } } return 1. - fraction; }