|
ACTS
Experiment-independent tracking
|
Algorithms and data structures for track fitting.
Namespaces | |
| namespace | Acts::Experimental::Gx2fConstants |
| namespace | Acts::GsfConstants |
Classes | |
| class | Acts::AtlasBetheHeitlerApprox |
| This class approximates the Bethe-Heitler distribution as a gaussian mixture. More... | |
| class | Acts::BetheHeitlerApprox |
| class | Acts::BetheHeitlerApproxSingleCmp |
| This class approximates the Bethe-Heitler with only one component. More... | |
| class | Acts::Experimental::Gx2Fitter< propagator_t, traj_t > |
| Global Chi Square fitter (GX2F) implementation. More... | |
| struct | Acts::Experimental::Gx2FitterExtensions< traj_t > |
| Extension struct which holds delegates to customise the GX2F behaviour. More... | |
| struct | Acts::Experimental::Gx2FitterOptions< traj_t > |
| Combined options for the Global-Chi-Square fitter. More... | |
| struct | Acts::Experimental::Gx2FitterResult< traj_t > |
| struct | Acts::Experimental::Gx2fSystem |
| A container to manage all properties of a gx2f system. More... | |
| struct | Acts::Experimental::ScatteringProperties |
| A container to store scattering properties for each material surface. More... | |
| class | Acts::GainMatrixSmoother |
| Kalman trajectory smoother based on gain matrix formalism. More... | |
| class | Acts::GainMatrixUpdater |
| Kalman update step using the gain matrix formalism. More... | |
| struct | Acts::GaussianSumFitter< propagator_t, traj_t > |
| Gaussian Sum Fitter implementation. More... | |
| struct | Acts::GsfComponent |
| Encapsulates a component of a Gaussian mixture as used by the GSF. More... | |
| struct | Acts::GsfExtensions< traj_t > |
| The extensions needed for the GSF. More... | |
| struct | Acts::GsfOptions< traj_t > |
| class | Acts::KalmanFitter< propagator_t, traj_t > |
| Kalman fitter implementation. More... | |
| struct | Acts::KalmanFitterExtensions< traj_t > |
| Extension struct which holds delegates to customise the KF behavior. More... | |
| struct | Acts::KalmanFitterOptions< traj_t > |
| Combined options for the Kalman fitter. More... | |
| struct | Acts::KalmanFitterResult< traj_t > |
| class | Acts::MbfSmoother |
| Kalman trajectory smoother based on the Modified Bryson–Frazier (mBF) smoother. More... | |
Enumerations | |
| enum class | Acts::ComponentMergeMethod { Acts::ComponentMergeMethod::eMean , Acts::ComponentMergeMethod::eMaxWeight } |
| Available reduction methods for the reduction of a Gaussian mixture. More... | |
Functions | |
| template<typename track_state_t> | |
| void | Acts::Experimental::addMaterialToGx2fSums (Gx2fSystem &extendedSystem, const std::size_t nMaterialsHandled, const std::unordered_map< GeometryIdentifier, ScatteringProperties > &scatteringMap, const track_state_t &trackState, const Logger &logger) |
| Process material and fill the aMatrix and bVector. | |
| template<std::size_t kMeasDim, typename track_state_t> | |
| void | Acts::Experimental::addMeasurementToGx2fSums (Gx2fSystem &extendedSystem, const std::vector< BoundMatrix > &jacobianFromStart, const track_state_t &trackState, const Logger &logger) |
| Process measurements and fill the aMatrix and bVector. | |
| void | Acts::Experimental::addMeasurementToGx2fSumsBackend (Gx2fSystem &extendedSystem, const std::vector< BoundMatrix > &jacobianFromStart, const Eigen::MatrixXd &covarianceMeasurement, const BoundVector &predicted, const Eigen::VectorXd &measurement, const Eigen::MatrixXd &projector, const Logger &logger) |
| Adds a measurement to the GX2F equation system in a modular backend function. | |
| Eigen::VectorXd | Acts::Experimental::computeGx2fDeltaParams (const Gx2fSystem &extendedSystem) |
| Solve the gx2f system to get the delta parameters for the update. | |
| template<TrackProxyConcept track_proxy_t> | |
| std::size_t | Acts::Experimental::countMaterialStates (const track_proxy_t track, const std::unordered_map< GeometryIdentifier, ScatteringProperties > &scatteringMap, const Logger &logger) |
| Count the valid material states in a track for scattering calculations. | |
| template<TrackProxyConcept track_proxy_t> | |
| void | Acts::Experimental::fillGx2fSystem (const track_proxy_t track, Gx2fSystem &extendedSystem, const bool multipleScattering, const std::unordered_map< GeometryIdentifier, ScatteringProperties > &scatteringMap, std::vector< GeometryIdentifier > &geoIdVector, const Logger &logger) |
| Fill the GX2F system with data from a track. | |
| void | Acts::reduceMixtureLargestWeights (std::vector< GsfComponent > &cmpCache, std::size_t maxCmpsAfterMerge, const Surface &surface) |
| Very simple mixture reduction method: Just removes the components with the smallest weight until the required number of components is reached. | |
| void | Acts::reduceMixtureWithKLDistance (std::vector< GsfComponent > &cmpCache, std::size_t maxCmpsAfterMerge, const Surface &surface) |
| Greedy component reduction algorithm. | |
| void | Acts::Experimental::updateGx2fCovarianceParams (BoundMatrix &fullCovariancePredicted, Gx2fSystem &extendedSystem) |
| Calculate and update the covariance of the fitted parameters. | |
| void | Acts::Experimental::updateGx2fParams (BoundTrackParameters ¶ms, const Eigen::VectorXd &deltaParamsExtended, const std::size_t nMaterialSurfaces, std::unordered_map< GeometryIdentifier, ScatteringProperties > &scatteringMap, const std::vector< GeometryIdentifier > &geoIdVector) |
| Update parameters (and scattering angles if applicable). | |
|
strong |
| void Acts::Experimental::addMaterialToGx2fSums | ( | Gx2fSystem & | extendedSystem, |
| const std::size_t | nMaterialsHandled, | ||
| const std::unordered_map< GeometryIdentifier, ScatteringProperties > & | scatteringMap, | ||
| const track_state_t & | trackState, | ||
| const Logger & | logger ) |
Process material and fill the aMatrix and bVector.
The function processes each material for the GX2F Actor fitting process. It extracts the information from the track state and adds it to aMatrix, bVector, and chi2sum.
| track_state_t | The type of the track state |
| extendedSystem | All parameters of the current equation system |
| nMaterialsHandled | How many materials we already handled. Used for the offset. |
| scatteringMap | The scattering map, containing all scattering angles and covariances |
| trackState | The track state to analyse |
| logger | A logger instance |
| void Acts::Experimental::addMeasurementToGx2fSums | ( | Gx2fSystem & | extendedSystem, |
| const std::vector< BoundMatrix > & | jacobianFromStart, | ||
| const track_state_t & | trackState, | ||
| const Logger & | logger ) |
Process measurements and fill the aMatrix and bVector.
The function processes each measurement for the GX2F Actor fitting process. It extracts the information from the track state and adds it to aMatrix, bVector, and chi2sum.
| kMeasDim | Number of dimensions of the measurement |
| track_state_t | The type of the track state |
| extendedSystem | All parameters of the current equation system to update |
| jacobianFromStart | The Jacobian matrix from start to the current state |
| trackState | The track state to analyse |
| logger | A logger instance |
| void Acts::Experimental::addMeasurementToGx2fSumsBackend | ( | Gx2fSystem & | extendedSystem, |
| const std::vector< BoundMatrix > & | jacobianFromStart, | ||
| const Eigen::MatrixXd & | covarianceMeasurement, | ||
| const BoundVector & | predicted, | ||
| const Eigen::VectorXd & | measurement, | ||
| const Eigen::MatrixXd & | projector, | ||
| const Logger & | logger ) |
Adds a measurement to the GX2F equation system in a modular backend function.
This function processes measurement data and integrates it into the GX2F system.
| extendedSystem | All parameters of the current equation system to update. |
| jacobianFromStart | The Jacobian matrix from the start to the current state. |
| covarianceMeasurement | The covariance matrix of the measurement. |
| predicted | The predicted state vector based on the track state. |
| measurement | The measurement vector. |
| projector | The projection matrix. |
| logger | A logger instance. |
| Eigen::VectorXd Acts::Experimental::computeGx2fDeltaParams | ( | const Gx2fSystem & | extendedSystem | ) |
Solve the gx2f system to get the delta parameters for the update.
This function computes the delta parameters for the GX2F Actor fitting process by solving the linear equation system [a] * delta = b. It uses the column-pivoting Householder QR decomposition for numerical stability.
| extendedSystem | All parameters of the current equation system |
| std::size_t Acts::Experimental::countMaterialStates | ( | const track_proxy_t | track, |
| const std::unordered_map< GeometryIdentifier, ScatteringProperties > & | scatteringMap, | ||
| const Logger & | logger ) |
Count the valid material states in a track for scattering calculations.
This function counts the valid material surfaces encountered in a track by examining each track state. The count is based on the presence of material flags and the availability of scattering information for each surface.
| track_proxy_t | The type of the track proxy |
| track | A constant track proxy to inspect |
| scatteringMap | Map of geometry identifiers to scattering properties, containing scattering angles and validation status |
| logger | A logger instance |
| void Acts::Experimental::fillGx2fSystem | ( | const track_proxy_t | track, |
| Gx2fSystem & | extendedSystem, | ||
| const bool | multipleScattering, | ||
| const std::unordered_map< GeometryIdentifier, ScatteringProperties > & | scatteringMap, | ||
| std::vector< GeometryIdentifier > & | geoIdVector, | ||
| const Logger & | logger ) |
Fill the GX2F system with data from a track.
This function processes a track proxy and updates the aMatrix, bVector, and chi2 values for the GX2F fitting system. It considers material only if multiple scattering is enabled.
| track_proxy_t | The type of the track proxy |
| track | A constant track proxy to inspect |
| extendedSystem | All parameters of the current equation system |
| multipleScattering | Flag to consider multiple scattering in the calculation |
| scatteringMap | Map of geometry identifiers to scattering properties, containing scattering angles and validation status |
| geoIdVector | A vector to store geometry identifiers for tracking processed elements |
| logger | A logger instance |
| void Acts::reduceMixtureLargestWeights | ( | std::vector< GsfComponent > & | cmpCache, |
| std::size_t | maxCmpsAfterMerge, | ||
| const Surface & | surface ) |
Very simple mixture reduction method: Just removes the components with the smallest weight until the required number of components is reached.
| cmpCache | the component collection |
| maxCmpsAfterMerge | the number of components we want to reach |
| surface | the surface type on which the components are (unused here) |
| void Acts::reduceMixtureWithKLDistance | ( | std::vector< GsfComponent > & | cmpCache, |
| std::size_t | maxCmpsAfterMerge, | ||
| const Surface & | surface ) |
Greedy component reduction algorithm.
Reduces the components with the minimal symmetric KL-distance (applied only to the q/p-dimension) until the required number of components is reached.
| cmpCache | the component collection |
| maxCmpsAfterMerge | the number of components we want to reach |
| surface | the surface type on which the components are |
| void Acts::Experimental::updateGx2fCovarianceParams | ( | BoundMatrix & | fullCovariancePredicted, |
| Gx2fSystem & | extendedSystem ) |
Calculate and update the covariance of the fitted parameters.
This function calculates the covariance of the fitted parameters using cov = inv([a]) It then updates the first square block of size ndfSystem. This ensures, that we only update the covariance for fitted parameters. (In case of no qop/time fit)
| fullCovariancePredicted | The covariance matrix to update |
| extendedSystem | All parameters of the current equation system |
| void Acts::Experimental::updateGx2fParams | ( | BoundTrackParameters & | params, |
| const Eigen::VectorXd & | deltaParamsExtended, | ||
| const std::size_t | nMaterialSurfaces, | ||
| std::unordered_map< GeometryIdentifier, ScatteringProperties > & | scatteringMap, | ||
| const std::vector< GeometryIdentifier > & | geoIdVector ) |
Update parameters (and scattering angles if applicable).
| params | Parameters to be updated |
| deltaParamsExtended | Delta parameters for bound parameter and scattering angles |
| nMaterialSurfaces | Number of material surfaces in the track |
| scatteringMap | Map of geometry identifiers to scattering properties, containing all scattering angles and covariances |
| geoIdVector | Vector of geometry identifiers corresponding to material surfaces |