10#include <klm/modules/KLMClustersReconstructor/KLMClustersReconstructorModule.h>
13#include <klm/dataobjects/KLMElementNumbers.h>
14#include <klm/dataobjects/bklm/BKLMElementNumbers.h>
15#include <klm/dataobjects/eklm/EKLMElementNumbers.h>
18#include <Math/VectorUtil.h>
23#include <unordered_set>
35 setDescription(
"Unified BKLM/EKLM module for the reconstruction of KLMClusters.");
40 "Vertex position mode: 'FullAverage', 'FirstLayer', 'FirstTwoLayers', "
41 "or 'SuccessiveTwoLayers'.",
42 std::string(
"FirstTwoLayers"));
44 "Clusterization mode ('AnyHit' or 'FirstHit').",
45 std::string(
"AnyHit"));
47 "If true, iteratively remove angular outliers after clustering "
48 "(dropped hits re-enter the pool). If false, behavior matches the "
49 "original algorithm.",
52 "Used only if removeOutlierHits: floor (minimum) angular threshold (rad). "
53 "Effective cut is max(OutlierTrimAngle, k * MAD(residuals)) with fixed k=3.",
75 B2FATAL(
"Incorrect PositionMode argument.");
81 B2FATAL(
"Incorrect ClusterMode argument.");
93static KLMHit2d* hitWithMinR(
const std::vector<KLMHit2d*>& hits)
97 for (
size_t k = 1; k < hits.size(); k++) {
98 double r = hits[k]->getPosition().R();
108static ROOT::Math::XYZVector unitDirection(
const ROOT::Math::XYZVector& v)
112 return ROOT::Math::XYZVector{0, 0, 0};
113 return v * (1.0 / r);
118static ROOT::Math::XYZVector referenceDirection(
const std::vector<KLMHit2d*>& inliers,
119 bool seedByInnermost)
122 return unitDirection(hitWithMinR(inliers)->
getPosition());
123 ROOT::Math::XYZVector sum{0, 0, 0};
125 sum = sum + unitDirection(h->getPosition());
127 return unitDirection(hitWithMinR(inliers)->
getPosition());
128 return sum * (1.0 / sum.R());
132static double adaptiveThreshold(
const ROOT::Math::XYZVector& ref,
133 const std::vector<KLMHit2d*>& hits,
134 double floorAngle,
double madFactor)
136 std::vector<double> residuals;
137 residuals.reserve(hits.size());
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);
153 std::vector<KLMHit2d*>& poolHits)
157 constexpr double kMADFactor = 3.0;
158 constexpr int kMaxOutlierIterations = 5;
159 constexpr double kMinInlierFraction = 0.3;
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())));
168 std::vector<KLMHit2d*> inliers = before;
169 for (
int iter = 0; iter < kMaxOutlierIterations; iter++) {
170 const ROOT::Math::XYZVector ref =
171 referenceDirection(inliers, iter == 0);
172 const double threshold =
174 std::vector<KLMHit2d*> next;
175 next.reserve(before.size());
177 if (ROOT::Math::VectorUtil::Angle(h->getPosition(), ref) <= threshold)
180 if (next.size() < minInliers)
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());
187 inliers = std::move(next);
191 inliers = std::move(next);
194 if (inliers.size() < minInliers)
196 if (inliers.size() == before.size())
201 std::unordered_set<KLMHit2d*> inlierSet(inliers.begin(), inliers.end());
202 std::vector<KLMHit2d*> outliers;
203 outliers.reserve(before.size() - inliers.size());
205 if (inlierSet.find(h) == inlierSet.end())
206 outliers.push_back(h);
208 clusterHits.swap(inliers);
209 poolHits.insert(poolHits.end(), outliers.begin(), outliers.end());
210 sort(poolHits.begin(), poolHits.end(), compareDistance);
216 int i, nLayers, innermostLayer, nHits;
219 int* layerHitsBKLM, *layerHitsEKLM;
222 std::vector<KLMHit2d*> klmHit2ds, klmClusterHits;
223 std::vector<KLMHit2d*>::iterator it, it0, it2;
225 layerHitsBKLM =
new int[nLayersBKLM];
226 layerHitsEKLM =
new int[nLayersEKLM];
229 for (i = 0; i < nHits; i++) {
235 sort(klmHit2ds.begin(), klmHit2ds.end(), compareDistance);
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();
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);
258 if (ROOT::Math::VectorUtil::Angle(
259 (*it)->getPosition(), (*it2)->getPosition()) <
261 klmClusterHits.push_back(*it);
262 it = klmHit2ds.erase(it);
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;
279 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
280 if (minTime < 0 || (*it)->getTime() < minTime)
281 minTime = (*it)->getTime();
283 layerHitsBKLM[(*it)->getLayer() - 1]++;
285 layerHitsEKLM[(*it)->getLayer() - 1]++;
287 it0 = klmClusterHits.begin();
290 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
291 clusterPosition = clusterPosition + (*it)->getPosition();
299 selSubdet[0] = (*it0)->getSubdetector();
300 selLayer[0] = (*it0)->getLayer();
303 for (i = 0; i < nLayersBKLM && selCount < 2; i++) {
304 if (layerHitsBKLM[i] > 0) {
306 selLayer[selCount] = i + 1;
310 for (i = 0; i < nLayersEKLM && selCount < 2; i++) {
311 if (layerHitsEKLM[i] > 0) {
313 selLayer[selCount] = i + 1;
318 selSubdet[0] = (*it0)->getSubdetector();
319 selLayer[0] = (*it0)->getLayer();
323 bool foundPair =
false;
324 for (i = 0; i < nLayersBKLM - 1; i++) {
325 if (layerHitsBKLM[i] > 0 && layerHitsBKLM[i + 1] > 0) {
336 for (i = 0; i < nLayersEKLM - 1; i++) {
337 if (layerHitsEKLM[i] > 0 && layerHitsEKLM[i + 1] > 0) {
349 selSubdet[0] = (*it0)->getSubdetector();
350 selLayer[0] = (*it0)->getLayer();
354 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it) {
355 int sd = (*it)->getSubdetector();
356 int ly = (*it)->getLayer();
358 for (i = 0; i < selCount; i++) {
359 if (sd == selSubdet[i] && ly == selLayer[i]) {
365 clusterPosition = clusterPosition + (*it)->getPosition();
370 B2FATAL(
"KLMClustersReconstructor: no hits for vertex position.");
372 clusterPosition = clusterPosition * (1.0 / nHits);
376 for (i = 0; i < nLayersBKLM; i++) {
377 if (layerHitsBKLM[i] > 0) {
379 if (innermostLayer < 0)
380 innermostLayer = i + 1;
383 for (i = 0; i < nLayersEKLM; i++) {
384 if (layerHitsEKLM[i] > 0) {
386 if (innermostLayer < 0)
387 innermostLayer = i + 1;
396 p = klmClusterHits.size() * 0.215;
406 clusterPosition.X(), clusterPosition.Y(), clusterPosition.Z(), minTime, nLayers,
408 for (it = klmClusterHits.begin(); it != klmClusterHits.end(); ++it)
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();
420 nDigits += klmdigits.size();
425 delete[] layerHitsBKLM;
426 delete[] layerHitsEKLM;
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
Store one reconstructed BKLM 1D hit as a ROOT object.
static constexpr int getMaximalLayerNumber()
Get maximal layer number.
void setKLMnDigits(int nDigits)
Set number of KLM digits in the cluster.
~KLMClustersReconstructorModule()
Destructor.
enum ClusterMode m_ClusterMode
Clusterization mode.
void initialize() override
Initializer.
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.
std::string m_ClusterModeString
Clusterization mode.
@ c_FullAverage
Full average.
@ c_SuccessiveTwoLayers
Innermost adjacent layer pair with hits; else FirstLayer.
@ c_FirstLayer
Innermost hit layer only.
@ c_FirstTwoLayers
Two innermost layers with hits (BKLM then EKLM).
KLMClustersReconstructorModule()
Constructor.
void beginRun() override
Called when entering a new run.
double m_ClusteringAngle
Clustering angle.
void applyOutlierRemoval(std::vector< KLMHit2d * > &clusterHits, std::vector< KLMHit2d * > &poolHits)
Optional post-cluster hit filtering.
StoreArray< KLMCluster > m_KLMClusters
KLM clusters.
@ c_FirstHit
Angle from first hit.
@ c_AnyHit
Angle from any hit.
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).
ROOT::Math::XYZVector getPosition() const
Get hit global position.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
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 ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
ExpRunEvt getPosition(const std::vector< Evt > &events, double tEdge)
Get the exp-run-evt number from the event time [hours].
Abstract base class for different kinds of events.