CMS 3D CMS Logo

SpikeAndDoubleSpikeCleaner.cc
Go to the documentation of this file.
2 #include <cmath>
3 
4 namespace {
5  std::pair<double, double> dCrack(double phi, double eta) {
6  constexpr double oneOverCrystalSize = 1.0 / 0.0175;
7  constexpr double pi = M_PI;
8  constexpr double twopi = 2 * pi;
9  // the result below is from unrolling
10  // the lazy-eval loop in PFClusterAlgo::dCrack
11  constexpr double cPhi[18] = {2.97025,
12  2.621184149601134,
13  2.272118299202268,
14  1.9230524488034024,
15  1.5739865984045365,
16  1.2249207480056705,
17  0.8758548976068048,
18  0.5267890472079388,
19  0.1777231968090729,
20  -0.17134265358979306,
21  -0.520408503988659,
22  -0.8694743543875245,
23  -1.2185402047863905,
24  -1.5676060551852569,
25  -1.9166719055841224,
26  -2.265737755982988,
27  -2.6148036063818543,
28  -2.9638694567807207};
29  constexpr double cEta[9] = {
30  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};
31  // shift for eta < 0
32  constexpr double delta_cPhi = 0.00638;
33  // let's calculate dphi
34  double defi = 0;
35  if (eta < 0)
36  phi += delta_cPhi;
37  if (phi >= -pi && phi <= pi) {
38  //the problem of the extrema
39  if (phi < cPhi[17] || phi >= cPhi[0]) {
40  if (phi < 0)
41  phi += 2 * pi;
42  defi = std::min(std::abs(phi - cPhi[0]), std::abs(phi - cPhi[17] - twopi));
43  } else { // between these extrema
44  bool OK = false;
45  unsigned i = 16;
46  while (!OK) {
47  if (phi < cPhi[i]) {
48  defi = std::min(std::abs(phi - cPhi[i + 1]), std::abs(phi - cPhi[i]));
49  OK = true;
50  } else {
51  i -= 1;
52  }
53  } // end while
54  }
55  } else { // if there's a problem assume we're in a crack
56  defi = 0;
57  }
58  // let's calculate deta
59  double deta = 999.0;
60  for (const double etaGap : cEta) {
61  deta = std::min(deta, std::abs(eta - etaGap));
62  }
63  defi *= oneOverCrystalSize;
64  deta *= oneOverCrystalSize;
65  return std::make_pair(defi, deta);
66  }
67 } // namespace
68 
71  _layerMap({{"PS2", (int)PFLayer::PS2},
72  {"PS1", (int)PFLayer::PS1},
73  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
74  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
75  {"NONE", (int)PFLayer::NONE},
76  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
77  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
78  // hack to deal with ring1 in HO
79  {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
80  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
81  {"HF_EM", (int)PFLayer::HF_EM},
82  {"HF_HAD", (int)PFLayer::HF_HAD}}) {
83  const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("cleaningByDetector");
84  for (const auto& pset : thresholds) {
85  spike_cleaning info;
86  const std::string& det = pset.getParameter<std::string>("detector");
87  info._minS4S1_a = pset.getParameter<double>("minS4S1_a");
88  info._minS4S1_b = pset.getParameter<double>("minS4S1_b");
89  info._doubleSpikeS6S2 = pset.getParameter<double>("doubleSpikeS6S2");
90  info._eneThreshMod = pset.getParameter<double>("energyThresholdModifier");
91  info._fracThreshMod = pset.getParameter<double>("fractionThresholdModifier");
92  info._doubleSpikeThresh = pset.getParameter<double>("doubleSpikeThresh");
93  info._singleSpikeThresh = pset.getParameter<double>("singleSpikeThresh");
94  auto entry = _layerMap.find(det);
95  if (entry == _layerMap.end()) {
96  throw cms::Exception("InvalidDetectorLayer") << "Detector layer : " << det << " is not in the list of recognized"
97  << " detector layers!";
98  }
99  _thresholds.emplace(_layerMap.find(det)->second, info);
100  }
101 }
102 
104  //need to run over energy sorted rechits
105  auto const& hits = *input;
106  std::vector<unsigned> ordered_hits(hits.size());
107  for (unsigned i = 0; i < hits.size(); ++i)
108  ordered_hits[i] = i;
109  std::sort(ordered_hits.begin(), ordered_hits.end(), [&](unsigned i, unsigned j) {
110  return hits[i].energy() > hits[j].energy();
111  });
112 
113  for (const auto& idx : ordered_hits) {
114  const unsigned i = idx;
115  if (!mask[i])
116  continue; // don't need to re-mask things :-)
117  const reco::PFRecHit& rechit = hits[i];
118  int hitlayer = (int)rechit.layer();
119  if (hitlayer == PFLayer::HCAL_BARREL2 && std::abs(rechit.positionREP().eta()) > 0.34) {
120  hitlayer *= 100;
121  }
122  const spike_cleaning& clean = _thresholds.find(hitlayer)->second;
123  if (rechit.energy() < clean._singleSpikeThresh)
124  continue;
125 
126  //Fix needed for HF. Here, we find the (up to) five companion rechits
127  //to work in conjunction with the neighbours4() implementation below for the full HF surrounding energy
128  float compsumE = 0.0;
129  if ((hitlayer == PFLayer::HF_EM || hitlayer == PFLayer::HF_HAD)) {
130  int comp = 1;
131  if (hitlayer == PFLayer::HF_EM)
132  comp = 2;
133  const HcalDetId& detid = (HcalDetId)rechit.detId();
134  int heta = detid.ieta();
135  int hphi = detid.iphi();
136 
137  //At eta>39, phi segmentation changes
138  int predphi = 2;
139  if (std::abs(heta) > 39)
140  predphi = 4;
141 
142  int curphiL = hphi - predphi;
143  int curphiH = hphi + predphi;
144 
145  //HcalDetId valid phi range (0-72)
146  while (curphiL < 0)
147  curphiL += 72;
148  while (curphiH > 72)
149  curphiH -= 72;
150 
151  std::pair<std::vector<int>, std::vector<int>> phietas({heta, heta + 1, heta - 1, heta, heta},
152  {hphi, hphi, hphi, curphiL, curphiH});
153 
154  std::vector<uint32_t> rawDetIds;
155  for (unsigned in = 0; in < phietas.first.size(); in++) {
156  HcalDetId tempID(HcalForward, phietas.first[in], phietas.second[in], comp);
157  rawDetIds.push_back(tempID.rawId());
158  }
159 
160  for (const auto& jdx : ordered_hits) {
161  const unsigned j = jdx;
162  const reco::PFRecHit& matchrechit = hits[j];
163  for (const auto& iID : rawDetIds)
164  if (iID == matchrechit.detId())
165  compsumE += matchrechit.energy();
166  }
167  }
168  //End of fix needed for HF
169 
170  const double rhenergy = rechit.energy();
171  // single spike cleaning
172  auto const& neighbours4 = rechit.neighbours4();
173  double surroundingEnergy = compsumE;
174  for (auto k : neighbours4) {
175  if (!mask[k])
176  continue;
177  auto const& neighbour = hits[k];
178  const double sum = neighbour.energy(); //energyUp is just rechit energy?
179  surroundingEnergy += sum;
180  }
181  // wannaBeSeed.energyUp()/wannaBeSeed.energy() : 1.;
182  // Fraction 1 is the balance between the hit and its neighbours
183  // from both layers
184  const double fraction1 = surroundingEnergy / rhenergy;
185  // removed spurious comments from old pfcluster algo...
186  // look there if you want more history
187  const double f1Cut = (clean._minS4S1_a * std::log10(rechit.energy()) + clean._minS4S1_b);
188  if (fraction1 < f1Cut) {
189  const double eta = rechit.positionREP().eta();
190  const double aeta = std::abs(eta);
191  const double phi = rechit.positionREP().phi();
192  std::pair<double, double> dcr = dCrack(phi, eta);
193  const double dcrmin = (rechit.layer() == PFLayer::ECAL_BARREL ? std::min(dcr.first, dcr.second) : dcr.second);
194  if (aeta < 5.0 && ((aeta < 2.85 && dcrmin > 1.0) || (rhenergy > clean._eneThreshMod * clean._singleSpikeThresh &&
195  fraction1 < f1Cut / clean._fracThreshMod))) {
196  mask[i] = false;
197  }
198  } //if initial fraction cut (single spike)
199  // double spike removal
200  if (mask[i] && rhenergy > clean._doubleSpikeThresh) {
201  //Determine energy surrounding the seed and the most energetic neighbour
202  double surroundingEnergyi = 0.0;
203  double enmax = -999.0;
204  unsigned int mostEnergeticNeighbour = 0;
205  auto const& neighbours4i = rechit.neighbours4();
206  for (auto k : neighbours4i) {
207  if (!mask[k])
208  continue;
209  auto const& neighbour = hits[k];
210  const double nenergy = neighbour.energy();
211  surroundingEnergyi += nenergy;
212  if (nenergy > enmax) {
213  enmax = nenergy;
214  mostEnergeticNeighbour = k;
215  }
216  }
217  // is there an energetic neighbour
218  if (enmax > 0.0) {
219  double surroundingEnergyj = 0.0;
220  auto const& neighbours4j = hits[mostEnergeticNeighbour].neighbours4();
221  for (auto k : neighbours4j) {
222  //if( !mask[k] && k != i) continue; // leave out?
223  surroundingEnergyj += hits[k].energy();
224  }
225  // The energy surrounding the double spike candidate
226  const double surroundingEnergyFraction =
227  (surroundingEnergyi + surroundingEnergyj) / (rechit.energy() + hits[mostEnergeticNeighbour].energy()) - 1.;
228  if (surroundingEnergyFraction < clean._doubleSpikeS6S2) {
229  const double eta = rechit.positionREP().eta();
230  const double aeta = std::abs(eta);
231  const double phi = rechit.positionREP().phi();
232  std::pair<double, double> dcr = dCrack(phi, eta);
233  const double dcrmin = (rechit.layer() == PFLayer::ECAL_BARREL ? std::min(dcr.first, dcr.second) : dcr.second);
234  if (aeta < 5.0 && ((aeta < 2.85 && dcrmin > 1.0) ||
235  (rhenergy > clean._eneThreshMod * clean._doubleSpikeThresh &&
236  surroundingEnergyFraction < clean._doubleSpikeS6S2 / clean._fracThreshMod))) {
237  mask[i] = false;
238  mask[mostEnergeticNeighbour] = false;
239  }
240  }
241  } // was there an energetic neighbour ?
242  } // if double spike thresh
243  } // rechit loop
244 }
mps_fire.i
i
Definition: mps_fire.py:428
input
static const std::string input
Definition: EdmProvDump.cc:48
particleFlowZeroSuppressionECAL_cff.thresholds
thresholds
Definition: particleFlowZeroSuppressionECAL_cff.py:31
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
reco::PFRecHit::energy
float energy() const
rechit energy
Definition: PFRecHit.h:99
HcalDetId::iphi
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
min
T min(T a, T b)
Definition: MathUtil.h:58
mps_splice.entry
entry
Definition: mps_splice.py:68
PFLayer::HCAL_ENDCAP
Definition: PFLayer.h:37
SpikeAndDoubleSpikeCleaner::_thresholds
std::unordered_map< int, spike_cleaning > _thresholds
Definition: SpikeAndDoubleSpikeCleaner.h:29
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition: FWCollectionSummaryWidget.cc:153
edm::Handle
Definition: AssociativeIterator.h:50
reco::OK
std::pair< int, edm::FunctionWithDict > OK
Definition: findMethod.cc:126
AlCaHLTBitMon_QueryRunRegistry.comp
comp
Definition: AlCaHLTBitMon_QueryRunRegistry.py:249
PFLayer::ECAL_BARREL
Definition: PFLayer.h:33
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
SpikeAndDoubleSpikeCleaner.h
PFLayer::PS1
Definition: PFLayer.h:31
PFLayer::HCAL_BARREL2
Definition: PFLayer.h:36
PFLayer::HF_EM
Definition: PFLayer.h:38
PVValHelper::eta
Definition: PVValidationHelpers.h:69
PFLayer::HCAL_BARREL1
Definition: PFLayer.h:35
PFLayer::NONE
Definition: PFLayer.h:34
dqmdumpme.k
k
Definition: dqmdumpme.py:60
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PFLayer::HF_HAD
Definition: PFLayer.h:39
HcalDetId::ieta
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
SpikeAndDoubleSpikeCleaner::clean
void clean(const edm::Handle< reco::PFRecHitCollection > &input, std::vector< bool > &mask) override
Definition: SpikeAndDoubleSpikeCleaner.cc:103
edm::ParameterSet
Definition: ParameterSet.h:47
PVValHelper::phi
Definition: PVValidationHelpers.h:68
SpikeAndDoubleSpikeCleaner::SpikeAndDoubleSpikeCleaner
SpikeAndDoubleSpikeCleaner(const edm::ParameterSet &conf)
Definition: SpikeAndDoubleSpikeCleaner.cc:69
SpikeAndDoubleSpikeCleaner::spike_cleaning
Definition: SpikeAndDoubleSpikeCleaner.h:11
recoMuon::in
Definition: RecoMuonEnumerators.h:6
HcalDetId
Definition: HcalDetId.h:12
createfilelist.int
int
Definition: createfilelist.py:10
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
RhoEtaPhi::eta
float eta() const
momentum pseudorapidity
Definition: PtEtaPhiMass.h:52
reco::PFRecHit::neighbours4
Neighbours neighbours4() const
Definition: PFRecHit.h:81
reco::PFRecHit::detId
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
HcalForward
Definition: HcalAssistant.h:36
DDAxes::phi
RhoEtaPhi::phi
float phi() const
momentum azimuthal angle
Definition: PtEtaPhiMass.h:54
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Exception
Definition: hltDiff.cc:246
reco::PFRecHit
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
reco::PFRecHit::layer
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
RecHitTopologicalCleanerBase
Definition: RecHitTopologicalCleanerBase.h:12
pi
const Double_t pi
Definition: trackSplitPlot.h:36
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
PFLayer::ECAL_ENDCAP
Definition: PFLayer.h:32
reco::PFRecHit::positionREP
RhoEtaPhi const & positionREP() const
Definition: PFRecHit.h:119
PFLayer::PS2
Definition: PFLayer.h:30
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27