Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-17 08:01:31

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