File indexing completed on 2026-04-09 07:50:08
0001
0002
0003
0004 #include <DD4hep/VolumeManager.h>
0005 #include <Evaluator/DD4hepUnits.h>
0006 #include <Math/GenVector/Cartesian3D.h>
0007 #include <Math/GenVector/DisplacementVector3D.h>
0008 #include <algorithms/geo.h>
0009 #include <edm4eic/Cov6f.h>
0010 #include <edm4eic/EDM4eicVersion.h>
0011 #include <edm4eic/MCRecoTrackParticleAssociationCollection.h>
0012 #include <edm4eic/MCRecoTrackerHitAssociationCollection.h>
0013 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0014 #include <edm4eic/MCRecoTrackerHitLinkCollection.h>
0015 #endif
0016 #include <edm4eic/Measurement2DCollection.h>
0017 #include <edm4eic/RawTrackerHit.h>
0018 #include <edm4eic/TrackCollection.h>
0019 #include <edm4eic/TrackerHit.h>
0020 #include <edm4hep/MCParticle.h>
0021 #include <edm4hep/SimTrackerHit.h>
0022 #include <edm4hep/Vector2f.h>
0023 #include <edm4hep/Vector3d.h>
0024 #include <edm4hep/Vector3f.h>
0025 #include <edm4hep/utils/vector_utils.h>
0026 #include <podio/LinkNavigator.h>
0027 #include <podio/RelationRange.h>
0028 #include <podio/detail/Link.h>
0029 #include <Eigen/Geometry>
0030 #include <Eigen/Householder>
0031 #include <Eigen/Jacobi>
0032 #include <Eigen/SVD>
0033 #include <algorithm>
0034 #include <cmath>
0035 #include <cstddef>
0036 #include <cstdint>
0037 #include <memory>
0038 #include <new>
0039 #include <unordered_map>
0040 #include <utility>
0041
0042 #include "FarDetectorLinearTracking.h"
0043 #include "algorithms/fardetectors/FarDetectorLinearTrackingConfig.h"
0044
0045 namespace eicrecon {
0046
0047 void FarDetectorLinearTracking::init() {
0048
0049
0050 m_layerWeights = Eigen::VectorXd::Constant(m_cfg.n_layer, 1);
0051
0052 for (std::size_t i = 0; i < std::min(m_cfg.layer_weights.size(), m_cfg.n_layer); i++) {
0053 m_layerWeights(i) = m_cfg.layer_weights[i];
0054 }
0055
0056
0057 m_optimumDirection = Eigen::Vector3d::UnitZ();
0058 m_optimumDirection =
0059 Eigen::AngleAxisd(m_cfg.optimum_theta, Eigen::Vector3d::UnitY()) * m_optimumDirection;
0060 m_optimumDirection =
0061 Eigen::AngleAxisd(m_cfg.optimum_phi, Eigen::Vector3d::UnitZ()) * m_optimumDirection;
0062
0063 m_cellid_converter = algorithms::GeoSvc::instance().cellIDPositionConverter();
0064 }
0065
0066 void FarDetectorLinearTracking::process(const FarDetectorLinearTracking::Input& input,
0067 const FarDetectorLinearTracking::Output& output) const {
0068
0069 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0070 const auto [inputhits, hitLinks, assocHits] = input;
0071 auto [outputTracks, trackLinks, assocTracks] = output;
0072 #else
0073 const auto [inputhits, assocHits] = input;
0074 auto [outputTracks, assocTracks] = output;
0075 #endif
0076
0077
0078 std::size_t nCollections = inputhits.size();
0079 if (nCollections != m_cfg.n_layer) {
0080 error("Wrong number of input collections passed to algorithm");
0081 return;
0082 }
0083
0084 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0085
0086 const bool do_assoc = hitLinks != nullptr && !hitLinks->empty();
0087 if (!do_assoc) {
0088 debug("Provided MCRecoTrackerHitLink collection is empty. No truth associations "
0089 "will be performed.");
0090 }
0091
0092 std::optional<podio::LinkNavigator<edm4eic::MCRecoTrackerHitLinkCollection>> link_nav;
0093 if (do_assoc) {
0094 link_nav.emplace(*hitLinks);
0095 }
0096 #else
0097 const bool do_assoc = assocHits != nullptr && !assocHits->empty();
0098 if (!do_assoc) {
0099 debug("Provided MCRecoTrackerHitAssociation collection is empty. No truth associations "
0100 "will be performed.");
0101 }
0102 #endif
0103
0104 std::vector<std::vector<Eigen::Vector3d>> convertedHits;
0105 std::vector<std::vector<edm4hep::MCParticle>> assocParts;
0106
0107
0108
0109
0110 for (const auto& layerHits : inputhits) {
0111 if ((*layerHits).size() > m_cfg.layer_hits_max) {
0112 info("Too many hits in layer");
0113 return;
0114 }
0115 if ((*layerHits).empty()) {
0116 trace("No hits in layer");
0117 return;
0118 }
0119 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0120 ConvertClusters(*layerHits, *link_nav, *assocHits, convertedHits, assocParts);
0121 #else
0122 ConvertClusters(*layerHits, *assocHits, convertedHits, assocParts);
0123 #endif
0124 }
0125
0126
0127 Eigen::MatrixXd hitMatrix(3, m_cfg.n_layer);
0128
0129
0130 std::vector<std::size_t> layerHitIndex(m_cfg.n_layer, 0);
0131
0132 int layer = 0;
0133
0134
0135 while (true) {
0136 hitMatrix.col(layer) << convertedHits[layer][layerHitIndex[layer]];
0137
0138 bool isValid = true;
0139
0140 if (layer > 0 && m_cfg.restrict_direction) {
0141 isValid = checkHitPair(hitMatrix.col(layer - 1), hitMatrix.col(layer));
0142 }
0143
0144
0145 if (isValid) {
0146 if (layer == static_cast<long>(m_cfg.n_layer) - 1) {
0147
0148 checkHitCombination(&hitMatrix, outputTracks,
0149 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0150 trackLinks,
0151 #endif
0152 assocTracks, inputhits, assocParts, layerHitIndex);
0153 } else {
0154 layer++;
0155 continue;
0156 }
0157 }
0158
0159
0160 layerHitIndex[layer]++;
0161
0162 bool doBreak = false;
0163
0164 while (layerHitIndex[layer] >= convertedHits[layer].size()) {
0165 layerHitIndex[layer] = 0;
0166 if (layer == 0) {
0167 doBreak = true;
0168 break;
0169 }
0170 layer--;
0171
0172 layerHitIndex[layer]++;
0173 }
0174 if (doBreak) {
0175 break;
0176 }
0177 }
0178 }
0179
0180 void FarDetectorLinearTracking::checkHitCombination(
0181 Eigen::MatrixXd* hitMatrix, edm4eic::TrackCollection* outputTracks,
0182 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0183 edm4eic::MCRecoTrackParticleLinkCollection* trackLinks,
0184 #endif
0185 edm4eic::MCRecoTrackParticleAssociationCollection* assocTracks,
0186 const std::vector<gsl::not_null<const edm4eic::Measurement2DCollection*>>& inputHits,
0187 const std::vector<std::vector<edm4hep::MCParticle>>& assocParts,
0188 const std::vector<std::size_t>& layerHitIndex) const {
0189
0190 Eigen::Vector3d weightedAnchor = (*hitMatrix) * m_layerWeights / (m_layerWeights.sum());
0191
0192 auto localMatrix = (*hitMatrix).colwise() - weightedAnchor;
0193
0194 Eigen::BDCSVD<Eigen::MatrixXd> svd(localMatrix.transpose(),
0195 Eigen::ComputeThinU | Eigen::ComputeThinV);
0196
0197 auto V = svd.matrixV();
0198
0199
0200 auto rotatedMatrix = localMatrix.transpose() * V;
0201 auto residuals = rotatedMatrix.rightCols(2);
0202 double chi2 = (residuals.array() * residuals.array()).sum() / (2 * m_cfg.n_layer);
0203
0204 if (chi2 > m_cfg.chi2_max) {
0205 return;
0206 }
0207
0208 edm4hep::Vector3d outPos = weightedAnchor.data();
0209 edm4hep::Vector3d outVec = V.col(0).data();
0210
0211
0212 if (outVec.z > 0) {
0213 outVec = outVec * -1;
0214 }
0215
0216 int32_t type{0};
0217 edm4hep::Vector3f position(outPos.x, outPos.y, outPos.z);
0218 edm4hep::Vector3f momentum(outVec.x, outVec.y, outVec.z);
0219 edm4eic::Cov6f positionMomentumCovariance;
0220 float time{0};
0221 float timeError{0};
0222 float charge{-1};
0223 int32_t ndf{static_cast<int32_t>(m_cfg.n_layer) - 1};
0224 int32_t pdg{11};
0225
0226
0227 auto track = (*outputTracks)
0228 .create(type, position, momentum, positionMomentumCovariance, time, timeError,
0229 charge, chi2, ndf, pdg);
0230
0231
0232 std::unordered_map<const edm4hep::MCParticle*, int> particleCount;
0233 for (std::size_t layer = 0; layer < layerHitIndex.size(); layer++) {
0234 track.addToMeasurements((*inputHits[layer])[layerHitIndex[layer]]);
0235 const auto& assocParticle = assocParts[layer][layerHitIndex[layer]];
0236 particleCount[&assocParticle]++;
0237 }
0238
0239
0240 for (const auto& [particle, count] : particleCount) {
0241 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0242 auto trackLink = trackLinks->create();
0243 trackLink.setFrom(track);
0244 trackLink.setTo(*particle);
0245 trackLink.setWeight(count / static_cast<double>(m_cfg.n_layer));
0246 #endif
0247 auto trackAssoc = assocTracks->create();
0248 trackAssoc.setRec(track);
0249 trackAssoc.setSim(*particle);
0250 trackAssoc.setWeight(count / static_cast<double>(m_cfg.n_layer));
0251 }
0252 }
0253
0254
0255 bool FarDetectorLinearTracking::checkHitPair(const Eigen::Vector3d& hit1,
0256 const Eigen::Vector3d& hit2) const {
0257
0258 Eigen::Vector3d hitDiff = hit2 - hit1;
0259 hitDiff.normalize();
0260
0261 double angle = std::acos(hitDiff.dot(m_optimumDirection));
0262
0263 debug("Vector: x={}, y={}, z={}", hitDiff.x(), hitDiff.y(), hitDiff.z());
0264 debug("Optimum: x={}, y={}, z={}", m_optimumDirection.x(), m_optimumDirection.y(),
0265 m_optimumDirection.z());
0266 debug("Angle: {}, Tolerance {}", angle, m_cfg.step_angle_tolerance);
0267
0268 return angle <= m_cfg.step_angle_tolerance;
0269 }
0270
0271
0272 void FarDetectorLinearTracking::ConvertClusters(
0273 const edm4eic::Measurement2DCollection& clusters,
0274 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0275 const podio::LinkNavigator<edm4eic::MCRecoTrackerHitLinkCollection>& link_nav,
0276 #endif
0277 [[maybe_unused]] const edm4eic::MCRecoTrackerHitAssociationCollection& assoc_hits,
0278 std::vector<std::vector<Eigen::Vector3d>>& pointPositions,
0279 std::vector<std::vector<edm4hep::MCParticle>>& assoc_parts) const {
0280
0281
0282 const dd4hep::VolumeManagerContext* context =
0283 m_cellid_converter->findContext(clusters[0].getSurface());
0284
0285 std::vector<Eigen::Vector3d> layerPositions;
0286 std::vector<edm4hep::MCParticle> assocParticles;
0287
0288 for (auto cluster : clusters) {
0289
0290 auto globalPos = context->localToWorld({cluster.getLoc()[0], cluster.getLoc()[1], 0});
0291 layerPositions.emplace_back(globalPos.x() / dd4hep::mm, globalPos.y() / dd4hep::mm,
0292 globalPos.z() / dd4hep::mm);
0293
0294
0295
0296 float maxWeight = 0;
0297 std::size_t maxIndex = cluster.getWeights().size();
0298 for (std::size_t i = 0; i < cluster.getWeights().size(); ++i) {
0299 if (cluster.getWeights()[i] > maxWeight) {
0300 maxWeight = cluster.getWeights()[i];
0301 maxIndex = i;
0302 }
0303 }
0304 if (maxIndex == cluster.getWeights().size()) {
0305
0306 continue;
0307 }
0308 auto maxHit = cluster.getHits()[maxIndex];
0309
0310 auto rawHit = maxHit.getRawHit();
0311
0312 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0313 const auto sim_hits = link_nav.getLinked(rawHit);
0314 if (!sim_hits.empty()) {
0315 auto particle = sim_hits[0].o.getParticle();
0316 assocParticles.push_back(particle);
0317 }
0318 #else
0319
0320 for (const auto& assoc : assoc_hits) {
0321 if (assoc.getRawHit() == rawHit) {
0322 auto particle = assoc.getSimHit().getParticle();
0323 assocParticles.push_back(particle);
0324 break;
0325 }
0326 }
0327 #endif
0328 }
0329
0330 pointPositions.push_back(layerPositions);
0331 assoc_parts.push_back(assocParticles);
0332 }
0333
0334 }