Skip to content

Commit

Permalink
Merge pull request #13117 from aminnj/update-hcal-strip-80X
Browse files Browse the repository at this point in the history
Update hcal strip 80 x
  • Loading branch information
davidlange6 committed Feb 1, 2016
2 parents 873bfaf + a468287 commit cecda2a
Show file tree
Hide file tree
Showing 5 changed files with 54 additions and 6 deletions.
6 changes: 6 additions & 0 deletions DataFormats/METReco/interface/HcalHaloData.h
Original file line number Diff line number Diff line change
Expand Up @@ -14,13 +14,19 @@
struct HaloTowerStrip {
std::vector<std::pair<uint8_t, CaloTowerDetId> > cellTowerIds;
float hadEt;
float energyRatio;
float emEt;
HaloTowerStrip() {
hadEt = 0.0;
energyRatio = -1.0;
emEt = 0.0;
cellTowerIds.clear();
}
HaloTowerStrip(const HaloTowerStrip& strip) {
cellTowerIds = strip.cellTowerIds;
hadEt = strip.hadEt;
energyRatio = strip.energyRatio;
emEt = strip.emEt;
}
};

Expand Down
3 changes: 2 additions & 1 deletion DataFormats/METReco/src/classes_def.xml
Original file line number Diff line number Diff line change
Expand Up @@ -172,7 +172,8 @@
</class>
<class name="edm::Wrapper<reco::HcalHaloData>"/>

<class name="HaloTowerStrip" ClassVersion="2">
<class name="HaloTowerStrip" ClassVersion="3">
<version ClassVersion="3" checksum="3149281836"/>
<version ClassVersion="2" checksum="665188928"/>
</class>
<class name="std::vector<HaloTowerStrip>"/>
Expand Down
39 changes: 36 additions & 3 deletions RecoMET/METAlgorithms/src/HcalHaloAlgo.cc
Original file line number Diff line number Diff line change
@@ -1,4 +1,5 @@
#include "RecoMET/METAlgorithms/interface/HcalHaloAlgo.h"
#include <map>

