Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-15 09:42:16

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/Algebra.hpp"
0012 #include "Acts/EventData/CompositeSpacePoint.hpp"
0013 #include "Acts/EventData/CompositeSpacePointCalibrator.hpp"
0014 #include "Acts/Seeding/detail/CompSpacePointAuxiliaries.hpp"
0015 #include "Acts/Utilities/CalibrationContext.hpp"
0016 #include "Acts/Utilities/Logger.hpp"
0017 #include "Acts/Utilities/RangeXD.hpp"
0018 #include "Acts/Utilities/Result.hpp"
0019 
0020 #include <memory>
0021 #include <vector>
0022 
0023 namespace Acts::Experimental::detail {
0024 
0025 ///  @brief The FastStrawLineFitter fits a straight line to a set of straw measurements
0026 ///         The space points  passed to the fitter need to be all straw space
0027 ///         points and their straw wires need to be approximately parallel to
0028 ///         each other. Under this assumption, each straw measurements can be
0029 ///         projected onto the plane which is perpendicular to its straw wire
0030 ///         and the residual to the k-th straw is written as
0031 ///
0032 ///                 R_{k} = T_{z,k} * sin \theta - (T_{y,k} -  y_{0})* cos \theta - sign_{k} * r_{k}
0033 ///
0034 ///         where, \theta and y_{0} are the two line parameters to fit, T_{z,k} and T_{y,k} are
0035 ///         the coordinates of the straw-tube in the plane, sign_{k} fixes the
0036 ///         left-right ambiguity, and r_{k} is the drift radius of the
0037 ///         measurement. The z & y coordinates are given by the dot product of
0038 ///         the straw's local position with its planeNormal & toNextSensor
0039 ///         vector, respectively. Assuming that the drift radius does not need
0040 ///         to be calibrated during the fit, the entire problem reduces to the
0041 ///         minimization of a polynomial having the form
0042 ///
0043 ///                   A * cos(\theta)  + B * sin(\theta)  + C * cos(2\theta) + D * sin(2\theta)
0044 ///         where A, B, C, D are fit constants depending on the configuration of
0045 ///         the tubes to fit. They are documented in
0046 ///         https://gitlab.cern.ch/atlas-nextgen/work-package-2.5/analyticalsegment
0047 ///
0048 ///         In general, the primary source of the straw drift radius is a record
0049 ///         of the time-of-arrival of an electronics signal which is then
0050 ///         translated using a r-t relation which is assumed to be twice
0051 ///         differentiable. Just because a particle might arrive a bit later,
0052 ///         the drift radius then becomes a function of a generic time offset,
0053 ///         t_{0}:
0054 ///
0055 ///                       r(t) = r(t_{time of record} - t_{0})
0056 ///
0057 ///         The \chi^{2} then becomes a function of the two parameters (\theta, t_{0})
0058 ///         which need to be minimized by varying these two parameters.
0059 
0060 class FastStrawLineFitter {
0061  public:
0062   /// @brief abrivation of the residual indices to fetch the proper covariance
0063   using ResidualIdx = CompSpacePointAuxiliaries::ResidualIdx;
0064   /// @brief Vector type
0065   using Vector = CompSpacePointAuxiliaries::Vector;
0066   /// @brief Configuration object
0067   struct Config {
0068     /// @brief Number of maximum iterations
0069     std::size_t maxIter{10};
0070     /// @brief Cutoff to define the fit to be converged if the parameter change is below threshold
0071     double precCutOff{1.e-9};
0072   };
0073   /// @brief Constructor of the fast straw line fitter
0074   /// @param cfg: Reference to the fitter's configuration object
0075   /// @param logger: Optional overwrite of the logging object
0076   explicit FastStrawLineFitter(const Config& cfg,
0077                                std::unique_ptr<const Logger> logger =
0078                                    getDefaultLogger("FastStrawLineFitter",
0079                                                     Logging::Level::INFO));
0080 
0081   /// @brief Helper struct to pack the result of the straw line fit
0082   struct FitResult {
0083     virtual ~FitResult() = default;
0084     /// @brief Printer method
0085     virtual void print(std::ostream& ostr) const;
0086     /// @brief Ostream operator
0087     friend std::ostream& operator<<(std::ostream& ostr, const FitResult& x) {
0088       x.print(ostr);
0089       return ostr;
0090     }
0091     /// @brief Fitted inclination angle
0092     double theta{0.};
0093     /// @brief Uncertainty on the fitted angle
0094     double dTheta{0.};
0095     /// @brief Fitted line intercept
0096     double y0{0.};
0097     /// @brief Uncertainty on the intercept
0098     double dY0{0.};
0099     /// @brief Evaluated chi2 postfit
0100     double chi2{0.};
0101     /// @brief Number of degrees of freedom
0102     std::size_t nDoF{0};
0103     /// @brief Number of iterations to converge
0104     std::size_t nIter{0};
0105   };
0106 
0107   /// @brief Expansion of the FitResult to also return information
0108   ///        on the fitted time offset t0
0109   struct FitResultT0 : public FitResult {
0110     void print(std::ostream& ostr) const override;
0111     /// @brief Fitted time offset t0
0112     double t0{0.};
0113     /// @brief Uncertainty on the time offset
0114     double dT0{0.};
0115   };
0116   /// @brief Fit a straight line to a set of straw measurements and return
0117   ///        the theta angle & the intercept as well as the associated
0118   ///        uncertainties.
0119   /// @param measuremments: Collection of straw measurements to fit
0120   /// @param signs: Convention of whether the straw-line shall be left or right
0121   ///               of a particular straw measurements. Needs to have the same
0122   ///               length as the list of measurements.
0123   template <CompositeSpacePointContainer StrawCont_t>
0124   std::optional<FitResult> fit(const StrawCont_t& measurements,
0125                                const std::vector<int>& signs) const;
0126 
0127   template <CompositeSpacePointContainer StripCont_t>
0128   std::optional<FitResult> fit(const StripCont_t& measurements,
0129                                const ResidualIdx projection) const;
0130 
0131   /// @brief Fit a straight line to a set of straw measurements taking the
0132   ///        time offset into account. The floating t0 parameter shrinks /
0133   ///        blows-up the drift radii of each straw measurement
0134   /// @param ctx: Reference to the experiment specific calibration context to
0135   ///             calculate the updated straw drift radius, velocity &
0136   ///             acceleration
0137   /// @param calibrator: Reference to the straw calibrator estimating the drift radius, velocity
0138   ///                    & acceleration.
0139   /// @param measuremments: Collection of straw measurements to fit
0140   /// @param signs: Convention of whether the straw-line shall be left or right
0141   ///               of a particular straw measurements. Needs to have the same
0142   ///               length as the list of measurements.
0143   /// @param startT0: Optional start estimate on the t0 parameter.
0144   template <CompositeSpacePointContainer StrawCont_t,
0145             CompositeSpacePointFastCalibrator<
0146                 Acts::RemovePointer_t<typename StrawCont_t::value_type>>
0147                 Calibrator_t>
0148   std::optional<FitResultT0> fit(
0149       const Acts::CalibrationContext& ctx, const Calibrator_t& calibrator,
0150       const StrawCont_t& measurements, const std::vector<int>& signs,
0151       std::optional<double> startT0 = std::nullopt) const;
0152 
0153  private:
0154   /// @brief Index of the drift circle covariance inside the straw
0155   ///        space-point's covariance array
0156   static constexpr auto s_covIdx = toUnderlying(ResidualIdx::bending);
0157 
0158   ///@brief Auxiliary struct to calculate the fast-fit constants
0159   struct FitAuxiliaries {
0160     /// @brief Default destructor
0161     virtual ~FitAuxiliaries() = default;
0162     /// @brief default constructor
0163     FitAuxiliaries() = default;
0164     /// @brief move constructor
0165     FitAuxiliaries(FitAuxiliaries&& other) = default;
0166     /// @brief Move assignment
0167     FitAuxiliaries& operator=(FitAuxiliaries&& other) = default;
0168     /// @brief Printer method
0169     virtual void print(std::ostream& ostr) const;
0170     /// @brief Ostream operator
0171     friend std::ostream& operator<<(std::ostream& ostr,
0172                                     const FitAuxiliaries& x) {
0173       x.print(ostr);
0174       return ostr;
0175     }
0176     ///@brief Tube position center - y weighted with inverse covariances
0177     double centerY{0.};
0178     ///@brief Tube position center - z weighted with inverse covariances
0179     double centerZ{0.};
0180     /// @brief Inverse covariances per straw measurement
0181     std::vector<double> invCovs{};
0182     ///@brief One over the sum of the inverse straw measurement covariances
0183     double covNorm{0.};
0184     ///@brief Expectation value of T_{z}^{2} - T_{y}^{2}
0185     double T_zzyy{0.};
0186     ///@brief Expectation value of T_{y} * T_{z}
0187     double T_yz{0.};
0188     ///@brief Expectation value of T_{z} * r
0189     double T_rz{0.};
0190     ///@brief Expectation value of T_{y} * r
0191     double T_ry{0.};
0192     ///@brief Prediced y0 given as the expectation value of the radii
0193     ///         divided by the inverse covariance sum.
0194     double fitY0{0.};
0195     /// @brief Number of degrees of freedom
0196     std::size_t nDoF{0};
0197   };
0198   ///@brief Extension of the auxiliary fit constants needed for the
0199   ///         seed refinement when T0 is floating
0200   struct FitAuxiliariesWithT0 : public FitAuxiliaries {
0201     ///@brief Constructor
0202     explicit FitAuxiliariesWithT0(FitAuxiliaries&& parent)
0203         : FitAuxiliaries{std::move(parent)} {}
0204     FitAuxiliariesWithT0(FitAuxiliariesWithT0&& other) = default;
0205 
0206     FitAuxiliariesWithT0& operator=(FitAuxiliariesWithT0&& other) = default;
0207     void print(std::ostream& ostr) const override;
0208     ///@brief Expectation value of T_{y} * v
0209     double T_vy{0.};
0210     ///@brief Expectation value of T_{z} * v
0211     double T_vz{0.};
0212     ///@brief Expectation value of T_{y} * a
0213     double T_ay{0.};
0214     ///@brief Expectation value of T_{z} * a
0215     double T_az{0.};
0216     ///@brief Expectation value of r * v
0217     double R_vr{0.};
0218     ///@brief Expectation value of v * v
0219     double R_vv{0.};
0220     ///@brief Expectation value of r * a
0221     double R_ar{0.};
0222     ///@brief First derivative of the fitted Y0
0223     double R_v{0.};
0224     ///@brief Second derivative of the ftted Y0
0225     double R_a{0.};
0226   };
0227   /// @brief Small Helper struct to calculate sin & cos of theta & 2 theta
0228   struct TrigonomHelper {
0229     /// @brief Constructor passing the theta angle
0230     explicit TrigonomHelper(const double theta)
0231         : cosTheta{std::cos(theta)},
0232           sinTheta{std::sin(theta)},
0233           cosTwoTheta{std::cos(2. * theta)},
0234           sinTwoTheta{std::sin(2. * theta)} {}
0235     double cosTheta{0.};
0236     double sinTheta{0.};
0237     double cosTwoTheta{0.};
0238     double sinTwoTheta{0.};
0239   };
0240   /// @brief Calculate the pure angular derivatives of the chi2 w.r.t theta
0241   /// @param angles: Wrapper object to pass sin & cos of theta
0242   /// @param fitPars: Constants of the current fit configuration
0243   /// @param thetaPrime: Mutable reference to which the first derivative is written
0244   /// @param thetaTwoPrime: Mutable reference to which the second derivative is written
0245   void calcAngularDerivatives(const TrigonomHelper& angles,
0246                               const FitAuxiliaries& fitPars, double& thetaPrime,
0247                               double& thetaTwoPrime) const;
0248   /// @brief Calculate the fit constants for a set of straw circles
0249   /// @param measurements: List of straw measurements
0250   /// @param signs: List of left/right signs to be annotated with each straw circle
0251   template <CompositeSpacePointContainer StrawCont_t>
0252   FitAuxiliaries fillAuxiliaries(const StrawCont_t& measurements,
0253                                  const std::vector<int>& signs) const;
0254   /// @brief Evaluate the chi2 after the fit has converged.
0255   /// @param measurements: List of straw measurements that have been used in the fit
0256   /// @param result: Mutable reference to the FitResult object. The object is used
0257   ///                to fetch the fit parameters and to store the chi2 result
0258   template <CompositeSpacePointContainer StrawCont_t>
0259   void calcPostFitChi2(const StrawCont_t& measurements,
0260                        FitResult& result) const;
0261   /// @brief Evaluate the chi2 after the fit with t0 has converged
0262   /// @param measurements: List of straw measurements
0263   /// @param calibrator: Reference to the calibrator used to retrieve an updated
0264   ///                    drift radius taking the fitted t0 into account
0265   /// @param result: Mutable reference to the FitResult object. The object is used
0266   ///                to fetch the fit parameters and to store the chi2 result
0267   template <CompositeSpacePointContainer StrawCont_t,
0268             CompositeSpacePointFastCalibrator<
0269                 Acts::RemovePointer_t<typename StrawCont_t::value_type>>
0270                 Calibrator_t>
0271   void calcPostFitChi2(const Acts::CalibrationContext& ctx,
0272                        const StrawCont_t& measurements,
0273                        const Calibrator_t& calibrator,
0274                        FitResultT0& result) const;
0275   /// @brief Calculate the chi2 term from a straw measurement given the known theta, y0
0276   /// @param angle: Helper struct caching the cosTheta & sinTheta
0277   /// @param y0: Intercept of the fitted line
0278   /// @param strawMeas: Reference to the straw measurement of interest
0279   /// @param r: Optional updated calibrated straw radius
0280   template <CompositeSpacePoint Point_t>
0281   double chi2Term(const TrigonomHelper& angle, const double y0,
0282                   const Point_t& strawMeas,
0283                   std::optional<double> r = std::nullopt) const;
0284   /// @brief Fit the track inclination angle and calculate the intercept
0285   ///        afterwards
0286   /// @param fitPars: Constants of the current measurement configuration
0287   std::optional<FitResult> fit(const FitAuxiliaries& fitPars) const;
0288   /// @brief Calculate the starting parameters on theta from the fit constants
0289   static double startTheta(const FitAuxiliaries& fitPars);
0290   /// @brief Calculate the time gradient of the chi2 w.r.t theta
0291   /// @param angles: Wrapper object to pass sin & cos of theta
0292   /// @param pars: Fit constants of the current fit configuration
0293   static double calcTimeGrad(const TrigonomHelper& angles,
0294                              const FitAuxiliariesWithT0& pars);
0295   /// @brief Fill the y0 and the uncertainties on theta and y0 in the result
0296   /// @param fitPars: Fit constants from the straw measurements
0297   /// @param thetaTwoPrime: Second derivative of the chi2 w.r.t theta
0298   /// @param result: Mutable reference to the FitResult object. The updated parameter are written
0299   ///                to this object
0300   void completeResult(const FitAuxiliaries& fitPars, const double thetaTwoPrime,
0301                       FitResult& result) const;
0302   /// @brief Enumeration to describe the outcome of the fit iteration
0303   enum class UpdateStatus {
0304     Converged,  ///< The fit converged
0305     Exceeded,   ////< Maximum number of iterations exceeded
0306     GoodStep,   ///< The fit did not converge yet
0307   };
0308   /// @brief Update the straw circle parameters for a fit with ts0 & theta
0309   /// @param fitPars: Fit constants from the straw measurements
0310   /// @param fitResult: Mutable reference to the FitResult object.
0311   ///                   The updated parameter are written to this object if the
0312   ///                   step was successful
0313   UpdateStatus updateIteration(const FitAuxiliariesWithT0& fitPars,
0314                                FitResultT0& fitResult) const;
0315   ///  @brief Calculate the extension of the fit constants needed for the simultaneous theta - t0 fit
0316   /// @param ctx: Reference to the experiment specific calibration context to
0317   ///             calculate the updated straw drift radius, velocity &
0318   ///             acceleration
0319   /// @param calibrator: Reference to the straw calibrator estimating the drift radius, velocity
0320   ///                    & acceleration.
0321   /// @param measuremments: Collection of straw measurements to fit
0322   /// @param measurements: List of straw measurements
0323   /// @param signs: List of left/right signs to be annotated with each straw circle
0324   /// @param t0: Current estimate on the time offset
0325   template <CompositeSpacePointContainer StrawCont_t,
0326             CompositeSpacePointFastCalibrator<
0327                 Acts::RemovePointer_t<typename StrawCont_t::value_type>>
0328                 Calibrator_t>
0329   FitAuxiliariesWithT0 fillAuxiliaries(const CalibrationContext& ctx,
0330                                        const Calibrator_t& calibrator,
0331                                        const StrawCont_t& measurements,
0332                                        const std::vector<int>& signs,
0333                                        const double t0) const;
0334 
0335   const Config m_cfg{};
0336   std::unique_ptr<const Acts::Logger> m_logger{};
0337 
0338   const Acts::Logger& logger() const;
0339 };
0340 
0341 }  // namespace Acts::Experimental::detail
0342 #include "Acts/Seeding/detail/FastStrawLineFitter.ipp"