6 #include <unordered_map> 27 const std::unordered_map<std::string, int>
_layerMap;
34 std::pair<double, double> dCrack(
double phi,
double eta) {
35 constexpr double oneOverCrystalSize = 1.0 / 0.0175;
59 0.0, 4.44747e-01, -4.44747e-01, 7.92824e-01, -7.92824e-01, 1.14090e+00, -1.14090e+00, 1.47464e+00, -1.47464e+00};
68 if (
phi < cPhi[17] ||
phi >= cPhi[0]) {
89 for (
const double etaGap : cEta) {
92 defi *= oneOverCrystalSize;
93 deta *= oneOverCrystalSize;
94 return std::make_pair(defi, deta);
112 const std::vector<edm::ParameterSet>&
thresholds = conf.getParameterSetVector(
"cleaningByDetector");
116 info._minS4S1_a =
pset.getParameter<
double>(
"minS4S1_a");
117 info._minS4S1_b =
pset.getParameter<
double>(
"minS4S1_b");
118 info._doubleSpikeS6S2 =
pset.getParameter<
double>(
"doubleSpikeS6S2");
119 info._eneThreshMod =
pset.getParameter<
double>(
"energyThresholdModifier");
120 info._fracThreshMod =
pset.getParameter<
double>(
"fractionThresholdModifier");
121 info._doubleSpikeThresh =
pset.getParameter<
double>(
"doubleSpikeThresh");
122 info._singleSpikeThresh =
pset.getParameter<
double>(
"singleSpikeThresh");
123 auto entry = _layerMap.find(det);
124 if (
entry == _layerMap.end()) {
125 throw cms::Exception(
"InvalidDetectorLayer") <<
"Detector layer : " << det <<
" is not in the list of recognized" 126 <<
" detector layers!";
128 _thresholds.emplace(_layerMap.find(det)->second,
info);
135 std::vector<unsigned> ordered_hits(
hits.size());
136 for (
unsigned i = 0;
i <
hits.size(); ++
i)
138 std::sort(ordered_hits.begin(), ordered_hits.end(), [&](
unsigned i,
unsigned j) {
142 for (
const auto&
idx : ordered_hits) {
143 const unsigned i =
idx;
157 float compsumE = 0.0;
163 int heta =
detid.ieta();
164 int hphi =
detid.iphi();
171 int curphiL = hphi - predphi;
172 int curphiH = hphi + predphi;
180 std::pair<std::vector<int>, std::vector<int>> phietas({heta, heta + 1, heta - 1, heta, heta},
181 {hphi, hphi, hphi, curphiL, curphiH});
183 std::vector<uint32_t> rawDetIds;
184 for (
unsigned in = 0;
in < phietas.first.size();
in++) {
186 rawDetIds.push_back(tempID.
rawId());
189 for (
const auto& jdx : ordered_hits) {
190 const unsigned j = jdx;
192 for (
const auto& iID : rawDetIds)
193 if (iID == matchrechit.
detId())
194 compsumE += matchrechit.
energy();
199 const double rhenergy = rechit.
energy();
202 double surroundingEnergy = compsumE;
203 for (
auto k : neighbours4) {
206 auto const& neighbour =
hits[
k];
207 const double sum = neighbour.energy();
208 surroundingEnergy += sum;
213 const double fraction1 = surroundingEnergy / rhenergy;
216 const double f1Cut = (
clean._minS4S1_a * std::log10(rechit.
energy()) +
clean._minS4S1_b);
217 if (fraction1 < f1Cut) {
221 std::pair<double, double> dcr = dCrack(
phi,
eta);
223 if (aeta < 5.0 && ((aeta < 2.85 && dcrmin > 1.0) || (rhenergy >
clean._eneThreshMod *
clean._singleSpikeThresh &&
224 fraction1 < f1Cut /
clean._fracThreshMod))) {
229 if (
mask[
i] && rhenergy >
clean._doubleSpikeThresh) {
231 double surroundingEnergyi = 0.0;
232 double enmax = -999.0;
233 unsigned int mostEnergeticNeighbour = 0;
235 for (
auto k : neighbours4i) {
238 auto const& neighbour =
hits[
k];
239 const double nenergy = neighbour.energy();
240 surroundingEnergyi += nenergy;
241 if (nenergy > enmax) {
243 mostEnergeticNeighbour =
k;
248 double surroundingEnergyj = 0.0;
249 auto const& neighbours4j =
hits[mostEnergeticNeighbour].neighbours4();
250 for (
auto k : neighbours4j) {
252 surroundingEnergyj +=
hits[
k].energy();
255 const double surroundingEnergyFraction =
256 (surroundingEnergyi + surroundingEnergyj) / (rechit.
energy() +
hits[mostEnergeticNeighbour].energy()) - 1.;
257 if (surroundingEnergyFraction <
clean._doubleSpikeS6S2) {
261 std::pair<double, double> dcr = dCrack(
phi,
eta);
263 if (aeta < 5.0 && ((aeta < 2.85 && dcrmin > 1.0) ||
264 (rhenergy >
clean._eneThreshMod *
clean._doubleSpikeThresh &&
265 surroundingEnergyFraction <
clean._doubleSpikeS6S2 /
clean._fracThreshMod))) {
267 mask[mostEnergeticNeighbour] =
false;
Neighbours neighbours4() const
uint32_t cc[maxCellsPerHit]
double _doubleSpikeThresh
const std::unordered_map< std::string, int > _layerMap
PFLayer::Layer layer() const
rechit layer
void clean(const edm::Handle< reco::PFRecHitCollection > &input, std::vector< bool > &mask) override
static std::string const input
unsigned detId() const
rechit detId
double _singleSpikeThresh
SpikeAndDoubleSpikeCleaner(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
float phi() const
momentum azimuthal angle
Abs< T >::type abs(const T &t)
RhoEtaPhi const & positionREP() const
std::pair< int, edm::FunctionWithDict > OK
constexpr uint32_t rawId() const
get the raw id
float eta() const
momentum pseudorapidity
std::unordered_map< int, spike_cleaning > _thresholds
float energy() const
rechit energy
SpikeAndDoubleSpikeCleaner & operator=(const SpikeAndDoubleSpikeCleaner &)=delete
#define DEFINE_EDM_PLUGIN(factory, type, name)