|
|
|||
Warning, file /acts/Core/include/Acts/Seeding/CompositeSpacePointLineSeeder.hpp was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
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/EventData/CompositeSpacePoint.hpp" 0012 #include "Acts/Seeding/detail/CompSpacePointAuxiliaries.hpp" 0013 #include "Acts/Utilities/CalibrationContext.hpp" 0014 #include "Acts/Utilities/Delegate.hpp" 0015 0016 namespace Acts::Experimental { 0017 namespace detail { 0018 /// @brief Concept for a utility class to fill the space points from 0019 /// one container to another. The filler is templated over the 0020 /// SpacePoint type from the uncalibrated source container over 0021 /// the target space point container type, which may differ 0022 /// The filler is responsible for the conversion between the two 0023 /// space point formats and to apply a re-calibation if needed 0024 template <typename SeedFiller_t, typename UnCalibSp_t, typename CalibCont_t> 0025 concept CompositeSpacePointSeedFiller = 0026 CompositeSpacePoint<UnCalibSp_t> && 0027 CompositeSpacePointContainer<CalibCont_t> && 0028 requires(const SeedFiller_t& filler, const CalibrationContext& cctx, 0029 const Vector3& pos, const Vector3& dir, const double t0, 0030 const UnCalibSp_t& testSp, CalibCont_t& seedContainer, 0031 const std::size_t lowerLayer, const std::size_t upperLayer) { 0032 /// @brief Utility function to choose the good straw space points 0033 /// for seeding 0034 /// @param testSp: Reference to the straw-type measurement to test 0035 { filler.goodCandidate(testSp) } -> std::same_as<bool>; 0036 /// @brief Calculates the pull of the space point w.r.t. to the 0037 /// candidate seed line. To improve the pull's precision 0038 /// the function may call the calibrator in the backend 0039 /// @param cctx: Reference to the calibration context to pipe 0040 /// the hook for conditions access to the caller 0041 /// @param pos: Position of the cancidate seed line 0042 /// @param dir: Direction of the candidate seed line 0043 /// @param t0: Offse in the time of arrival of the particle 0044 /// @param testSp: Reference to the straw space point to test 0045 { 0046 filler.candidatePull(cctx, pos, dir, t0, testSp) 0047 } -> std::same_as<double>; 0048 /// @brief Creates a new empty container to construct a new segment seed 0049 /// @param cctx: Calibration context in case that the container shall be 0050 /// part of a predefined memory block 0051 { filler.newContainer(cctx) } -> std::same_as<CalibCont_t>; 0052 /// @brief Appends the candidate candidate space point to the segment 0053 /// seed container & optionally calibrates the parameters 0054 /// @param cctx: Reference to the calibration context to pipe 0055 /// the hook for conditions access to the caller 0056 /// @param pos: Position of the cancidate seed line 0057 /// @param dir: Direction of the candidate seed line 0058 /// @param t0: Offse in the time of arrival of the particle 0059 /// @param testSp: Reference to the straw space point to test 0060 /// @param seedContainer: Reference to the target container to 0061 /// which the space point is appended to 0062 { 0063 filler.append(cctx, pos, dir, t0, testSp, seedContainer) 0064 } -> std::same_as<void>; 0065 /// @brief Helper method to send a stop signal to the line seeder, if for instance, 0066 /// the two layers are too close to each other. The method is 0067 /// called after every layer update. If true is returned no seeds 0068 /// are produced further 0069 /// @param lowerLayer: Index of the current lower straw layer 0070 /// @param upperLayer: Index of the current upper straw layer 0071 { filler.stopSeeding(lowerLayer, upperLayer) } -> std::same_as<bool>; 0072 }; 0073 /// @brief Define the concept of the space point measurement sorter. The sorter shall take a collection 0074 /// of station space points and split them into straw && strip hits. 0075 /// Hits 0076 /// from each category are then subdivided further into the particular detector 0077 /// layers. 0078 /// 0079 /// A possible implementation of the CompositeSpacePointSorter needs to 0080 /// have 0081 /// the following attributes 0082 /// 0083 /// using SpVec_t = Standard container satisfiyng the 0084 /// CompositeSpacePointContainer concept 0085 /// 0086 /// 0087 /// const std::vector<SpVec_t>& strawHits(); 0088 /// const std::vector<SpVec_t>& stripHits(); 0089 /// Each SpVec_t contains all measurements from a particular detector layer 0090 template <typename Splitter_t, typename SpacePointCont_t> 0091 concept CompositeSpacePointSorter = 0092 CompositeSpacePointContainer<SpacePointCont_t> && 0093 requires(const Splitter_t& sorter) { 0094 /// @brief Return the straw-hit space point sorted by straw layer 0095 { 0096 sorter.strawHits() 0097 } -> std::same_as<const std::vector<SpacePointCont_t>&>; 0098 /// @brief Return the strip-hit space points sorted by detector layer 0099 { 0100 sorter.stripHits() 0101 } -> std::same_as<const std::vector<SpacePointCont_t>&>; 0102 }; 0103 0104 template <typename SeedAuxiliary_t, typename UnCalibCont_t, 0105 typename CalibCont_t> 0106 concept CompSpacePointSeederDelegate = 0107 CompositeSpacePointSeedFiller< 0108 SeedAuxiliary_t, RemovePointer_t<typename UnCalibCont_t::value_type>, 0109 CalibCont_t> && 0110 CompositeSpacePointSorter<SeedAuxiliary_t, UnCalibCont_t>; 0111 0112 } // namespace detail 0113 0114 /// @brief Initial line parameters from a pattern recognition like 0115 /// the Hough transform are often not suitable for a line fit 0116 /// as the resolution of the hough bins usually exceeds the size 0117 /// of the straws. 0118 /// The CompositeSpacePointLineSeeder refines the parameters 0119 /// and the selected measurements such that both become 0120 /// candidates for a stright line fit. The user needs to 0121 /// split the straw measurements per logical straw layer. 0122 /// Further, the interface needs to provide some auxiliary 0123 /// methods to interact with an empty space point container 0124 /// & to calculate a calbrated candidate pull. 0125 /// From these ingredients, the `CompositeSpacePointLineSeeder` 0126 /// iterates from the outermost layers at both ends and tries 0127 /// to construct new candidates. Tangent lines are constructed 0128 /// to a pair of circles from each seeding layer and then straw 0129 /// measurements from the other layers are tried to be added 0130 /// onto the line. If the number of straws exceed the threshold, 0131 /// compatible strip measurements from each strip layer are added. 0132 class CompositeSpacePointLineSeeder { 0133 public: 0134 /// @brief Configuration of the cuts to sort out generated 0135 /// seeds with poor quality. 0136 struct Config { 0137 /// @brief Cut on the theta angle 0138 std::array<double, 2> thetaRange{0, 0}; 0139 /// @brief Cut on the intercept range 0140 std::array<double, 2> interceptRange{0, 0}; 0141 0142 /// @brief Upper cut on the hit chi2 w.r.t. seed in order to be associated to the seed 0143 double hitPullCut{5.}; 0144 /// @brief How many drift circles may be on a layer to be used for seeding 0145 std::size_t busyLayerLimit{2}; 0146 /// @brief Layers may contain measurements with bad hits and hence the 0147 bool busyLimitCountGood{true}; 0148 0149 /// @brief How many drift circle hits needs the seed to contain in order to be valid 0150 std::size_t nStrawHitCut{3}; 0151 /// @brief Hit cut based on the fraction of collected tube layers. 0152 /// The seed must pass the tighter of the two requirements. 0153 double nStrawLayHitCut{2. / 3.}; 0154 /// @brief Once a seed with even more than initially required hits is found, 0155 /// reject all following seeds with less hits 0156 bool tightenHitCut{true}; 0157 /// @brief Check whether a new seed candidate shares the same left-right solution with already accepted ones 0158 /// Reject the seed if it has the same amount of hits 0159 bool overlapCorridor{true}; 0160 }; 0161 /// @brief Class constructor 0162 /// @param cfg Reference to the seeder configuration object 0163 /// @param logger Logger object used for debug print out 0164 explicit CompositeSpacePointLineSeeder( 0165 const Config& cfg, 0166 std::unique_ptr<const Logger> logger = getDefaultLogger( 0167 "CompositeSpacePointLineSeeder", Logging::Level::INFO)); 0168 /// @brief Return the configuration object of the seeder 0169 const Config& config() const { return m_cfg; } 0170 /// @brief Use the assignment of the parameter indices from the CompSpacePointAuxiliaries 0171 using ParIdx = detail::CompSpacePointAuxiliaries::FitParIndex; 0172 /// @brief Use the vector from the CompSpacePointAuxiliaires 0173 using Vector = detail::CompSpacePointAuxiliaries::Vector; 0174 /// @brief Vector containing the 5 straight segment line parameters 0175 using SeedParam_t = std::array<double, toUnderlying(ParIdx::nPars)>; 0176 /// @brief Abrivation of the straight line. The first element is the 0177 /// reference position and the second element is the direction 0178 using Line_t = std::pair<Vector, Vector>; 0179 0180 /// @brief Enumeration to pick one of the four tangent lines to 0181 /// the straw circle pair. 0182 enum class TangentAmbi : std::uint8_t { 0183 RR = 0, //< Both circles are on the right side 0184 RL = 1, //< The top circle is on the right and the bottom circle on the 0185 // left < side 0186 LR = 2, //< The top circle is on the left and the bottom circle on the 0187 //< right side 0188 LL = 3, //< Both circles are on the left side 0189 }; 0190 0191 /// @brief Helper struct describing the line parameters that are 0192 /// tangential to a pair of straw measurements 0193 struct TwoCircleTangentPars { 0194 /// @brief Estimated angle 0195 double theta{0.}; 0196 /// @brief Estimated intercept 0197 double y0{0.}; 0198 /// @brief Uncertainty on the angle 0199 double dTheta{0.}; 0200 /// @brief Uncertainty on the intercept 0201 double dY0{0.}; 0202 /// @brief Flag indicating which solution is constructed 0203 TangentAmbi ambi{TangentAmbi::LL}; 0204 /// @brief Default destructor 0205 virtual ~TwoCircleTangentPars() = default; 0206 0207 /// @brief Definition of the print operator 0208 /// @param ostr: Mutable reference to the stream to print to 0209 /// @param pars: The parameters to be printed 0210 friend std::ostream& operator<<(std::ostream& ostr, 0211 const TwoCircleTangentPars& pars) { 0212 pars.print(ostr); 0213 return ostr; 0214 } 0215 0216 protected: 0217 /// @brief Actual implementation of the printing 0218 /// @param ostr: Mutable reference to the stream to print to 0219 virtual void print(std::ostream& ostr) const; 0220 }; 0221 0222 /// @brief Converts the line tangent ambiguity into a string 0223 static std::string toString(const TangentAmbi ambi); 0224 /// @brief Translate the combination of two drift signs into the proper 0225 /// tangent ambiguity enum value 0226 /// @param signTop: Left/right sign of the top straw tube 0227 /// @param signBottom: Left/right sign of the bottom straw tube 0228 static TangentAmbi encodeAmbiguity(const int signTop, const int signBottom); 0229 /// @brief Construct the line that is tangential to a pair of two straw circle measurements 0230 /// @param topHit: First straw hit 0231 /// @param bottomHit: Second straw hit 0232 /// @param ambi: Left right ambiguity of the bottom & top hit 0233 template <CompositeSpacePoint Sp_t> 0234 static TwoCircleTangentPars constructTangentLine(const Sp_t& topHit, 0235 const Sp_t& bottomHit, 0236 const TangentAmbi ambi); 0237 /// @brief Creates the direction vector from the reference hit used to 0238 /// construct the tangent seed and the result on theta 0239 /// @param refHit: Reference hit to define the local axes (Bottom hit) 0240 /// @param tanAngle: Theta value from the TwoCircleTangentPars 0241 template <CompositeSpacePoint Spt_t> 0242 static Vector makeDirection(const Spt_t& refHit, const double tanAngle); 0243 0244 private: 0245 /// @brief Cache object of a constructed & valid seed solution. 0246 /// It basically consists out of the generated parameters. 0247 /// the straw hits contributing to the seed & the left/right 0248 /// ambiguities given the parameters of the solutions. 0249 /// To avoid the copy of the memory, the hits are encoded as a 0250 /// pair of indices representing the straw layer & hit number. 0251 template <CompositeSpacePointContainer UnCalibCont_t, 0252 detail::CompositeSpacePointSorter<UnCalibCont_t> Splitter_t> 0253 struct SeedSolution : public TwoCircleTangentPars { 0254 /// @brief Constructor taking the constructed tangential parameters & 0255 /// the pointer to the splitter to associate the hits to the seed 0256 /// @param pars: Theta & intercept describing the tangential line 0257 /// @param layerSorter: Pointer to the sorter object carrying a sorted 0258 /// collection of hits that are split per logical layer 0259 explicit SeedSolution(const TwoCircleTangentPars& pars, 0260 const Splitter_t& layerSorter) 0261 : TwoCircleTangentPars{pars}, m_splitter{layerSorter} {} 0262 0263 /// @brief Abrivation of the underlying space point reference 0264 using SpacePoint_t = ConstDeRef_t<typename UnCalibCont_t::value_type>; 0265 /// @brief Helper function to calculate the straw signs of the seed hits 0266 /// cached by this solution w.r.t. an external line 0267 /// @param seedPos: Reference point of the segment line 0268 /// @param seedDir: Direction of the segment line 0269 std::vector<int> leftRightAmbiguity(const Vector& seedPos, 0270 const Vector3& seedDir) const; 0271 0272 /// @brief Returns the number of hits cached by the seed 0273 std::size_t size() const { return m_seedHits.size(); } 0274 /// @brief Returns the i-th seed hit 0275 /// @param idx: Index of the hit to to return 0276 SpacePoint_t getHit(const std::size_t idx) const; 0277 /// @brief Appends a new seed hit to the solution 0278 void append(const std::size_t layIdx, const std::size_t hitIdx); 0279 /// @brief Vector of the associate left-rignt ambiguities 0280 std::vector<int> solutionSigns{}; 0281 ///@brief Number of good straw measurements 0282 std::size_t nStrawHits{0ul}; 0283 0284 private: 0285 /// @brief Pointer to the space point per layer splitter to gain access to the 0286 /// input space point container 0287 const Splitter_t& m_splitter; 0288 /// @brief Set of hits collected onto the seed. For each element 0289 /// the first index represents the layer & 0290 /// the second one the particular hit in that layer 0291 std::vector<std::pair<std::size_t, std::size_t>> m_seedHits{}; 0292 /// @brief Prints the seed solution to the screen 0293 /// @param ostr: Mutable reference to the stream to print to 0294 void print(std::ostream& ostr) const final; 0295 }; 0296 0297 public: 0298 /// @brief Helper struct to pack the parameters and the associated 0299 /// measurements into a common object. Returned by the 0300 /// central nextSeed method (cf. below) 0301 template <CompositeSpacePointContainer contType_t> 0302 struct SegmentSeed { 0303 /// @brief Constructor taking the seed parameters && 0304 /// a new hit container 0305 /// @param _pars: The seed line parameter 0306 /// @param _hits A new empty container to be filled 0307 explicit SegmentSeed(SeedParam_t _pars, contType_t&& _hits) noexcept 0308 : parameters{_pars}, hits{std::move(_hits)} {} 0309 /// @brief Seed line parameters 0310 SeedParam_t parameters; 0311 /// @brief Collection of hits 0312 contType_t hits; 0313 }; 0314 0315 /// @brief Central auxiliary struct to steer the seeding process. 0316 /// First, the user needs to implement the experiment specific 0317 /// CompSpacePointSeederDelegate over which the template of the 0318 /// SeedingState is then specified. The base class provides the 0319 /// straw hits split per logical layer and the SeedingState holds 0320 /// the indices to select iteratively a straw measurement each 0321 /// from the first and last layer. Also it keeps track of the 0322 /// previously constructed seeds. 0323 template <CompositeSpacePointContainer UncalibCont_t, 0324 CompositeSpacePointContainer CalibCont_t, 0325 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0326 Delegate_t> 0327 struct SeedingState : public Delegate_t { 0328 /// @brief Use all constructors from the base class for this one as well 0329 using Delegate_t::Delegate_t; 0330 /// @brief radius of the straw tubes used to reject hits outside the tube 0331 double strawRadius{0.}; 0332 /// @brief Try at the first time the external seed parameters as candidate 0333 bool startWithPattern{false}; 0334 /// @brief Estimated parameters from pattern 0335 SeedParam_t patternParams{}; 0336 /// @brief Stringstream output operator 0337 friend std::ostream& operator<<(std::ostream& ostr, 0338 const SeedingState& opts) { 0339 opts.print(ostr); 0340 return ostr; 0341 } 0342 /// @brief Return the number of generated seeds 0343 std::size_t nGenSeeds() const { return m_seenSolutions.size(); } 0344 /// @brief Grant the embedding class access to the private members 0345 friend CompositeSpacePointLineSeeder; 0346 0347 private: 0348 /// @brief List of straw measurement already constructed straw measurement seeds 0349 std::vector<SeedSolution<UncalibCont_t, Delegate_t>> m_seenSolutions{}; 0350 /// @brief @brief Index of the upper layer under consideration for the seeding 0351 std::optional<std::size_t> m_upperLayer{std::nullopt}; 0352 /// @brief Index of the lower layer under consideration for the seeding 0353 std::optional<std::size_t> m_lowerLayer{std::nullopt}; 0354 /// @brief Index of the hit in the lower layer under consideration for the seeding 0355 std::size_t m_lowerHitIndex{0ul}; 0356 /// @brief Index of the hit in the upper layer under consideration for the seeding 0357 std::size_t m_upperHitIndex{0ul}; 0358 /// @brief Index of the sign combination under consideration for the seeding 0359 std::size_t m_signComboIndex{0ul}; 0360 /// @brief Number of minimum straw hits a seed must have 0361 std::size_t m_nStrawCut{0ul}; 0362 /// @brief Flag toggling whether the upper of the lower layer shall be moved 0363 bool m_moveUpLayer{true}; 0364 /// @brief Prints the seed solution to the screen 0365 void print(std::ostream& ostr) const; 0366 }; 0367 /// @brief Main interface method provided by the SeederClass. The user instantiates 0368 /// a SeedingState object containing all the straw hit candidates from 0369 /// which the seed shall be constructed. Then, the nextSeed() returns 0370 /// the next best seed candidate which can then be fitted. The user 0371 /// continues to call the method until a nullopt is returned. 0372 /// @param cctx: Experiment specific calibration context to be piped back to the 0373 /// caller such that the space points may be calibrated during 0374 /// the seeding process. 0375 /// @param state: Mutable reference to the SeedingState object from which all the 0376 /// segment seeds are constructed. 0377 template <CompositeSpacePointContainer UncalibCont_t, 0378 CompositeSpacePointContainer CalibCont_t, 0379 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0380 Delegate_t> 0381 std::optional<SegmentSeed<CalibCont_t>> nextSeed( 0382 const CalibrationContext& cctx, 0383 SeedingState<UncalibCont_t, CalibCont_t, Delegate_t>& state) const; 0384 0385 private: 0386 /// @brief Reference to the logger object 0387 const Logger& logger() const { return *m_logger; } 0388 /// @brief Abrivation of the selector delegate to skip invalid straw hits in the seed 0389 template <CompositeSpacePointContainer Cont_t> 0390 using Selector_t = Delegate<bool(ConstDeRef_t<typename Cont_t::value_type>)>; 0391 /// @brief Abrivation of the split hit containers 0392 template <CompositeSpacePointContainer Cont_t> 0393 using StrawLayers_t = std::vector<Cont_t>; 0394 /// @brief Counts the number of hits inside the container. Depending on whether 0395 /// the busyLimitCountGood flag is true, bad hits are not considered 0396 /// @param container: Reference to the container which size is to be evaluated 0397 /// @param selector: Delegate method to skip bad bad hits 0398 template <CompositeSpacePointContainer Cont_t> 0399 std::size_t countHits(const Cont_t& container, 0400 const Selector_t<Cont_t>& selector) const; 0401 /// @brief Moves to the hit index to the next good hit inside the layer. 0402 /// The index is incremented until the underlying hit is accepted 0403 /// by the selector or all hits in the container were tried 0404 /// @param hitVec: Reference to the straw hits inside the layer 0405 /// @param selector: Delegate method to skip bad bad hits 0406 /// @param hitIdx: Mutable reference to the index that incremented 0407 template <CompositeSpacePointContainer UnCalibCont_t> 0408 bool moveToNextHit(const UnCalibCont_t& hitVec, 0409 const Selector_t<UnCalibCont_t>& selector, 0410 std::size_t& hitIdx) const; 0411 /// @brief Sets the parsed index to the first good hit inside the straw layer. 0412 /// @param hitVec: Reference to the straw hits inside the layer 0413 /// @param selector: Delegate method to skip bad bad hits 0414 /// @param hitIdx: Mutable reference to the index that incremented 0415 template <CompositeSpacePointContainer UnCalibCont_t> 0416 bool firstGoodHit(const UnCalibCont_t& hitVec, 0417 const Selector_t<UnCalibCont_t>& selector, 0418 std::size_t& hitIdx) const; 0419 /// @brief Move the layer index towards the possible value or,if the 0420 /// layer index is not yet initializes the lyaer index to 0421 // the next possible value. 0422 /// @param strawLayers: List of all straw hits split into the particular layers 0423 /// @param selector: Delegate method to skip bad hits 0424 /// @param boundary: Boundary value that the layer index must not cross 0425 /// @param layerIndex: Mutable reference to the layer index that needs to be moved 0426 /// @param hitIdx: Mutable reference to the associated hit index inside the layer 0427 /// @param moveForward: Flag toggling whether the layer index shall be incremented or 0428 /// decremented. 0429 template <CompositeSpacePointContainer UnCalibCont_t> 0430 bool nextLayer(const StrawLayers_t<UnCalibCont_t>& strawLayers, 0431 const Selector_t<UnCalibCont_t>& selector, 0432 const std::size_t boundary, 0433 std::optional<std::size_t>& layerIndex, std::size_t& hitIdx, 0434 bool moveForward) const; 0435 /// @brief Move the layer and hit indices inside the state towards the next candidate 0436 /// pair. First, the L-R ambiguities are incremented, then it is 0437 /// searched for the next pair inside the lower && upper layer pair. 0438 /// Finally, the indices are moved towards the next layer 0439 /// @param selector: Delegate method to skip bad hits 0440 /// @param state: Mutable reference to the SeedingState object carring the state indices 0441 template <CompositeSpacePointContainer UncalibCont_t, 0442 CompositeSpacePointContainer CalibCont_t, 0443 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0444 Delegate_t> 0445 void moveToNextCandidate( 0446 const Selector_t<UncalibCont_t>& selector, 0447 SeedingState<UncalibCont_t, CalibCont_t, Delegate_t>& state) const; 0448 /// @brief Attempts to construct the next seed from the given configuration of 0449 /// seed circles. The seed needs to contain a minimum number of other 0450 /// straw hits and there must be no other previously constructed seed 0451 /// with the same Left-Right solution 0452 /// @param cctx: Calibration context to be piped to the experiment's implementation 0453 /// such that conditions data access becomes possible 0454 /// @param selector: Delegate method to skip bad hits 0455 /// @param state: Mutable reference to the SeedingState object carring the state indices 0456 template <CompositeSpacePointContainer UncalibCont_t, 0457 CompositeSpacePointContainer CalibCont_t, 0458 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0459 Delegate_t> 0460 std::optional<SegmentSeed<CalibCont_t>> buildSeed( 0461 const CalibrationContext& cctx, const Selector_t<UncalibCont_t>& selector, 0462 SeedingState<UncalibCont_t, CalibCont_t, Delegate_t>& state) const; 0463 /// @brief Checks whether the new seed candidate passes the quality cuts on 0464 /// the number of good straw hits and whether it is not within the 0465 /// same overlap corridor as previously produced seeds 0466 /// @param tangentSeed: Pair of reference position & direction constructed 0467 /// from the two line tangent seed 0468 /// @param newSolution: The new seed solution that's to be tested 0469 /// @param state: The cache carrying the already produced solutions 0470 template <CompositeSpacePointContainer UncalibCont_t, 0471 CompositeSpacePointContainer CalibCont_t, 0472 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0473 Delegate_t> 0474 bool passSeedCuts( 0475 const Line_t& tangentSeed, 0476 SeedSolution<UncalibCont_t, Delegate_t>& newSolution, 0477 SeedingState<UncalibCont_t, CalibCont_t, Delegate_t>& state) const; 0478 /// @brief Converts the accepted seed solution to the segment seed returned by 0479 /// nextSeed and adds the strip measurements to the seed. The solution 0480 /// is then appended to the state 0481 /// @param cctx: Calibration context to be piped to the experiment's implementation 0482 /// such that conditions data access becomes possible 0483 /// @param tangentSeed: Position and direction constructed from the current tangent seed 0484 /// @param state: Mutable reference to the state from which the strip measurements are drawn 0485 /// and to which the newSolution is then appended 0486 /// @param newSolution: Current tangent seed solution object holding the straw measurements 0487 /// to be put onto the seed. 0488 template <CompositeSpacePointContainer UncalibCont_t, 0489 CompositeSpacePointContainer CalibCont_t, 0490 detail::CompSpacePointSeederDelegate<UncalibCont_t, CalibCont_t> 0491 Delegate_t> 0492 SegmentSeed<CalibCont_t> consructSegmentSeed( 0493 const CalibrationContext& cctx, const Line_t& tangentSeed, 0494 SeedingState<UncalibCont_t, CalibCont_t, Delegate_t>& state, 0495 SeedSolution<UncalibCont_t, Delegate_t>&& newSolution) const; 0496 /// @brief Construct the final seed parameters by combining the initial 0497 /// pattern parameters with the parameter from two circle tangent 0498 /// @param tangentSeed: Pair of reference position & direction constructed 0499 /// from the two line tangent seed 0500 /// @param patternParams: Parameter estimate from the hit pattern 0501 SeedParam_t combineWithPattern(const Line_t& tangentSeed, 0502 const SeedParam_t& patternParams) const; 0503 /// @brief Constructs a line from the parsed seed parameters. The 0504 /// first element is the reference point && the second one 0505 /// is the direction 0506 /// @param pars: Reference to the line parameters from which the line is created 0507 Line_t makeLine(const SeedParam_t& pars) const; 0508 /// @brief Check whether the generated seed parameters are within the ranges defined by the used 0509 /// @param tangentPars: Reference to the seed parameters to check 0510 bool isValidLine(const TwoCircleTangentPars& tangentPars) const; 0511 /// @brief Configuration object 0512 Config m_cfg{}; 0513 /// @brief Logger instance 0514 std::unique_ptr<const Logger> m_logger{}; 0515 /// @brief Array encoding the four possible left right solutions. 0516 /// The first (second) index encodes the ambiguity of the 0517 /// bottom (top) straw measurement. The order of the index 0518 /// pairs is coupled to the TangentAmbi index 0519 static constexpr std::array<std::array<int, 2>, 4> s_signCombo{ 0520 std::array{1, 1}, std::array{1, -1}, std::array{-1, 1}, 0521 std::array{-1, -1}}; 0522 }; 0523 } // namespace Acts::Experimental 0524 #include "Acts/Seeding/CompositeSpacePointLineSeeder.ipp"
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|