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