Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Spike cleaner fix #17859

Merged
merged 7 commits into from
Aug 4, 2017
Merged
Show file tree
Hide file tree
Changes from 6 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Original file line number Diff line number Diff line change
Expand Up @@ -117,19 +117,53 @@ clean(const edm::Handle<reco::PFRecHitCollection>& input,
}
const spike_cleaning& clean = _thresholds.find(hitlayer)->second;
if( rechit.energy() < clean._singleSpikeThresh ) continue;

//Fix needed for HF. Here, we find the (up to) five companion rechits
//to work in conjunction with the neighbours4() implementation below for the full HF surrounding energy
float compsumE = 0.0;
if ((hitlayer==PFLayer::HF_EM || hitlayer==PFLayer::HF_HAD)) {
int comp = 1;
if (hitlayer==PFLayer::HF_EM) comp = 2;
const HcalDetId& detid = (HcalDetId)rechit.detId();
int heta = detid.ieta();
int hphi = detid.iphi();

//At eta>39, phi segmentation changes
int predphi =2;
if (std::abs(heta)>39) predphi=4;

int curphiL = hphi-predphi;
int curphiH = hphi+predphi;

//HcalDetId valid phi range (1-72)
while (curphiL<0) curphiL+=72;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

hi @knash - what about curphiL=0? given the comment, 0 is not handled correctly by this code.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Will change the comment in the next version. I double checked the code and phi=0 indeed should be included

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Indeed, it is a bit unconfortable having those "0", "72" and "39" hardcoded here...

while (curphiH>72) curphiH-=72;

std::pair<std::vector<int>, std::vector<int>> phietas({heta,heta+1,heta-1,heta,heta},{hphi,hphi,hphi,curphiL,curphiH});

std::vector<uint32_t> rawDetIds;
for(int in=0;in<5;in++) {
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

please receive the length of the vector (5) - but better yet, don't hardwire 5, ask the vectors in the pair how many elements they have.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The next version will include this

HcalDetId tempID (HcalForward, phietas.first[in], phietas.second[in], comp);
rawDetIds.push_back(tempID.rawId());
}

for( const auto& jdx : ordered_hits ) {
const unsigned j = jdx;
const reco::PFRecHit& matchrechit = hits[j];
for( const auto& iID : rawDetIds ) if (iID==matchrechit.detId())compsumE+=matchrechit.energy();
}
}
//End of fix needed for HF

const double rhenergy = rechit.energy();
// single spike cleaning
auto const & neighbours4 = rechit.neighbours4();
double surroundingEnergy = rechit.energy();
double neighbourEnergy = 0.0;
double layerEnergy = 0.0;
double surroundingEnergy = compsumE;
for( auto k : neighbours4 ) {
if( !mask[k] ) continue;
const auto & neighbour = hits[k];
const double sum = neighbour.energy(); //energyUp is just rechit energy?
surroundingEnergy += sum;
neighbourEnergy += sum;
layerEnergy += neighbour.energy();
}
// wannaBeSeed.energyUp()/wannaBeSeed.energy() : 1.;
// Fraction 1 is the balance between the hit and its neighbours
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@
#define __SpikeAndDoubleSpikeCleaner_H__

#include "RecoParticleFlow/PFClusterProducer/interface/RecHitTopologicalCleanerBase.h"
#include "DataFormats/HcalRecHit/interface/HFRecHit.h"

#include <unordered_map>

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -104,7 +104,7 @@
particleFlowClusterHF = cms.EDProducer(
"PFClusterProducer",
recHitsSource = cms.InputTag("particleFlowRecHitHF"),
recHitCleaners = cms.VPSet(_spikeAndDoubleSpikeCleaner_HF),
recHitCleaners = cms.VPSet(),
seedFinder = _localMaxSeeds_HF,
initialClusteringStep = _topoClusterizer_HF,
pfClusterBuilder = _pfClusterizer_HF,
Expand Down