|
|
|||
File indexing completed on 2025-12-16 09:22:20
0001 // This file is part of the ACTS project. 0002 // 0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project 0004 // 0005 // This Source Code Form is subject to the terms of the Mozilla Public 0006 // License, v. 2.0. If a copy of the MPL was not distributed with this 0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/. 0008 0009 #pragma once 0010 0011 #include "Acts/Definitions/TrackParametrization.hpp" 0012 #include "Acts/EventData/MultiTrajectory.hpp" 0013 #include "Acts/EventData/MultiTrajectoryBackendConcept.hpp" 0014 #include "Acts/EventData/ParticleHypothesis.hpp" 0015 #include "Acts/EventData/TrackContainerBackendConcept.hpp" 0016 #include "Acts/EventData/TrackParameters.hpp" 0017 #include "Acts/EventData/TrackProxyConcept.hpp" 0018 #include "Acts/EventData/TrackStatePropMask.hpp" 0019 #include "Acts/Utilities/HashedString.hpp" 0020 #include "Acts/Utilities/TypeTraits.hpp" 0021 #include "Acts/Utilities/UnitVectors.hpp" 0022 0023 #include <iterator> 0024 0025 namespace Acts { 0026 0027 namespace detail_tp { 0028 inline constexpr HashedString kTipIndexKey = hashString("tipIndex"); 0029 inline constexpr HashedString kStemIndexKey = hashString("stemIndex"); 0030 inline constexpr HashedString kMeasurementsKey = hashString("nMeasurements"); 0031 inline constexpr HashedString kHolesKey = hashString("nHoles"); 0032 inline constexpr HashedString kOutliersKey = hashString("nOutliers"); 0033 inline constexpr HashedString kSharedHitsKey = hashString("nSharedHits"); 0034 inline constexpr HashedString kChi2Key = hashString("chi2"); 0035 inline constexpr HashedString kNdfKey = hashString("ndf"); 0036 inline constexpr HashedString kNextKey = hashString("next"); 0037 } // namespace detail_tp 0038 0039 template <TrackContainerBackend track_container_t, 0040 CommonMultiTrajectoryBackend traj_t, 0041 template <typename> class holder_t> 0042 class TrackContainer; 0043 0044 /// Proxy class representing a single track. 0045 /// This class provides a **view** into an associated @ref TrackContainer, and 0046 /// has **reference semantics**. You can think of it as a pointer to a vector 0047 /// of tracks, which exposes an object-oriented interface for accessing the 0048 /// track properties. 0049 /// 0050 /// @tparam track_container_t the container backend 0051 /// @tparam trajectory_t the track state container backend 0052 /// @tparam holder_t ownership management class for the backend 0053 /// @tparam read_only true if this track container is not mutable 0054 template <typename track_container_t, typename trajectory_t, 0055 template <typename> class holder_t, bool read_only = true> 0056 class TrackProxy { 0057 public: 0058 /// Indicates whether this track proxy is read-only or if it can be modified 0059 static constexpr bool ReadOnly = read_only; 0060 0061 /// The track container backend given as a template parameter 0062 using Container = track_container_t; 0063 0064 /// The track state container backend given as a template parameter 0065 using Trajectory = trajectory_t; 0066 0067 /// The index type of the track container 0068 using IndexType = typename Container::IndexType; 0069 0070 /// Sentinel value that indicates an invalid index 0071 static constexpr IndexType kInvalid = Container::kInvalid; 0072 0073 /// Alias for the mutable version of this track proxy, with the same backends 0074 using MutableTrackProxy = 0075 TrackProxy<track_container_t, trajectory_t, holder_t, false>; 0076 0077 /// Alias for the const version of this track proxy, with the same backends 0078 using ConstTrackProxy = 0079 TrackProxy<track_container_t, trajectory_t, holder_t, true>; 0080 0081 /// Alias for an associated const track proxy, with the same backends 0082 using ConstProxyType = ConstTrackProxy; 0083 0084 /// Alias for an associated mutable track state proxy, with the same backends 0085 using TrackStateProxy = typename Trajectory::TrackStateProxy; 0086 0087 /// Alias for an associated const track state proxy, with the same backends 0088 using ConstTrackStateProxy = typename Trajectory::ConstTrackStateProxy; 0089 0090 /// Map-type for a bound parameter vector. This has reference semantics, i.e. 0091 /// points at a matrix by an internal pointer. 0092 using Parameters = 0093 typename detail_lt::FixedSizeTypes<eBoundSize, false>::CoefficientsMap; 0094 0095 /// Same as @ref Parameters, but with const semantics 0096 using ConstParameters = 0097 typename detail_lt::FixedSizeTypes<eBoundSize, true>::CoefficientsMap; 0098 0099 /// Map-type for a bound covariance. This has reference semantics, i.e. 0100 /// points at a matrix by an internal pointer. 0101 using Covariance = 0102 typename detail_lt::FixedSizeTypes<eBoundSize, false>::CovarianceMap; 0103 0104 /// Same as @ref Covariance, but with const semantics 0105 using ConstCovariance = 0106 typename detail_lt::FixedSizeTypes<eBoundSize, true>::CovarianceMap; 0107 0108 #ifndef DOXYGEN 0109 friend TrackContainer<Container, Trajectory, holder_t>; 0110 friend MutableTrackProxy; 0111 friend ConstTrackProxy; 0112 // Track proxies are friends, not food! 0113 template <typename A, typename B, template <typename> class H, bool R> 0114 friend class TrackProxy; 0115 #endif 0116 0117 /// @anchor track_proxy_construct 0118 /// @name Constructors and assignment operator 0119 /// 0120 /// Public constructors and assignment operators for @c TrackProxy only 0121 /// allow construction from another @c TrackProxy. You should generally 0122 /// not have to construct @c TrackProxy manually. 0123 /// 0124 /// @{ 0125 0126 /// Copy constructor: const to const or mutable to mutable 0127 /// @param other the other track proxy 0128 TrackProxy(const TrackProxy& other) = default; 0129 0130 /// Copy assignment operator: const to const or mutable to mutable 0131 /// @param other the other track proxy 0132 /// @return reference to this track proxy 0133 TrackProxy& operator=(const TrackProxy& other) = default; 0134 0135 /// Constructor from mutable track proxy 0136 /// @note Only available if the track proxy is read-only 0137 /// @param other the other track proxy 0138 explicit TrackProxy(const MutableTrackProxy& other) 0139 requires ReadOnly 0140 : m_container{other.m_container}, m_index{other.m_index} {} 0141 0142 /// Copy assignment operator from mutable track proxy 0143 /// @note Only available if the track proxy is read-only 0144 /// @param other the other track proxy 0145 /// @return reference to this track proxy 0146 TrackProxy& operator=(const MutableTrackProxy& other) 0147 requires ReadOnly 0148 { 0149 m_container = other.m_container; 0150 m_index = other.m_index; 0151 return *this; 0152 } 0153 0154 /// @} 0155 0156 /// Equality operator with another track proxy 0157 /// Checks the container identity and the track index 0158 /// @param other Other track proxy to compare with 0159 /// @return True if the track proxies refer to the same track 0160 bool operator==(const TrackProxy& other) const { 0161 return &(*m_container) == &(*other.m_container) && m_index == other.m_index; 0162 } 0163 0164 /// @anchor track_proxy_props 0165 /// @name TrackProxy properties 0166 /// Methods that give access to the properties of a track represented by 0167 /// @c TrackProxy. 0168 /// 0169 /// Many of these methods come in a @c const and a non-@c const version. The 0170 /// non-@c const version is only available if you have an instance of 0171 /// @c TrackProxy that does not have the @c read_only template parameter set to 0172 /// @c true, even if you hold it as an lvalue. 0173 /// 0174 /// @{ 0175 0176 /// Get the tip index, i.e. the entry point into the track state container 0177 /// @return the tip index by value 0178 IndexType tipIndex() const { 0179 return component<IndexType>(detail_tp::kTipIndexKey); 0180 } 0181 0182 /// Index of the stem, i.e. the innermost track state of the track. 0183 /// This might be invalid, signifying that the track state is not 0184 /// forward-linked. 0185 /// @return the stem index 0186 IndexType stemIndex() const { 0187 return component<IndexType>(detail_tp::kStemIndexKey); 0188 } 0189 0190 /// Get a mutable reference to the tip index, i.e. the entry point into the 0191 /// track container 0192 /// @note Only available if the track proxy is not read-only 0193 /// @return mutable reference to the tip index 0194 IndexType& tipIndex() 0195 requires(!ReadOnly) 0196 { 0197 return component<IndexType>(detail_tp::kTipIndexKey); 0198 } 0199 0200 /// Index of the stem, i.e. the innermost track state of the track. 0201 /// This might be invalid, signifying that the track state is not 0202 /// forward-linked. 0203 /// @note Only available if the track proxy is not read-only 0204 /// @return mutable reference to the stem index 0205 IndexType& stemIndex() 0206 requires(!ReadOnly) 0207 { 0208 return component<IndexType>(detail_tp::kStemIndexKey); 0209 } 0210 0211 /// Get the reference surface of the track (e.g. the perigee) 0212 /// @return the reference surface 0213 const Surface& referenceSurface() const { 0214 return *m_container->container().referenceSurface_impl(m_index); 0215 } 0216 0217 // NOLINTBEGIN(performance-unnecessary-value-param) 0218 // looks like a false-positive. clang-tidy believes `srf` is not movable. 0219 /// Set a new reference surface for this track 0220 /// @param srf The surface to set 0221 void setReferenceSurface(std::shared_ptr<const Surface> srf) 0222 requires(!ReadOnly) 0223 { 0224 m_container->container().setReferenceSurface_impl(m_index, std::move(srf)); 0225 } 0226 // NOLINTEND(performance-unnecessary-value-param) 0227 0228 /// Returns whether the track has a reference surface or not 0229 /// @return whether a surface exists or not 0230 bool hasReferenceSurface() const { 0231 // @TODO: This could be more efficient 0232 return m_container->container().referenceSurface_impl(m_index) != nullptr; 0233 } 0234 0235 /// Get the parameters of the track at the reference surface (e.g. perigee). 0236 /// Const version 0237 /// @return Proxy vector for the parameters 0238 ConstParameters parameters() const { 0239 return m_container->parameters(m_index); 0240 } 0241 0242 /// Get the covariance of the track at the reference surface (e.g. perigee). 0243 /// Const version 0244 /// @return Proxy matrix for the covariance 0245 ConstCovariance covariance() const { 0246 return m_container->covariance(m_index); 0247 } 0248 0249 /// Get the parameters of the track at the reference surface (e.g. perigee). 0250 /// Mutable version 0251 /// @note Only available if the track proxy is not read-only 0252 /// @return Proxy vector for the parameters 0253 Parameters parameters() 0254 requires(!ReadOnly) 0255 { 0256 return m_container->parameters(m_index); 0257 } 0258 0259 /// Get the covariance of the track at the reference surface (e.g. perigee). 0260 /// Mutable version 0261 /// @note Only available if the track proxy is not read-only 0262 /// @return Proxy matrix for the covariance 0263 Covariance covariance() 0264 requires(!ReadOnly) 0265 { 0266 return m_container->covariance(m_index); 0267 } 0268 0269 /// Access the theta parameter of the track at the reference surface 0270 /// @return The theta parameter 0271 double theta() const { return parameters()[eBoundTheta]; } 0272 0273 /// Access the phi parameter of the track at the reference surface 0274 /// @return The phi parameter 0275 double phi() const { return parameters()[eBoundPhi]; } 0276 0277 /// Access the loc0 parameter of the track at the reference surface 0278 /// @return The loc0 parameter 0279 double loc0() const { return parameters()[eBoundLoc0]; } 0280 0281 /// Access the loc1 parameter of the track at the reference surface 0282 /// @return The loc1 parameter 0283 double loc1() const { return parameters()[eBoundLoc1]; } 0284 0285 /// Access the time parameter of the track at the reference surface 0286 /// @return The time parameter 0287 double time() const { return parameters()[eBoundTime]; } 0288 0289 /// Access the q/p (curvature) parameter of the track at the reference surface 0290 /// @return The q/p parameter 0291 double qOverP() const { return parameters()[eBoundQOverP]; } 0292 0293 /// Get the particle hypothesis 0294 /// @return the particle hypothesis 0295 ParticleHypothesis particleHypothesis() const { 0296 return m_container->container().particleHypothesis_impl(m_index); 0297 } 0298 0299 /// Set a new particle hypothesis for this track 0300 /// @note Only available if the track proxy is not read-only 0301 /// @param particleHypothesis The particle hypothesis to set 0302 void setParticleHypothesis(const ParticleHypothesis& particleHypothesis) 0303 requires(!ReadOnly) 0304 { 0305 m_container->container().setParticleHypothesis_impl(m_index, 0306 particleHypothesis); 0307 } 0308 0309 /// Get the charge of the tack 0310 /// @note this depends on the charge hypothesis 0311 /// @return The absolute track momentum 0312 double charge() const { return particleHypothesis().extractCharge(qOverP()); } 0313 0314 /// Get the absolute momentum of the tack 0315 /// @return The absolute track momentum 0316 double absoluteMomentum() const { 0317 return particleHypothesis().extractMomentum(qOverP()); 0318 } 0319 0320 /// Get the transverse momentum of the track 0321 /// @return The track transverse momentum value 0322 double transverseMomentum() const { 0323 return std::sin(theta()) * absoluteMomentum(); 0324 } 0325 0326 /// Get a unit vector along the track direction at the reference surface 0327 /// @return The direction unit vector 0328 Vector3 direction() const { 0329 return makeDirectionFromPhiTheta(phi(), theta()); 0330 } 0331 0332 /// Get the global momentum vector 0333 /// @return the global momentum vector 0334 Vector3 momentum() const { return absoluteMomentum() * direction(); } 0335 0336 /// Get the four-momentum vector: (px, py, pz, e) 0337 /// @return the four-momentum vector 0338 Vector4 fourMomentum() const { 0339 Vector4 p4 = Vector4::Zero(); 0340 0341 Vector3 p3 = momentum(); 0342 p4[eMom0] = p3[eMom0]; 0343 p4[eMom1] = p3[eMom1]; 0344 p4[eMom2] = p3[eMom2]; 0345 0346 float m = particleHypothesis().mass(); 0347 p4[eEnergy] = std::sqrt(m * m + p3.squaredNorm()); 0348 0349 return p4; 0350 } 0351 0352 /// Return the number of track states associated to this track 0353 /// @note This is calculated by iterating over the track states which is 0354 /// somewhat expensive. Consider caching this value if you need It 0355 /// more than once. 0356 /// @return The number of track states 0357 unsigned int nTrackStates() const { 0358 // @TODO: This should probably be cached, distance is expensive 0359 // without random access 0360 if (tipIndex() == kInvalid) { 0361 // no tip index -> no track states 0362 return 0; 0363 } 0364 auto tsRange = trackStatesReversed(); 0365 return std::distance(tsRange.begin(), tsRange.end()); 0366 } 0367 0368 /// Return a mutable reference to the number of measurements for the track. 0369 /// Mutable version 0370 /// @note Only available if the track proxy is not read-only 0371 /// @return The number of measurements 0372 unsigned int& nMeasurements() 0373 requires(!ReadOnly) 0374 { 0375 return component<unsigned int, detail_tp::kMeasurementsKey>(); 0376 } 0377 0378 /// Return the number of measurements for the track. Const version 0379 /// @return The number of measurements 0380 unsigned int nMeasurements() const { 0381 return component<unsigned int, detail_tp::kMeasurementsKey>(); 0382 } 0383 0384 /// Return a mutable reference to the number of holes for the track. 0385 /// Mutable version 0386 /// @note Only available if the track proxy is not read-only 0387 /// @return The number of holes 0388 unsigned int& nHoles() 0389 requires(!ReadOnly) 0390 { 0391 return component<unsigned int, detail_tp::kHolesKey>(); 0392 } 0393 0394 /// Return the number of measurements for the track. Const version 0395 /// @return The number of measurements 0396 unsigned int nHoles() const { 0397 return component<unsigned int, detail_tp::kHolesKey>(); 0398 } 0399 0400 /// Return a mutable reference to the number of outliers for the track. 0401 /// Mutable version 0402 /// @note Only available if the track proxy is not read-only 0403 /// @return The number of outliers 0404 unsigned int& nOutliers() 0405 requires(!ReadOnly) 0406 { 0407 return component<unsigned int, detail_tp::kOutliersKey>(); 0408 } 0409 0410 /// Return the number of outliers for the track. Const version 0411 /// @return The number of outliers 0412 unsigned int nOutliers() const { 0413 return component<unsigned int, detail_tp::kOutliersKey>(); 0414 } 0415 0416 /// Return a mutable reference to the number of shared hits for the track. 0417 /// Mutable version 0418 /// @note Only available if the track proxy is not read-only 0419 /// @return The number of shared hits 0420 unsigned int& nSharedHits() 0421 requires(!ReadOnly) 0422 { 0423 return component<unsigned int, detail_tp::kSharedHitsKey>(); 0424 } 0425 0426 /// Return the number of shared hits for the track. Const version 0427 /// @return The number of shared hits 0428 unsigned int nSharedHits() const { 0429 return component<unsigned int, detail_tp::kSharedHitsKey>(); 0430 } 0431 0432 /// Return a mutable reference to the chi squared 0433 /// Mutable version 0434 /// @note Only available if the track proxy is not read-only 0435 /// @return The chi squared 0436 float& chi2() 0437 requires(!ReadOnly) 0438 { 0439 return component<float, detail_tp::kChi2Key>(); 0440 } 0441 0442 /// Return the chi squared for the track. Const version 0443 /// @return The chi squared 0444 float chi2() const { return component<float, detail_tp::kChi2Key>(); } 0445 0446 /// Return a mutable reference to the number of degrees of freedom for the 0447 /// track. Mutable version 0448 /// @note Only available if the track proxy is not read-only 0449 /// @return The number of degrees of freedom 0450 unsigned int& nDoF() 0451 requires(!ReadOnly) 0452 { 0453 return component<unsigned int, detail_tp::kNdfKey>(); 0454 } 0455 0456 /// Return the number of degrees of freedom for the track. Const version 0457 /// @return The number of degrees of freedom 0458 unsigned int nDoF() const { 0459 return component<unsigned int, detail_tp::kNdfKey>(); 0460 } 0461 0462 /// Return the index of this track in the track container 0463 /// @note This is separate from the tip index 0464 /// @return the track index 0465 IndexType index() const { return m_index; } 0466 0467 /// @} 0468 0469 /// @anchor track_proxy_track_states 0470 /// @name TrackProxy track state access 0471 /// Methods that give access to the track states of a track represented by @c TrackProxy. 0472 /// @{ 0473 0474 /// Return a const track state proxy to the outermost track state 0475 /// @return The outermost track state proxy 0476 ConstTrackStateProxy outermostTrackState() const { 0477 return m_container->trackStateContainer().getTrackState(tipIndex()); 0478 } 0479 0480 /// Return a mutable track state proxy to the outermost track state 0481 /// @return The outermost track state proxy 0482 TrackStateProxy outermostTrackState() 0483 requires(!ReadOnly) 0484 { 0485 return m_container->trackStateContainer().getTrackState(tipIndex()); 0486 } 0487 0488 /// Return a const track state proxy to the innermost track state 0489 /// @note This is only available, if the track is forward linked 0490 /// @return The innermost track state proxy 0491 auto innermostTrackState() const { 0492 using proxy_t = decltype(m_container->trackStateContainer().getTrackState( 0493 std::declval<IndexType>())); 0494 0495 IndexType stem = component<IndexType, detail_tp::kStemIndexKey>(); 0496 if (stem == kInvalid) { 0497 return std::optional<proxy_t>{}; 0498 } else { 0499 return std::optional<proxy_t>{ 0500 m_container->trackStateContainer().getTrackState(stem)}; 0501 } 0502 } 0503 0504 /// Return a mutable track state proxy to the innermost track state 0505 /// @note This is only available, if the track is forward linked 0506 /// @note Only available if the track proxy is not read-only 0507 /// @return The innermost track state proxy 0508 auto innermostTrackState() 0509 requires(!ReadOnly) 0510 { 0511 using proxy_t = decltype(m_container->trackStateContainer().getTrackState( 0512 std::declval<IndexType>())); 0513 0514 IndexType stem = component<IndexType>(detail_tp::kStemIndexKey); 0515 if (stem == kInvalid) { 0516 return std::optional<proxy_t>{}; 0517 } else { 0518 return std::optional<proxy_t>{ 0519 m_container->trackStateContainer().getTrackState(stem)}; 0520 } 0521 } 0522 0523 /// Get a range over the track states of this track. Return value is 0524 /// compatible with range based for loop. Const version 0525 /// @note This range is from the outside inwards! 0526 /// @return Track state range to iterate over 0527 auto trackStatesReversed() const { 0528 return m_container->reverseTrackStateRange(m_index); 0529 } 0530 0531 /// Get a range over the track states of this track. Return value is 0532 /// compatible with range based for loop. Mutable version 0533 /// @note Only available if the track proxy is not read-only 0534 /// @note This range is from the outside inwards! 0535 /// @return Track state range to iterate over 0536 auto trackStatesReversed() 0537 requires(!ReadOnly) 0538 { 0539 return m_container->reverseTrackStateRange(m_index); 0540 } 0541 0542 /// Get a range over the track states of this track. Return value is 0543 /// compatible with range based for loop. This overload returns a const-only 0544 /// track state range, which means you cannot modify the track states obtained 0545 /// in the iteration. 0546 /// @note This range is from the inside out! 0547 /// @warning This access direction is only possible if the track states are 0548 /// **forward-linked**. 0549 /// @return Track state range to iterate over 0550 auto trackStates() const { 0551 return m_container->forwardTrackStateRange(m_index); 0552 } 0553 0554 /// Get a range over the track states of this track. 0555 /// Return value is compatible with range based for loop. 0556 /// This overload returns a mutable track state range, which means you 0557 /// can modify the track states obtained in the iteration. 0558 /// @note Only available if the track proxy is not read-only 0559 /// @note This range is from the inside out! 0560 /// @warning This access direction is only possible if the track states are 0561 /// **forward-linked**. 0562 /// @return Track state range to iterate over 0563 auto trackStates() 0564 requires(!ReadOnly) 0565 { 0566 return m_container->forwardTrackStateRange(m_index); 0567 } 0568 0569 /// @} 0570 0571 /// @anchor track_proxy_track_state_manipulation 0572 /// @name TrackProxy track state manipulation 0573 /// Methods that manipulate the track states of a track represented by @c 0574 /// TrackProxy. 0575 /// 0576 /// **Copy Methods Overview:** 0577 /// 0578 /// Three main copy methods are available with different behaviors: 0579 /// - @c copyFrom(): Deep copy including all track states (creates new track 0580 /// states) 0581 /// - @c copyFromWithoutStates(): Copy only track properties, invalidate 0582 /// track state indices 0583 /// - @c copyFromShallow(): Shallow copy sharing the same track states (copy 0584 /// indices only) 0585 /// 0586 /// Choose based on your needs: 0587 /// - Use @c copyFrom() for independent track copies with separate track 0588 /// states 0589 /// - Use @c copyFromWithoutStates() to update track metadata without 0590 /// affecting trajectories 0591 /// - Use @c copyFromShallow() for lightweight copies when track states can 0592 /// be shared 0593 /// @{ 0594 0595 /// Forward connect a track. 0596 /// This means setting indices from the inside out on all track states. 0597 /// @note Only available if the track proxy is not read-only 0598 void linkForward() 0599 requires(!ReadOnly) 0600 { 0601 IndexType last = kInvalid; 0602 for (auto ts : trackStatesReversed()) { 0603 ts.template component<IndexType>(detail_tp::kNextKey) = last; 0604 last = ts.index(); 0605 } 0606 stemIndex() = last; 0607 } 0608 0609 /// Append a track state to this track. 0610 /// This will modify the tip index to point at the newly created track state, 0611 /// which will be directly after the previous track state at tip index. 0612 /// @note Only available if the track proxy is not read-only 0613 /// @param mask The allocation prop mask for the new track state 0614 /// @return The newly added track state 0615 auto appendTrackState(TrackStatePropMask mask = TrackStatePropMask::All) 0616 requires(!ReadOnly) 0617 { 0618 auto& tsc = m_container->trackStateContainer(); 0619 auto ts = tsc.makeTrackState(mask, tipIndex()); 0620 tipIndex() = ts.index(); 0621 return ts; 0622 } 0623 0624 /// Create a complete deep copy of another track, including all track states. 0625 /// This creates new track states in the destination trajectory and copies 0626 /// all data from the source track states. The track state sequence order 0627 /// is preserved. 0628 /// 0629 /// **Implementation details:** 0630 /// - Track states are initially copied in reversed order for efficiency 0631 /// - The track state links are then updated using reverseTrackStates() 0632 /// - As a consequence, the resulting track is forward-linked 0633 /// 0634 /// **What gets copied:** 0635 /// - All track-level properties (parameters, covariance, particle hypothesis, 0636 /// etc.) 0637 /// - Reference surface (shared pointer is copied) 0638 /// - Track summary data (nMeasurements, nHoles, chi2, etc.) 0639 /// - All dynamic track columns 0640 /// - Complete sequence of track states with all their data 0641 /// - All dynamic track state columns 0642 /// 0643 /// **Result:** 0644 /// - The destination track will have newly created track states 0645 /// - tipIndex() and stemIndex() will point to the new track states 0646 /// - Track state indices will be different from the source 0647 /// - All track state data will be identical to the source 0648 /// - The track will be forward-linked (stemIndex() will be valid) 0649 /// 0650 /// @note Only available if the track proxy is not read-only 0651 /// @note Both track containers must have compatible dynamic columns 0652 /// @tparam track_proxy_t the other track proxy's type 0653 /// @param other The source track proxy to copy from 0654 template <TrackProxyConcept track_proxy_t> 0655 void copyFrom(const track_proxy_t& other) 0656 requires(!ReadOnly) 0657 { 0658 copyFromWithoutStates(other); 0659 0660 // append track states (cheap), but they're in the wrong order 0661 for (const auto& srcTrackState : other.trackStatesReversed()) { 0662 auto destTrackState = appendTrackState(srcTrackState.getMask()); 0663 destTrackState.copyFrom(srcTrackState, Acts::TrackStatePropMask::All, 0664 true); 0665 } 0666 0667 // reverse using standard linked list reversal algorithm 0668 reverseTrackStates(); 0669 } 0670 0671 /// Copy track-level properties from another track, but not the track states. 0672 /// This copies all track metadata and properties but leaves the track state 0673 /// sequence unchanged. Useful when you want to copy track properties to an 0674 /// existing track that may already have track states. 0675 /// 0676 /// **What gets copied:** 0677 /// - Track parameters at reference surface 0678 /// - Covariance matrix at reference surface 0679 /// - Particle hypothesis 0680 /// - Reference surface (shared pointer is copied) 0681 /// - Track summary data (nMeasurements, nHoles, nOutliers, nSharedHits, chi2, 0682 /// nDoF) 0683 /// - All dynamic track columns 0684 /// 0685 /// **What does NOT get copied:** 0686 /// - Track states (existing track states remain unchanged in the container) 0687 /// 0688 /// **Result:** 0689 /// - All track-level properties are updated to match the source 0690 /// - tipIndex() and stemIndex() are set to kInvalid (track states become 0691 /// inaccessible) 0692 /// - Existing track states remain in the container but are no longer linked 0693 /// to this track 0694 /// - nTrackStates() will return 0 due to invalid indices 0695 /// 0696 /// @note Only available if the track proxy is not read-only 0697 /// @note Both track containers must have compatible dynamic columns 0698 /// @tparam track_proxy_t the other track proxy's type 0699 /// @param other The source track proxy to copy properties from 0700 template <TrackProxyConcept track_proxy_t> 0701 void copyFromWithoutStates(const track_proxy_t& other) 0702 requires(!ReadOnly) 0703 { 0704 setParticleHypothesis(other.particleHypothesis()); 0705 0706 if (other.hasReferenceSurface()) { 0707 setReferenceSurface(other.referenceSurface().getSharedPtr()); 0708 parameters() = other.parameters(); 0709 covariance() = other.covariance(); 0710 } else { 0711 setReferenceSurface(nullptr); 0712 } 0713 0714 nMeasurements() = other.nMeasurements(); 0715 nHoles() = other.nHoles(); 0716 nOutliers() = other.nOutliers(); 0717 nSharedHits() = other.nSharedHits(); 0718 chi2() = other.chi2(); 0719 nDoF() = other.nDoF(); 0720 0721 m_container->copyDynamicFrom(m_index, other.m_container->container(), 0722 other.m_index); 0723 0724 tipIndex() = kInvalid; 0725 stemIndex() = kInvalid; 0726 } 0727 0728 /// Create a shallow copy from another track, sharing the same track states. 0729 /// This copies all track-level properties and makes the destination track 0730 /// point to the same track state sequence as the source. The track states 0731 /// themselves are not duplicated - both tracks will reference the same 0732 /// track state objects in memory. 0733 /// 0734 /// **What gets copied:** 0735 /// - All track-level properties (parameters, covariance, particle hypothesis, 0736 /// etc.) 0737 /// - Reference surface (shared pointer is copied) 0738 /// - Track summary data (nMeasurements, nHoles, chi2, etc.) 0739 /// - All dynamic track columns 0740 /// - tipIndex() and stemIndex() (track state linking information) 0741 /// 0742 /// **What gets shared (not duplicated):** 0743 /// - Track states (both tracks reference the same track state objects) 0744 /// 0745 /// **Result:** 0746 /// - The destination track will have the same nTrackStates() as the source 0747 /// - Both tracks will iterate over the same track state sequence 0748 /// - Modifications to track states will be visible in both tracks 0749 /// - Track state indices will be identical between tracks 0750 /// - The destination track will have a different track index than the source 0751 /// 0752 /// @warning Modifying track states through either track will affect both tracks 0753 /// since they share the same track state objects 0754 /// @warning It is the user's responsibility to ensure that the tip and stem 0755 /// indices from the source track are valid in the destination 0756 /// track's track state container. No validation is performed - 0757 /// invalid indices will lead to undefined behavior when accessing 0758 /// track states 0759 /// @note Only available if the track proxy is not read-only 0760 /// @note Both track containers must have compatible dynamic columns 0761 /// @tparam track_proxy_t the other track proxy's type 0762 /// @param other The source track proxy to create a shallow copy from 0763 template <TrackProxyConcept track_proxy_t> 0764 void copyFromShallow(const track_proxy_t& other) 0765 requires(!ReadOnly) 0766 { 0767 copyFromWithoutStates(other); 0768 tipIndex() = other.tipIndex(); 0769 stemIndex() = other.stemIndex(); 0770 } 0771 0772 /// Reverse the ordering of track states for this track 0773 /// Afterwards, the previous endpoint of the track state sequence will be 0774 /// the "innermost" track state 0775 /// @note Only available if the track proxy is not read-only 0776 /// @note This is dangerous with branching track state sequences, as it will break them 0777 /// @note This also automatically forward-links the track! 0778 /// @param invertJacobians Whether to invert the Jacobians of the track states 0779 void reverseTrackStates(bool invertJacobians = false) 0780 requires(!ReadOnly) 0781 { 0782 IndexType current = tipIndex(); 0783 IndexType next = kInvalid; 0784 IndexType prev = kInvalid; 0785 0786 stemIndex() = tipIndex(); 0787 0788 // @TODO: Maybe refactor to not need this variable if invertJacobians == false 0789 BoundMatrix nextJacobian; 0790 0791 while (current != kInvalid) { 0792 auto ts = m_container->trackStateContainer().getTrackState(current); 0793 prev = ts.previous(); 0794 ts.template component<IndexType>(detail_tp::kNextKey) = prev; 0795 ts.previous() = next; 0796 if (invertJacobians) { 0797 if (next != kInvalid) { 0798 BoundMatrix curJacobian = ts.jacobian(); 0799 ts.jacobian() = nextJacobian.inverse(); 0800 nextJacobian = curJacobian; 0801 } else { 0802 nextJacobian = ts.jacobian(); 0803 ts.jacobian().setZero(); 0804 } 0805 } 0806 next = current; 0807 tipIndex() = current; 0808 current = prev; 0809 } 0810 } 0811 0812 /// @} 0813 0814 /// @anchor track_proxy_generic_component 0815 /// @name TrackProxy generic component access 0816 /// Methods that give access to generic components of a track represented by 0817 /// @c TrackProxy. Internally, a compile-time hash of the component name is 0818 /// used to identify which component is being requested. Most of the named 0819 /// methods in @ref track_proxy_props "TrackProxy properties" use these 0820 /// methods to retrieve the actual data. 0821 /// 0822 /// A number of overloads exist, where you can either supply the 0823 /// @ref HashedString @c key as a template parameter or a runtime argument. The 0824 /// former has the advantage of being guaranteed to be evaluated at 0825 /// compile-time. 0826 /// 0827 /// @{ 0828 0829 /// Retrieve a mutable reference to a component 0830 /// @tparam T The type of the component to access 0831 /// @tparam key String key for the component to access 0832 /// @return Mutable reference to the component given by @p key 0833 template <typename T, HashedString key> 0834 constexpr T& component() 0835 requires(!ReadOnly) 0836 { 0837 return m_container->template component<T, key>(m_index); 0838 } 0839 0840 /// Retrieve a mutable reference to a component 0841 /// @tparam T The type of the component to access 0842 /// @param key String key for the component to access 0843 /// @return Mutable reference to the component given by @p key 0844 template <typename T> 0845 constexpr T& component(HashedString key) 0846 requires(!ReadOnly) 0847 { 0848 return m_container->template component<T>(key, m_index); 0849 } 0850 0851 /// Retrieve a mutable reference to a component 0852 /// @tparam T The type of the component to access 0853 /// @param key String key for the component to access 0854 /// @note This might hash the @p key at runtime instead of compile-time 0855 /// @return Mutable reference to the component given by @p key 0856 template <typename T> 0857 constexpr T& component(std::string_view key) 0858 requires(!ReadOnly) 0859 { 0860 return m_container->template component<T>(hashStringDynamic(key), m_index); 0861 } 0862 0863 /// Retrieve a const reference to a component 0864 /// @tparam T The type of the component to access 0865 /// @tparam key String key for the component to access 0866 /// @return Const reference to the component given by @p key 0867 template <typename T, HashedString key> 0868 constexpr const T& component() const { 0869 return m_container->template component<T, key>(m_index); 0870 } 0871 0872 /// Check whether a dynamic column exists 0873 /// @param key String key for the component to check 0874 /// @return whether the column exists 0875 bool hasColumn(HashedString key) const { return m_container->hasColumn(key); } 0876 0877 /// Retrieve a const reference to a component 0878 /// @tparam T The type of the component to access 0879 /// @param key String key for the component to access 0880 /// @return Const reference to the component given by @p key 0881 template <typename T> 0882 constexpr const T& component(HashedString key) const { 0883 return m_container->template component<T>(key, m_index); 0884 } 0885 0886 /// Retrieve a const reference to a component 0887 /// @tparam T The type of the component to access 0888 /// @param key String key for the component to access 0889 /// @note This might hash the @p key at runtime instead of compile-time 0890 /// @return Const reference to the component given by @p key 0891 template <typename T> 0892 constexpr const T& component(std::string_view key) const { 0893 return m_container->template component<T>(hashStringDynamic(key), m_index); 0894 } 0895 0896 /// @} 0897 0898 /// Return the track parameters at the reference surface 0899 /// @note The parameters are created on the fly 0900 /// @return the track parameters 0901 BoundTrackParameters createParametersAtReference() const { 0902 return BoundTrackParameters(referenceSurface().getSharedPtr(), parameters(), 0903 covariance(), particleHypothesis()); 0904 } 0905 0906 /// Convert a track state into track parameters 0907 /// @note The parameters are created on the fly 0908 /// @param trackState Track state to convert to parameters 0909 /// @return the track parameters 0910 BoundTrackParameters createParametersFromState( 0911 const ConstTrackStateProxy& trackState) const { 0912 return BoundTrackParameters(trackState.referenceSurface().getSharedPtr(), 0913 trackState.parameters(), 0914 trackState.covariance(), particleHypothesis()); 0915 } 0916 0917 /// Return a reference to the track container backend, mutable version. 0918 /// @note Only available if the track proxy is not read-only 0919 /// @return reference to the track container backend 0920 auto& container() 0921 requires(!ReadOnly) 0922 { 0923 return *m_container; 0924 } 0925 0926 /// Return a reference to the track container backend, const version. 0927 /// @return reference to the track container backend 0928 const auto& container() const { return *m_container; } 0929 0930 private: 0931 TrackProxy( 0932 const_if_t<ReadOnly, TrackContainer<Container, Trajectory, holder_t>>& 0933 container, 0934 IndexType itrack) 0935 : m_container{&container}, m_index{itrack} {} 0936 0937 detail_lt::TransitiveConstPointer< 0938 const_if_t<ReadOnly, TrackContainer<Container, Trajectory, holder_t>>> 0939 m_container; 0940 IndexType m_index; 0941 }; 0942 0943 } // namespace Acts
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|