Warning, /acts/docs/core/reconstruction/pattern_recognition/seeding.md is written in an unsupported language. File is not indexed.
0001 (seeding_core)=
0002
0003 # Seeding
0004
0005 To reduce the time needed to reconstruct particle tracks, a track seed
0006 (henceforth: seed) is created which serves as the initial direction for the track
0007 reconstruction algorithm (henceforth: the tracking). The tracking then tries to
0008 find all measurements belonging to a single particle in this direction in order
0009 to reconstruct the track. This means, if no seed exists for a particle, this
0010 particle will not be reconstructed. On the other hand, finding too many seeds
0011 which either correspond to particles with already existing seeds or
0012 which do not correspond to particles at all increases the time needed for
0013 tracking.
0014
0015 A good seeding algorithm, therefore, has the following properties:
0016
0017 * It finds at least one seed for each particle that should be found
0018 * It doesn't find many seeds which do NOT correspond to particles
0019 * It doesn't find many seeds per particle
0020
0021 The most typical way to create seeds is to combine measurements. In a homogeneous magnetic field, 3 measurements perfectly describe the helical path of a charged particle. One such triplet of measurements would then constitute a seed and defines, in close bounds, where the tracking needs to look for additional measurements to create a track spanning the whole detector. The difficulty is in choosing the correct measurements, as a helix can be fitted through any 3 measurements in a collision event with potentially tens of thousands of measurements. Therefore, many constraints or “cuts” are defined to reduce the number of candidates. Cuts may define where particles originate or the range of energy of particles to be found or otherwise restrict the combination of measurements for seed creation.
0022
0023 ## ACTS Implementation
0024
0025 The seeding implementation in `Core/include/Acts/Seeding/` was written with a focus on parallelism and maintainability and as detector agnostic as possible, only assuming a (near) homogeneous magnetic field with particles originating from the central detector region. Cuts are configurable and can be plugged in as an algorithm which is called by the seeding. The seeding works on measurements or “SpacePoints” (SP), which need to provide $(x,y,z)$ coordinates with the $z$ axis being along the magnetic field, and $x$ and $y$. For the seeding algorithm to function the particle point of origin has to have a radius smaller than the radius of the detector layer closest to the interaction region. In other words, this seeding algorithm is not suitable for secondary particles originating far from the interaction region.
0026
0027 :::{figure} figures/seeding/3Dcoordinates.svg
0028 :width: 550px
0029 :align: center
0030 Sketch of the detector with 3 layers. The interaction region is supposed to be located along the $z$-axis and has a size significantly smaller than the radius of the innermost detector layer.
0031 :::
0032
0033 :::{attention}
0034 Note that the seeding algorithm breaks down for particles with a particle
0035 track whose helix diameter is smaller than the detector radius until which
0036 seeds are to be created. This is due to ordering assumptions of SP
0037 locations as well as due to approximations which become inaccurate for
0038 lower energy particles.
0039 :::
0040
0041 ### SP Grid and Groups Formation
0042
0043 The SPs in each detector layer are projected on a rectangular grid of
0044 configurable granularity. The search for seed starts by selecting SP in the
0045 middle detector layer. Then matching SPs are searched in the inner and outer
0046 layers. Grouping of the SPs in the aforementioned grid allows to limit the
0047 search to neighbouring grid cells thus improving significantly algorithm
0048 performance (see {numref}`tripletsFormation`). The number of neighboring bins
0049 used in the SP search can be defined separately for the bottom and top layer
0050 SPs in the $z$ and $\phi$ directions.
0051
0052 (tripletsFormation)=
0053 :::{figure} figures/seeding/tripletsFormation.svg
0054 :width: 450px
0055 :align: center
0056 Representation of the search for triplet combinations in the $(r, z)$ plane. The bins used in the search are represented in different colours.
0057 :::
0058
0059 ### The Seed Finder
0060
0061 The {func}`Acts::SeedFinder::createSeedsForGroup` function receives three iterators
0062 over SPs constructed from detector layers of increasing radii. The seedfinder will
0063 then attempt to create seeds, with each seed containing exactly one SP returned by
0064 each of the three iterators. It starts by iterating over SPs in the middle layer
0065 (2nd iterator), and within this loop separately iterates once over the bottom SP
0066 and once over the top SP. Within each of the nested loops, SP pairs are tested for
0067 compatibility by applying a set of configurable cuts that can be tested with
0068 two SP only (pseudorapidity, origin along $z$-axis, distance in $r$ between SP,
0069 compatibility with interaction point).
0070
0071 :::{doxygenfunction} Acts::SeedFinder::createSeedsForGroup(const Acts::SeedFinderOptions &options, SeedingState &state, const grid_t &grid, container_t &outputCollection, const sp_range_t &bottomSPs, const std::size_t middleSPs, const sp_range_t &topSPs, const Acts::Range1D<float> &rMiddleSPRange) const
0072 :::
0073
0074 For all pairs passing the selection the triplets of bottom-middle-top SPs are formed.
0075 Each triplet is then confronted with the helix hypothesis. In order to perform calculations
0076 only once, the circle calculation is spread out over the three loops.
0077
0078 (x-yCoordinates)=
0079 :::{figure} figures/seeding/x-yCoordinates.svg
0080 :width: 400px
0081 :align: center
0082 The x-y projection of the detector with the charged particle helical track originating from the centre of the detector. Signals left by the passage of the track through the detector layers are marked with green crosses.
0083 :::
0084
0085 From the helix circle (see {numref}`x-yCoordinates`), particle energy and impact parameters can be estimated.
0086 To calculate the helix circle in the $x/y$ plane, the $x,y$ coordinates are
0087 transformed into a $u/v$ plane to calculate the circle with a linear equation
0088 instead of a quadratic equation for speed. The conformal transformation is given by:
0089
0090 \begin{equation*}
0091 u = \frac{x}{x^2+y^2}, \quad \quad v = \frac{y}{x^2+y^2} ,
0092 \end{equation*}
0093
0094 where the circle containing the three SPs are transformed into a line with equation $v = Au + B$.
0095 The angular coefficient $A$ can be evaluated by the slope of the linear function between the top and bottom layer SPs, after transforming the coordinates of these SPs from $x/y$ to $u/v$ using the previous equations:
0096
0097 \begin{equation*}
0098 A = \frac{v_t-v_b}{u_t-u_b} ,
0099 \end{equation*}
0100
0101 Then, $B$ can be obtained by inserting $A$ into the linear equation for the bottom layer SP:
0102 \begin{equation*}
0103 v_b = Au_b + B \rightarrow B = v_b - Au_B .
0104 \end{equation*}
0105
0106 Inserting the coefficients in the circle equation and assuming that the circle goes through the origin we obtain:
0107
0108 \begin{equation*}
0109 (2R)^2 = \frac{A^2+1}{B^2} .
0110 \end{equation*}
0111
0112 Now we can apply a cut on the estimate of the minimum helix diameter `minHelixDiameter2` without the extra overhead of conversions or computationally complex calculations. The seed is accepted if
0113
0114 \begin{equation*}
0115 \frac{A^2+1}{B^2} > (2 R^{min})^2 = \left ( \frac{2 \cdot p_T^{min}}{300 \cdot B_z} \right)^2 \equiv \textnormal{minHelixDiameter2} ,
0116 \end{equation*}
0117
0118 where $B_z$ is the magnetic field.
0119
0120 (r-zCoordinates)=
0121 :::{figure} figures/seeding/r-zCoordinates.svg
0122 :width: 500px
0123 :align: center
0124 The r-z projection of the detector with the same charged particle track. The track is depicted with the same colours as in the previous figure.
0125 :::
0126
0127 The track is not an ideal helix. At each detector layer (or any other material)
0128 scattering may occur making the helix approximate.
0129 The algorithm will check if the triplet forms a nearly straight line
0130 in the $r/z$ plane (see {numref}`r-zCoordinates`) as the particle path in the $r/z$ plane is
0131 unaffected by the magnetic field. This is split into two parts; the first test occurs before the calculation of the helix
0132 circle. Therefore, the deviation from a straight line is compared to the
0133 maximum allowed scattering at minimum $p_T$ scaled by the forward angle:
0134
0135 \begin{equation*}
0136 \left (\cot \theta_b - \cot \theta_t \right )^2 < \sigma^2_{p_T^{min}} + \sigma_f^2,
0137 \end{equation*}
0138
0139 This check takes into account the squared uncertainty in the difference between slopes ($\sigma_f^2$) and the
0140 scattering term `scatteringInRegion2` ($\sigma^2_{p_T^{min}}$), which is calculated from
0141 `sigmaScattering`, the configurable number of sigmas of scattering angle
0142 to be considered, and `maxScatteringAngle2`, which is evaluated from the
0143 Lynch & Dahl correction[^Tanabashi:2018]$^,$[^Lynch:1991] of the Highland equation assuming the lowest
0144 allowed $p_T$:
0145
0146 \begin{equation*}
0147 \Theta_0^{min} = \frac{13.6 \text{MeV}}{p_T^{min}} \cdot \sqrt{\frac{z_q^2 L}{\beta^2 L_0}} \left(1+0.038 \ln \left(\frac{z_q^2 L}{\beta^2 L_0}\right) \right),
0148 \end{equation*}
0149
0150 The second part check against the multiple scattering term `p2scatterSigma` ($\sigma^2_{p_T^{estimated}}$) assuming the seed $p_T$ estimation, instead of the minimum allowed $p_T$.
0151 This term is inversely proportional to the momentum and accounts for the curvature of the seed; smaller scattering angle is permitted for higher momentum.
0152
0153 \begin{equation*}
0154 \left (\cot \theta_b - \cot \theta_t \right ) ^2 < \sigma^2_{p_T^{estimated}} + \sigma_f^2,
0155 \end{equation*}
0156
0157 The last cut applied in this function is on the transverse impact parameter (or DCA -
0158 distance of closest approach), which is the distance of the perigee of a track from
0159 the interaction region in $mm$ of detector radius. It is calculated and cut on
0160 before storing all top SP compatible with both the current middle SP and current
0161 bottom SP.
0162
0163 (impactParameter)=
0164 :::{figure} figures/seeding/impactParameter.svg
0165 :width: 400px
0166 :align: center
0167 Helix representation in $x/y$ reference frame with central space-point (SP$_m$) in the origin.
0168 :::
0169
0170 Assuming the middle layer SP is in the origin of the $x/y$ frame, as in {numref}`impactParameter`.
0171 The distance between the centre of the helix and the interaction point (IP) is given by
0172 \begin{equation*}
0173 (x_0 + r_m)^2 + y_0^2 = (R + d_0)^2 \quad \xrightarrow{R^2 = x_0^2 + y_0^2} \quad \frac{d_0^2}{R^2} + 2 \frac{d_0}{R} = \frac{2 x_0 r_m + r_m^2}{R^2} .
0174 \end{equation*}
0175
0176 Considering that $d_0 << R$ (we can neglect the term proportional to $d_0^2$) and using the $u/v$ line equation calculated previously,
0177 the cut can now be estimated using a linear function in the $u/v$ plane instead of a quartic function:
0178
0179 \begin{equation*}
0180 d_0 \leq \left| \left( A - B \cdot r_M \right) \cdot r_M \right|
0181 \end{equation*}
0182
0183 ### The Seed Filter
0184
0185 After creating the potential seeds we apply a seed filter procedure that compares the seeds with other SPs compatible with the seed curvature.
0186 This process ranks the potential seeds based on certain quality criteria and selects the ones that are more likely to produce high-quality tracks
0187 The filter is divided into two functions {func}`Acts::SeedFilter::filterSeeds_2SpFixed` and {func}`Acts::SeedFilter::filterSeeds_1SpFixed`.
0188
0189 The first function compares the middle and bottom layer SPs of the seeds to other top layer SPs; seeds only differing in top SP are
0190 compatible if they have similar helix radius with the same sign (i.e. the same charge). The SPs must have a minimum distance in
0191 detector radius, such that SPs from the same layer cannot be considered compatible. The second function iterates over the seeds with
0192 only a common middle layer SP and selects the higher quality combinations.
0193
0194 :::{doxygenfunction} Acts::SeedFilter::filterSeeds_2SpFixed
0195 :outline:
0196 :::
0197
0198 This function assigns a weight (which should correspond to the likelihood that
0199 a seed is good) to all seeds and applies detector-specific selection of seeds based on weights.
0200 The weight is a “soft cut”, which means that it is only
0201 used to discard tracks if many seeds are created for the same middle SP.
0202 This process is important to improving computational
0203 performance and the quality of the final track collections by rejecting lower-quality seeds.
0204
0205 The weight is calculated by:
0206 \begin{equation*}
0207 w = (c_1 \cdot N_{t} - c_2 \cdot d_0 - c_3 |z_0| ) + \textnormal{detector specific cuts}.
0208 \end{equation*}
0209
0210 The transverse ($d_0$) and longitudinal ($z_0$) impact parameters are multiplied by a configured factor and subtracted from
0211 the weight, as seeds with higher impact parameters are assumed to be less
0212 likely to stem from a particle than another seed using the same middle SP with
0213 smaller impact parameters. The number of compatible seeds ($N_t$) is used to increase the weight, as a higher number of measurements
0214 will lead to higher quality tracks. Finally, the weight can also be affected by optional detector-specific cuts.
0215
0216 The {func}`Acts::SeedFilter::filterSeeds_2SpFixed` function also includes a configurable {struct}`Acts::SeedConfirmationRangeConfig` seed confirmation step that, when enabled,
0217 classifies higher quality seeds as "quality confined" seeds if they fall within a predefined range of parameters ($d_0$, $z_0$ and $N_t$) that also
0218 depends on the region of the detector (i.e., forward or central region). If the seed is not
0219 classified as "quality confined" seed, it will only be accepted if its weight is greater
0220 than a certain threshold and no other high-quality seed has been found.
0221
0222 The seed confirmation also sets a limit on the number of seeds produced for each middle SP,
0223 which retains only the higher quality seeds. If this limit is exceeded, the algorithm
0224 checks if there is any low-quality seed in the seed container of this middle SP that can be removed.
0225
0226
0227 :::{doxygenfunction} Acts::SeedFilter::filterSeeds_1SpFixed (Acts::SpacePointMutableData& , CandidatesForMiddleSp<const external_spacepoint_t>&, collection_t&) const
0228 :outline:
0229 :::
0230
0231 This function allows the detector-specific cuts to filter based on all
0232 seeds with a common middle SP and limits the number of seeds per middle SP to
0233 the configured limit. It sorts the seeds by weight and, to achieve a
0234 well-defined ordering in the rare case weights are equal, sorts them by
0235 location. The ordering by location is only done to make sure reimplementations
0236 (such as the GPU code) are comparable and return the bitwise exactly the same
0237 result.
0238
0239 [^Tanabashi:2018]: M. Tanabashi et al. (Particle Data Group), Passage of Particles Through Matter, Phys. Rev. D 98 0300001 (2018) 2.
0240
0241 [^Lynch:1991]: G.R. Lynch and O.I Dahl, Nucl. Instrum. Methods B58, Phys. Rev. D 98 6 (1991) 2.