File indexing completed on 2026-06-04 08:11:44
0001
0002
0003
0004
0005
0006
0007
0008
0009 #pragma once
0010
0011 #include "Acts/Definitions/TrackParametrization.hpp"
0012 #include "Acts/Surfaces/CylinderSurface.hpp"
0013 #include "Acts/TrackFitting/GsfComponent.hpp"
0014 #include "Acts/TrackFitting/GsfOptions.hpp"
0015 #include "Acts/Utilities/detail/periodic.hpp"
0016
0017 #include <iosfwd>
0018 #include <stdexcept>
0019 #include <tuple>
0020
0021 namespace Acts::detail::Gsf {
0022
0023
0024
0025
0026
0027
0028
0029
0030 std::tuple<BoundVector, BoundMatrix> mergeGaussianMixture(
0031 std::span<const GsfComponent> mixture, const Surface &surface,
0032 ComponentMergeMethod method);
0033
0034
0035
0036
0037
0038
0039
0040 GsfComponent mergeTwoComponents(const GsfComponent &a, const GsfComponent &b,
0041 const Surface &surface);
0042
0043
0044 template <BoundIndices Idx>
0045 struct CyclicAngle {
0046 constexpr static BoundIndices idx = Idx;
0047 constexpr static double constant = 1.0;
0048 };
0049
0050 template <BoundIndices Idx>
0051 struct CyclicRadiusAngle {
0052 constexpr static BoundIndices idx = Idx;
0053 double constant = 1.0;
0054 };
0055
0056
0057 template <Surface::SurfaceType type_t>
0058 struct AngleDescription {
0059 using Desc = std::tuple<CyclicAngle<eBoundPhi>>;
0060 };
0061
0062 template <>
0063 struct AngleDescription<Surface::Disc> {
0064 using Desc = std::tuple<CyclicAngle<eBoundLoc1>, CyclicAngle<eBoundPhi>>;
0065 };
0066
0067 template <>
0068 struct AngleDescription<Surface::Cylinder> {
0069 using Desc =
0070 std::tuple<CyclicRadiusAngle<eBoundLoc0>, CyclicAngle<eBoundPhi>>;
0071 };
0072
0073
0074
0075
0076
0077
0078 template <typename Callable>
0079 auto angleDescriptionSwitch(const Surface &surface, Callable &&callable) {
0080 switch (surface.type()) {
0081 case Surface::Cylinder: {
0082 AngleDescription<Surface::Cylinder>::Desc desc{};
0083 const auto &bounds =
0084 static_cast<const CylinderSurface &>(surface).bounds();
0085 std::get<0>(desc).constant = bounds.get(CylinderBounds::eR);
0086 return callable(desc);
0087 }
0088 case Surface::Disc: {
0089 AngleDescription<Surface::Disc>::Desc desc{};
0090 return callable(desc);
0091 }
0092 default: {
0093 AngleDescription<Surface::Plane>::Desc desc{};
0094 return callable(desc);
0095 }
0096 }
0097 }
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107 template <typename component_range_t, typename projector_t,
0108 typename angle_desc_t>
0109 BoundVector mergeGaussianMixtureMean(const component_range_t &cmps,
0110 const projector_t &projector,
0111 const angle_desc_t &angleDesc) {
0112
0113 if (cmps.size() == 1) {
0114 const auto [weight_l, pars_l] = projector(cmps.front());
0115 return pars_l;
0116 }
0117
0118
0119
0120
0121
0122 using CVec = Eigen::Matrix<std::complex<double>, eBoundSize, 1>;
0123 CVec cMean = CVec::Zero();
0124 double sumOfWeights = 0;
0125
0126 for (const auto &cmp : cmps) {
0127 const auto [weight_l, pars_l] = projector(cmp);
0128 CVec cPars_l = pars_l;
0129
0130 const auto setPolar = [&](const auto &desc) {
0131 cPars_l[desc.idx] = std::polar(1.0, pars_l[desc.idx] / desc.constant);
0132 };
0133 std::apply([&](auto... dsc) { (setPolar(dsc), ...); }, angleDesc);
0134
0135 sumOfWeights += weight_l;
0136 cMean += weight_l * cPars_l;
0137 }
0138
0139 cMean /= sumOfWeights;
0140
0141 BoundVector mean = cMean.real();
0142
0143 const auto getArg = [&](const auto &desc) {
0144 mean[desc.idx] = desc.constant * std::arg(cMean[desc.idx]);
0145 };
0146 std::apply([&](auto... dsc) { (getArg(dsc), ...); }, angleDesc);
0147
0148 return mean;
0149 }
0150
0151
0152
0153
0154
0155
0156
0157
0158 template <typename component_range_t, typename projector_t>
0159 BoundVector mergeGaussianMixtureMaxWeight(const component_range_t &cmps,
0160 const projector_t &projector) {
0161 const auto maxWeightIt =
0162 std::ranges::max_element(cmps, {}, [&](const auto &cmp) {
0163 const auto [weight_l, pars_l] = projector(cmp);
0164 return weight_l;
0165 });
0166 const auto [weight_max, pars_max] = projector(*maxWeightIt);
0167 return pars_max;
0168 }
0169
0170
0171
0172
0173
0174
0175
0176
0177
0178 template <typename component_range_t, typename projector_t,
0179 typename angle_desc_t>
0180 BoundMatrix mergeGaussianMixtureCov(const component_range_t &cmps,
0181 const projector_t &projector,
0182 const BoundVector &mean,
0183 const angle_desc_t &angleDesc) {
0184 if (cmps.size() == 1) {
0185 const auto [weight_l, pars_l, cov_l] = projector(cmps.front());
0186 return cov_l;
0187 }
0188
0189 BoundMatrix cov = BoundMatrix::Zero();
0190 double sumOfWeights = 0;
0191
0192 for (const auto &cmp : cmps) {
0193 const auto [weight_l, pars_l, cov_l] = projector(cmp);
0194
0195 cov += weight_l * cov_l;
0196
0197 BoundVector diff = pars_l - mean;
0198
0199
0200 const auto handleCyclicCov = [&l = pars_l, &m = mean,
0201 &diff = diff](const auto &desc) {
0202 diff[desc.idx] = difference_periodic(l[desc.idx] / desc.constant,
0203 m[desc.idx] / desc.constant,
0204 2 * std::numbers::pi) *
0205 desc.constant;
0206 };
0207 std::apply([&](auto... dsc) { (handleCyclicCov(dsc), ...); }, angleDesc);
0208
0209 sumOfWeights += weight_l;
0210 cov += weight_l * diff * diff.transpose();
0211 }
0212
0213 cov /= sumOfWeights;
0214
0215 return cov;
0216 }
0217
0218
0219
0220
0221
0222
0223
0224
0225
0226
0227 template <typename component_range_t, typename projector_t,
0228 typename angle_desc_t>
0229 std::tuple<BoundVector, BoundMatrix> mergeGaussianMixtureMeanCov(
0230 const component_range_t &cmps, const projector_t &projector,
0231 const angle_desc_t &angleDesc) {
0232 if (cmps.size() == 1) {
0233 const auto [weight_l, pars_l, cov_l] = projector(cmps.front());
0234 return {pars_l, cov_l};
0235 }
0236
0237 const BoundVector mean = mergeGaussianMixtureMean(
0238 cmps,
0239 [&](const auto &c) {
0240 const auto [weight_l, pars_l, cov_l] = projector(c);
0241 return std::tuple(weight_l, pars_l);
0242 },
0243 angleDesc);
0244
0245
0246 const BoundMatrix cov =
0247 mergeGaussianMixtureCov(cmps, projector, mean, angleDesc);
0248
0249
0250 return {mean, cov};
0251 }
0252
0253
0254
0255
0256
0257
0258
0259
0260
0261 template <typename component_range_t, typename projector_t,
0262 typename angle_desc_t>
0263 BoundVector mergeGaussianMixtureParams(const component_range_t &cmps,
0264 const projector_t &projector,
0265 const angle_desc_t &angleDesc,
0266 ComponentMergeMethod method) {
0267 if (method == ComponentMergeMethod::eMean) {
0268 return mergeGaussianMixtureMean(cmps, projector, angleDesc);
0269 } else if (method == ComponentMergeMethod::eMaxWeight) {
0270 return mergeGaussianMixtureMaxWeight(cmps, projector);
0271 } else {
0272 throw std::logic_error("Invalid component merge method");
0273 }
0274 }
0275
0276
0277
0278
0279
0280
0281
0282
0283
0284
0285 template <typename component_range_t, typename projector_t,
0286 typename angle_desc_t>
0287 std::tuple<BoundVector, BoundMatrix> mergeGaussianMixture(
0288 const component_range_t &cmps, const projector_t &projector,
0289 const angle_desc_t &angleDesc, ComponentMergeMethod method) {
0290 const auto [mean, cov] =
0291 mergeGaussianMixtureMeanCov(cmps, projector, angleDesc);
0292
0293 if (method == ComponentMergeMethod::eMean) {
0294 return {mean, cov};
0295 } else if (method == ComponentMergeMethod::eMaxWeight) {
0296 const BoundVector mode =
0297 mergeGaussianMixtureMaxWeight(cmps, [&](const auto &c) {
0298 const auto [weight_l, pars_l, cov_l] = projector(c);
0299 return std::tuple(weight_l, pars_l);
0300 });
0301 return {mode, cov};
0302 } else {
0303 throw std::logic_error("Invalid component merge method");
0304 }
0305 }
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315 template <typename component_range_t, typename projector_t>
0316 BoundVector mergeGaussianMixtureParams(const component_range_t &cmps,
0317 const projector_t &projector,
0318 const Surface &surface,
0319 ComponentMergeMethod method) {
0320 return angleDescriptionSwitch(surface, [&](const auto &desc) {
0321 return mergeGaussianMixtureParams(cmps, projector, desc, method);
0322 });
0323 }
0324
0325
0326
0327
0328
0329
0330
0331
0332
0333
0334 template <typename component_range_t, typename projector_t>
0335 std::tuple<BoundVector, BoundMatrix> mergeGaussianMixture(
0336 const component_range_t &cmps, const projector_t &projector,
0337 const Surface &surface, ComponentMergeMethod method) {
0338 return angleDescriptionSwitch(surface, [&](const auto &desc) {
0339 return mergeGaussianMixture(cmps, projector, desc, method);
0340 });
0341 }
0342
0343
0344 class SymmetricKLDistanceMatrix {
0345 using Array = Eigen::Array<double, Eigen::Dynamic, 1>;
0346 using Mask = Eigen::Array<bool, Eigen::Dynamic, 1>;
0347
0348 Array m_distances;
0349 Mask m_mask;
0350 std::vector<std::pair<std::size_t, std::size_t>> m_mapToPair;
0351 std::size_t m_numberComponents{};
0352
0353 template <typename array_t, typename setter_t>
0354 void setAssociated(std::size_t n, array_t &array, const setter_t &setter) {
0355 const std::size_t indexConst = (n - 1) * n / 2;
0356
0357 for (std::size_t i = 0ul; i < n; ++i) {
0358 array[indexConst + i] = setter(n, i);
0359 }
0360
0361 for (std::size_t i = n + 1; i < m_numberComponents; ++i) {
0362 array[(i - 1) * i / 2 + n] = setter(n, i);
0363 }
0364 }
0365
0366
0367
0368 static double computeSymmetricKlDivergence(const GsfComponent &a,
0369 const GsfComponent &b);
0370
0371 public:
0372 explicit SymmetricKLDistanceMatrix(std::span<const GsfComponent> cmps);
0373
0374 double at(std::size_t i, std::size_t j) const;
0375
0376 void recomputeAssociatedDistances(std::size_t n,
0377 std::span<const GsfComponent> cmps);
0378
0379 void maskAssociatedDistances(std::size_t n);
0380
0381 std::pair<std::size_t, std::size_t> minDistancePair() const;
0382
0383 std::ostream &toStream(std::ostream &os) const;
0384
0385 friend std::ostream &operator<<(std::ostream &os,
0386 const SymmetricKLDistanceMatrix &m) {
0387 return m.toStream(os);
0388 }
0389 };
0390
0391 }