Belle II Software prerelease-11-00-00a
KLMClustersReconstructorModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <klm/modules/KLMClustersReconstructor/KLMClustersReconstructorModule.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/KLMElementNumbers.h>
14#include <klm/dataobjects/bklm/BKLMElementNumbers.h>
15#include <klm/dataobjects/eklm/EKLMElementNumbers.h>
16
17/* ROOT headers. */
18#include <Math/VectorUtil.h>
19
20/* C++ headers. */
21#include <algorithm>
22#include <cmath>
23#include <unordered_set>
24
25using namespace Belle2;
26
27REG_MODULE(KLMClustersReconstructor);
28
34{
35 setDescription("Unified BKLM/EKLM module for the reconstruction of KLMClusters.");
37 addParam("ClusteringAngle", m_ClusteringAngle, "Clustering angle (rad).",
38 0.26);
39 addParam("PositionMode", m_PositionModeString,
40 "Vertex position mode: 'FullAverage', 'FirstLayer', 'FirstTwoLayers', "
41 "or 'SuccessiveTwoLayers'.",
42 std::string("FirstTwoLayers"));
43 addParam("ClusterMode", m_ClusterModeString,
44 "Clusterization mode ('AnyHit' or 'FirstHit').",
45 std::string("AnyHit"));
46 addParam("RemoveOutlierHits", m_RemoveOutlierHits,
47 "If true, iteratively remove angular outliers after clustering "
48 "(dropped hits re-enter the pool). If false, behavior matches the "
49 "original algorithm.",
50 false);
51 addParam("OutlierTrimAngle", m_OutlierTrimAngle,
52 "Used only if removeOutlierHits: floor (minimum) angular threshold (rad). "
53 "Effective cut is max(OutlierTrimAngle, k * MAD(residuals)) with fixed k=3.",
54 0.18);
55}
56
60
62{
63 m_KLMClusters.registerInDataStore();
64 m_Hit2ds.isRequired();
65 m_KLMClusters.registerRelationTo(m_Hit2ds);
66 if (m_PositionModeString == "FullAverage")
68 else if (m_PositionModeString == "FirstLayer")
70 else if (m_PositionModeString == "FirstTwoLayers")
72 else if (m_PositionModeString == "SuccessiveTwoLayers")
74 else
75 B2FATAL("Incorrect PositionMode argument.");
76 if (m_ClusterModeString == "AnyHit")
78 else if (m_ClusterModeString == "FirstHit")
80 else
81 B2FATAL("Incorrect ClusterMode argument.");
82}
83
87
88static bool compareDistance(const KLMHit2d* hit1, const KLMHit2d* hit2)
89{
90 return hit1->getPosition().R() < hit2->getPosition().R();
91}
92
93static KLMHit2d* hitWithMinR(const std::vector<KLMHit2d*>& hits)
94{
95 KLMHit2d* best = hits[0];
96 double bestR = best->getPosition().R();
97 for (size_t k = 1; k < hits.size(); k++) {
98 double r = hits[k]->getPosition().R();
99 if (r < bestR) {
100 bestR = r;
101 best = hits[k];
102 }
103 }
104 return best;
105}
106
107/* Unit vector from origin to the hit position (zero if hit is at origin). */
108static ROOT::Math::XYZVector unitDirection(const ROOT::Math::XYZVector& v)
109{
110 double r = v.R();
111 if (r <= 0)
112 return ROOT::Math::XYZVector{0, 0, 0};
113 return v * (1.0 / r);
114}
115
116/* Reference direction: innermost-R hit direction on first iteration,
117 * otherwise the normalized sum of unit directions of the current inliers. */
118static ROOT::Math::XYZVector referenceDirection(const std::vector<KLMHit2d*>& inliers,
119 bool seedByInnermost)
120{
121 if (seedByInnermost)
122 return unitDirection(hitWithMinR(inliers)->getPosition());
123 ROOT::Math::XYZVector sum{0, 0, 0};
124 for (KLMHit2d* h : inliers)
125 sum = sum + unitDirection(h->getPosition());
126 if (sum.R() <= 0)
127 return unitDirection(hitWithMinR(inliers)->getPosition());
128 return sum * (1.0 / sum.R());
129}
130
131/* Adaptive angular threshold: max(floor, k * MAD of angular residuals to ref). */
132static double adaptiveThreshold(const ROOT::Math::XYZVector& ref,
133 const std::vector<KLMHit2d*>& hits,
134 double floorAngle, double madFactor)
135{
136 std::vector<double> residuals;
137 residuals.reserve(hits.size());
138 for (KLMHit2d* h : hits)
139 residuals.push_back(ROOT::Math::VectorUtil::Angle(h->getPosition(), ref));
140 std::vector<double> sorted = residuals;
141 std::sort(sorted.begin(), sorted.end());
142 double median = sorted[sorted.size() / 2];
143 std::vector<double> absDev;
144 absDev.reserve(residuals.size());
145 for (double r : residuals)
146 absDev.push_back(std::fabs(r - median));
147 std::sort(absDev.begin(), absDev.end());
148 double mad = absDev[absDev.size() / 2];
149 return std::max(floorAngle, madFactor * mad);
150}
151
152void KLMClustersReconstructorModule::applyOutlierRemoval(std::vector<KLMHit2d*>& clusterHits,
153 std::vector<KLMHit2d*>& poolHits)
154{
155 /* Fixed tuning: standard k=3 for MAD scale; iteration cap (fixed-point exits early);
156 * min inlier fraction 0.3 avoids over-trimming a cluster to fewer than 30% the hits. */
157 constexpr double kMADFactor = 3.0;
158 constexpr int kMaxOutlierIterations = 5;
159 constexpr double kMinInlierFraction = 0.3;
160
161 if (!m_RemoveOutlierHits || clusterHits.size() <= 1)
162 return;
163
164 const std::vector<KLMHit2d*> before = clusterHits;
165 const std::size_t minInliers = std::max<std::size_t>(
166 2, static_cast<std::size_t>(std::ceil(kMinInlierFraction * before.size())));
167
168 std::vector<KLMHit2d*> inliers = before;
169 for (int iter = 0; iter < kMaxOutlierIterations; iter++) {
170 const ROOT::Math::XYZVector ref =
171 referenceDirection(inliers, /*seedByInnermost=*/iter == 0);
172 const double threshold =
173 adaptiveThreshold(ref, before, m_OutlierTrimAngle, kMADFactor);
174 std::vector<KLMHit2d*> next;
175 next.reserve(before.size());
176 for (KLMHit2d* h : before) {
177 if (ROOT::Math::VectorUtil::Angle(h->getPosition(), ref) <= threshold)
178 next.push_back(h);
179 }
180 if (next.size() < minInliers)
181 break;
182 if (next.size() == inliers.size()) {
183 std::vector<KLMHit2d*> a = next, b = inliers;
184 std::sort(a.begin(), a.end());
185 std::sort(b.begin(), b.end());
186 if (a == b) {
187 inliers = std::move(next);
188 break;
189 }
190 }
191 inliers = std::move(next);
192 }
193
194 if (inliers.size() < minInliers)
195 return;
196 if (inliers.size() == before.size())
197 return;
198 if (inliers.empty())
199 return;
200
201 std::unordered_set<KLMHit2d*> inlierSet(inliers.begin(), inliers.end());
202 std::vector<KLMHit2d*> outliers;
203 outliers.reserve(before.size() - inliers.size());
204 for (KLMHit2d* h : before) {
205 if (inlierSet.find(h) == inlierSet.end())
206 outliers.push_back(h);
207 }
208 clusterHits.swap(inliers);
209 poolHits.insert(poolHits.end(), outliers.begin(), outliers.end());
210 sort(poolHits.begin(), poolHits.end(), compareDistance);
211}
212
214{
215 //static double mass = Const::Klong.getMass();
216 int i, nLayers, innermostLayer, nHits;
219 int* layerHitsBKLM, *layerHitsEKLM;
220 float minTime = -1;
221 double p;//, v;
222 std::vector<KLMHit2d*> klmHit2ds, klmClusterHits;
223 std::vector<KLMHit2d*>::iterator it, it0, it2;
224 KLMCluster* klmCluster;
225 layerHitsBKLM = new int[nLayersBKLM];
226 layerHitsEKLM = new int[nLayersEKLM];
227 /* Fill vector of 2d hits. */
228 nHits = m_Hit2ds.getEntries();
229 for (i = 0; i < nHits; i++) {
230 if (m_Hit2ds[i]->isOutOfTime())
231 continue;
232 klmHit2ds.push_back(m_Hit2ds[i]);
233 }
234 /* Sort by the distance from center. */
235 sort(klmHit2ds.begin(), klmHit2ds.end(), compareDistance);
236 /* Clustering. */
237 while (klmHit2ds.size() > 0) {
238 klmClusterHits.clear();
239 it = klmHit2ds.begin();
240 klmClusterHits.push_back(*it);
241 it = klmHit2ds.erase(it);
242 while (it != klmHit2ds.end()) {
243 it2 = klmClusterHits.begin();
244 switch (m_ClusterMode) {
245 case c_AnyHit:
246 while (it2 != klmClusterHits.end()) {
247 if (ROOT::Math::VectorUtil::Angle(
248 (*it)->getPosition(), (*it2)->getPosition()) <
250 klmClusterHits.push_back(*it);
251 it = klmHit2ds.erase(it);
252 goto clusterFound;
253 } else
254 ++it2;
255 }
256 break;
257 case c_FirstHit:
258 if (ROOT::Math::VectorUtil::Angle(
259 (*it)->getPosition(), (*it2)->getPosition()) <
261 klmClusterHits.push_back(*it);
262 it = klmHit2ds.erase(it);
263 goto clusterFound;
264 }
265 break;
266 }
267 ++it;
268clusterFound:;
269 }
271 applyOutlierRemoval(klmClusterHits, klmHit2ds);
272 ROOT::Math::XYZVector clusterPosition{0, 0, 0};
273 for (i = 0; i < nLayersBKLM; i++)
274 layerHitsBKLM[i] = 0;
275 for (i = 0; i < nLayersEKLM; i++)
276 layerHitsEKLM[i] = 0;
277 minTime = -1;
278 /* Minimal time and per-layer hit counts. */
279 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
280 if (minTime < 0 || (*it)->getTime() < minTime)
281 minTime = (*it)->getTime();
282 if ((*it)->getSubdetector() == KLMElementNumbers::c_BKLM)
283 layerHitsBKLM[(*it)->getLayer() - 1]++;
284 else
285 layerHitsEKLM[(*it)->getLayer() - 1]++;
286 }
287 it0 = klmClusterHits.begin();
288 nHits = 0;
290 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
291 clusterPosition = clusterPosition + (*it)->getPosition();
292 nHits++;
293 }
294 } else {
295 int selSubdet[2];
296 int selLayer[2];
297 int selCount = 0;
299 selSubdet[0] = (*it0)->getSubdetector();
300 selLayer[0] = (*it0)->getLayer();
301 selCount = 1;
302 } else if (m_PositionMode == c_FirstTwoLayers) {
303 for (i = 0; i < nLayersBKLM && selCount < 2; i++) {
304 if (layerHitsBKLM[i] > 0) {
305 selSubdet[selCount] = KLMElementNumbers::c_BKLM;
306 selLayer[selCount] = i + 1;
307 selCount++;
308 }
309 }
310 for (i = 0; i < nLayersEKLM && selCount < 2; i++) {
311 if (layerHitsEKLM[i] > 0) {
312 selSubdet[selCount] = KLMElementNumbers::c_EKLM;
313 selLayer[selCount] = i + 1;
314 selCount++;
315 }
316 }
317 if (selCount == 0) {
318 selSubdet[0] = (*it0)->getSubdetector();
319 selLayer[0] = (*it0)->getLayer();
320 selCount = 1;
321 }
323 bool foundPair = false;
324 for (i = 0; i < nLayersBKLM - 1; i++) {
325 if (layerHitsBKLM[i] > 0 && layerHitsBKLM[i + 1] > 0) {
326 selSubdet[0] = KLMElementNumbers::c_BKLM;
327 selLayer[0] = i + 1;
328 selSubdet[1] = KLMElementNumbers::c_BKLM;
329 selLayer[1] = i + 2;
330 selCount = 2;
331 foundPair = true;
332 break;
333 }
334 }
335 if (!foundPair) {
336 for (i = 0; i < nLayersEKLM - 1; i++) {
337 if (layerHitsEKLM[i] > 0 && layerHitsEKLM[i + 1] > 0) {
338 selSubdet[0] = KLMElementNumbers::c_EKLM;
339 selLayer[0] = i + 1;
340 selSubdet[1] = KLMElementNumbers::c_EKLM;
341 selLayer[1] = i + 2;
342 selCount = 2;
343 foundPair = true;
344 break;
345 }
346 }
347 }
348 if (!foundPair) {
349 selSubdet[0] = (*it0)->getSubdetector();
350 selLayer[0] = (*it0)->getLayer();
351 selCount = 1;
352 }
353 }
354 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
355 int sd = (*it)->getSubdetector();
356 int ly = (*it)->getLayer();
357 bool use = false;
358 for (i = 0; i < selCount; i++) {
359 if (sd == selSubdet[i] && ly == selLayer[i]) {
360 use = true;
361 break;
362 }
363 }
364 if (use) {
365 clusterPosition = clusterPosition + (*it)->getPosition();
366 nHits++;
367 }
368 }
369 if (nHits == 0)
370 B2FATAL("KLMClustersReconstructor: no hits for vertex position.");
371 }
372 clusterPosition = clusterPosition * (1.0 / nHits);
373 /* Find innermost layer. */
374 nLayers = 0;
375 innermostLayer = -1;
376 for (i = 0; i < nLayersBKLM; i++) {
377 if (layerHitsBKLM[i] > 0) {
378 nLayers++;
379 if (innermostLayer < 0)
380 innermostLayer = i + 1;
381 }
382 }
383 for (i = 0; i < nLayersEKLM; i++) {
384 if (layerHitsEKLM[i] > 0) {
385 nLayers++;
386 if (innermostLayer < 0)
387 innermostLayer = i + 1;
388 }
389 }
390 /* Calculate energy. */
391 //if (it0->inBKLM()) {
392 /*
393 * TODO: The constant is from BKLM K0L reconstructor,
394 * it must be recalculated.
395 */
396 p = klmClusterHits.size() * 0.215;
397 /* FIXME: Reimplement time calculation after completion of time calibration.
398 } else {
399 v = clusterPosition.R() / minTime / Const::speedOfLight;
400 if (v < 0.999999)
401 p = mass * v / sqrt(1.0 - v * v);
402 else
403 p = 0;
404 }*/
405 klmCluster = m_KLMClusters.appendNew(
406 clusterPosition.X(), clusterPosition.Y(), clusterPosition.Z(), minTime, nLayers,
407 innermostLayer, p);
408 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it)
409 klmCluster->addRelationTo(*it);
410
411 /* number of KLM digits in the cluster (BKLM via BKLMHit1d and EKLM directly). */
412 int nDigits = 0;
413 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
414 auto bklmhit1ds = (*it)->getRelationsTo<BKLMHit1d>();
415 for (const auto& bklmhit1d : bklmhit1ds) {
416 auto klmdigits = bklmhit1d.getRelationsTo<KLMDigit>();
417 nDigits += klmdigits.size();
418 }
419 auto klmdigits = (*it)->getRelationsTo<KLMDigit>();
420 nDigits += klmdigits.size();
421 }
422 klmCluster->setKLMnDigits(nDigits);
423 }
424
425 delete[] layerHitsBKLM;
426 delete[] layerHitsEKLM;
427}
428
432
436
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition BKLMHit1d.h:30
static constexpr int getMaximalLayerNumber()
Get maximal layer number.
KLM cluster data.
Definition KLMCluster.h:29
void setKLMnDigits(int nDigits)
Set number of KLM digits in the cluster.
Definition KLMCluster.h:276
enum ClusterMode m_ClusterMode
Clusterization mode.
void event() override
This method is called for each event.
void endRun() override
This method is called if the current run ends.
void terminate() override
This method is called at the end of the event processing.
double m_OutlierTrimAngle
Floor (minimum) angular threshold in rad; effective cut is max(floor, k*MAD).
std::string m_PositionModeString
Vertex position calculation mode.
enum PositionMode m_PositionMode
Vertex position calculation mode.
@ c_SuccessiveTwoLayers
Innermost adjacent layer pair with hits; else FirstLayer.
@ c_FirstTwoLayers
Two innermost layers with hits (BKLM then EKLM).
void beginRun() override
Called when entering a new run.
void applyOutlierRemoval(std::vector< KLMHit2d * > &clusterHits, std::vector< KLMHit2d * > &poolHits)
Optional post-cluster hit filtering.
StoreArray< KLMCluster > m_KLMClusters
KLM clusters.
bool m_RemoveOutlierHits
If true, drop angular outliers after clustering and re-queue them for other clusters.
StoreArray< KLMHit2d > m_Hit2ds
Two-dimensional hits.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition KLMDigit.h:29
KLM 2d hit.
Definition KLMHit2d.h:33
ROOT::Math::XYZVector getPosition() const
Get hit global position.
Definition KLMHit2d.h:315
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition Module.cc:208
Module()
Constructor.
Definition Module.cc:30
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition Module.h:80
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
ExpRunEvt getPosition(const std::vector< Evt > &events, double tEdge)
Get the exp-run-evt number from the event time [hours].
Definition Splitter.h:341
Abstract base class for different kinds of events.