/*
[class]: HcalHaloAlgo
Expand Down Expand Up @@ -128,24 +129,38 @@ HcalHaloData HcalHaloAlgo::Calculate(const CaloGeometry& TheCaloGeometry, edm::H
}
}


// Don't use HF.
int maxAbsIEta = 29;


std::map<int, float> iPhiHadEtMap;
std::vector<const CaloTower*> sortedCaloTowers;
for(CaloTowerCollection::const_iterator tower = TheCaloTowers->begin(); tower != TheCaloTowers->end(); tower++) {
if(abs(tower->ieta()) <= maxAbsIEta && tower->numProblematicHcalCells() > 0)
sortedCaloTowers.push_back(&(*tower));
if(abs(tower->ieta()) > maxAbsIEta) continue;

int iPhi = tower->iphi();
if(!iPhiHadEtMap.count(iPhi)) iPhiHadEtMap[iPhi] = 0.0;
iPhiHadEtMap[iPhi] += tower->hadEt();

if(tower->numProblematicHcalCells() > 0) sortedCaloTowers.push_back(&(*tower));

}


// Sort towers such that lowest iphi and ieta are first, highest last, and towers
// with same iphi value are consecutive. Then we can do everything else in one loop.
std::sort(sortedCaloTowers.begin(), sortedCaloTowers.end(), CompareTowers);

HaloTowerStrip strip;


int prevIEta = -99, prevIPhi = -99;
float prevHadEt = 0.;
float prevEmEt = 0.;
std::pair<uint8_t, CaloTowerDetId> prevPair, towerPair;
bool wasContiguous = true;

// Loop through and store a vector of pairs (problematicCells, DetId) for each contiguous strip we find
for(unsigned int i = 0; i < sortedCaloTowers.size(); i++) {
const CaloTower* tower = sortedCaloTowers[i];
Expand All @@ -162,22 +177,40 @@ HcalHaloData HcalHaloAlgo::Calculate(const CaloGeometry& TheCaloGeometry, edm::H
strip.cellTowerIds.push_back(prevPair);
strip.cellTowerIds.push_back(towerPair);
strip.hadEt += prevHadEt + tower->hadEt();
strip.emEt += prevEmEt + tower->emEt();
}

if(wasContiguous && isContiguous) {
strip.cellTowerIds.push_back(towerPair);
strip.hadEt += tower->hadEt();
strip.emEt += tower->emEt();
}

if((wasContiguous && !isContiguous) || i == sortedCaloTowers.size()-1) { //ended the strip, so flush it
if(strip.cellTowerIds.size() > 2) {

if(strip.cellTowerIds.size() > 3) {

int iPhi = strip.cellTowerIds.at(0).second.iphi();
int iPhiLower = (iPhi == 1) ? 72 : iPhi - 1;
int iPhiUpper = (iPhi == 72) ? 1 : iPhi + 1;

float energyRatio = 0.0;
if(iPhiHadEtMap.count(iPhiLower)) energyRatio += iPhiHadEtMap[iPhiLower];
if(iPhiHadEtMap.count(iPhiUpper)) energyRatio += iPhiHadEtMap[iPhiUpper];
iPhiHadEtMap[iPhi] = max(iPhiHadEtMap[iPhi], 0.001F);

energyRatio /= iPhiHadEtMap[iPhi];
strip.energyRatio = energyRatio;

TheHcalHaloData.getProblematicStrips().push_back( strip );

}
strip = HaloTowerStrip();
}

wasContiguous = isContiguous;
prevPair = towerPair;
prevEmEt = tower->emEt();
prevIPhi = tower->iphi();
prevIEta = tower->ieta();
prevHadEt = tower->hadEt();
Expand Down
8 changes: 7 additions & 1 deletion RecoMET/METFilters/plugins/HcalStripHaloFilter.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,12 +18,16 @@ class HcalStripHaloFilter : public edm::global::EDFilter<> {

const bool taggingMode_;
const int maxWeightedStripLength_;
const double maxEnergyRatio_;
const double minHadEt_;
edm::EDGetTokenT<reco::BeamHaloSummary> beamHaloSummaryToken_;
};

HcalStripHaloFilter::HcalStripHaloFilter(const edm::ParameterSet & iConfig)
: taggingMode_ (iConfig.getParameter<bool> ("taggingMode"))
, maxWeightedStripLength_ (iConfig.getParameter<int>("maxWeightedStripLength"))
, maxEnergyRatio_ (iConfig.getParameter<double>("maxEnergyRatio"))
, minHadEt_ (iConfig.getParameter<double>("minHadEt"))
, beamHaloSummaryToken_(consumes<reco::BeamHaloSummary>(edm::InputTag("BeamHaloSummary")))
{
produces<bool>();
Expand All @@ -42,7 +46,9 @@ bool HcalStripHaloFilter::filter(edm::StreamID iID, edm::Event & iEvent, const e
for (unsigned int iTower = 0; iTower < problematicStrip.cellTowerIds.size(); iTower++) {
numContiguousCells += (int)problematicStrip.cellTowerIds[iTower].first;
}
if(numContiguousCells > maxWeightedStripLength_) {
if( numContiguousCells > maxWeightedStripLength_
&& problematicStrip.energyRatio < maxEnergyRatio_
&& problematicStrip.hadEt > minHadEt_ ) {
pass = false;
break;
}
Expand Down
4 changes: 3 additions & 1 deletion RecoMET/METFilters/python/HcalStripHaloFilter_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,5 +3,7 @@
HcalStripHaloFilter = cms.EDFilter(
"HcalStripHaloFilter",
taggingMode = cms.bool(False),
maxWeightedStripLength = cms.int32(9) # values higher than this are rejected
maxWeightedStripLength = cms.int32(7),
maxEnergyRatio = cms.double(0.15),
minHadEt = cms.double(100.0)
)

0 comments on commit cecda2a

Please sign in to comment.