ATLAS Offline Software
Loading...
Searching...
No Matches
SegmentLineFitter.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Compile this file assuming that FP operations may trap.
6// Prevents spurious FPEs in the clang build.
9
12
15
16#include <ActsInterop/Logger.h>
19
20
21#include <format>
22
25
26
27namespace MuonR4::SegmentFit{
28 using namespace Acts;
29 using namespace Acts::UnitLiterals;
30
34
35 namespace {
38 bool isGoodHit(const MuonR4::CalibratedSpacePoint& hit) {
40 return hit.fitState() == Valid;
41 }
43 bool isPrecisionHit(const MuonR4::CalibratedSpacePoint& hit) {
44 using enum xAOD::UncalibMeasType;
45 return isGoodHit(hit) && (
47 hit.type() == MdtDriftCircleType || hit.type() == MMClusterType ||
49 (hit.type() == sTgcStripType &&
50 static_cast<const xAOD::sTgcMeasurement*>(hit.spacePoint()->primaryMeasurement())->channelType() ==
52 );
53 }
56 inline unsigned countPrecHits(const HitVec_t& hits) {
57 return std::ranges::count_if(hits, [](const Hit_t& hit){
58 return isPrecisionHit(*hit);
59 });
60 }
62 inline unsigned countPhiHits(const HitVec_t& hits) {
63 return std::ranges::count_if(hits, [](const Hit_t& hit){
64 return isGoodHit(*hit) && hit->measuresPhi();
65 });
66 }
70 inline void removeBeamSpot(HitVec_t& hits){
71 hits.erase(std::remove_if(hits.begin(), hits.end(),
72 [](const Hit_t& a){
73 return a->type() == xAOD::UncalibMeasType::Other;
74 }), hits.end());
75 }
77 HitVec_t copy(const HitVec_t& hits) {
78 HitVec_t copied{};
79 copied.reserve(hits.size());
80 std::ranges::transform(hits, std::back_inserter(copied),
81 [](const Hit_t& hit) {
82 return std::make_unique<CalibratedSpacePoint>(*hit);
83 });
84 return copied;
85 }
86 Result_t copy(const Result_t& toCopy) {
87 Result_t toRet{};
88 using FitPars_t = SegmentLineFitter::FitPars_t;
89 static_cast<FitPars_t&>(toRet) = toCopy;
90 toRet.measurements = copy(toCopy.measurements);
91 return toRet;
92 }
93 }
94 SegmentLineFitter::Config::RangeArray
96 RangeArray rng{};
97 constexpr double spatRang = 10._m;
98 constexpr double timeTange = 50._ns;
99 using enum ParamDefs;
100 rng[toUnderlying(y0)] = std::array{-spatRang, spatRang};
101 rng[toUnderlying(x0)] = std::array{-spatRang, spatRang};
102 rng[toUnderlying(phi)] = std::array{-179._degree, 179._degree};
103 rng[toUnderlying(theta)] = std::array{0._degree, 175._degree};
104 rng[toUnderlying(t0)] = std::array{-timeTange, timeTange};
105 return rng;
106 }
108 AthMessaging{name},
110 m_cfg{config} {
111 m_goodHitSel.connect<isGoodHit>();
112 }
113 Result_t SegmentLineFitter::callLineFit(const Acts::CalibrationContext& cctx,
114 const Parameters& startPars,
115 const Amg::Transform3D& localToGlobal,
116 HitVec_t&& calibHits) const {
118 if (m_cfg.doBeamSpot && countPhiHits(calibHits) > 0) {
119 const Amg::Transform3D globToLoc{localToGlobal.inverse()};
120 Amg::Vector3D beamSpot{globToLoc.translation()};
121 Amg::Vector3D beamLine{globToLoc.linear().col(2)};
122 SpacePoint::Cov_t covariance{};
123 covariance[toUnderlying(AxisDefs::etaCov)] = square(m_cfg.beamSpotRadius);
124 covariance[toUnderlying(AxisDefs::phiCov)] = square(m_cfg.beamSpotLength);
126 auto beamSpotSP = std::make_unique<CalibratedSpacePoint>(nullptr, std::move(beamSpot));
127 beamSpotSP->setBeamDirection(std::move(beamLine));
128 beamSpotSP->setCovariance(std::move(covariance));
129 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__<<": Beam spot constraint "
130 <<Amg::toString(beamSpotSP->localPosition())<<", "<<beamSpotSP->covariance());
131 calibHits.push_back(std::move(beamSpotSP));
132 }
133
134 if (msgLvl(MSG::VERBOSE)) {
135 const auto [pos, dir] = makeLine(startPars);
136 std::stringstream hitStream{};
137 for (const Hit_t& hit : calibHits) {
138 hitStream<<" **** "<< (*hit)<<", pull: "
139 <<std::sqrt(SeedingAux::chi2Term(pos, dir, *hit)) <<std::endl;
140 }
141 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": Start segment fit with parameters "
142 <<toString(startPars) <<", plane location: "<<Amg::toString(localToGlobal)<<std::endl
143 <<hitStream.str());
144 }
145 FitOpts_t fitOpts{};
147 fitOpts.calibContext = cctx;
148 fitOpts.calibrator = m_cfg.calibrator;
149 fitOpts.selector = m_goodHitSel;
150 //check the degrees of freedom before try the fit
151 const auto dOF = m_fitter.countDoF(calibHits, fitOpts.selector);
152 if((dOF.bending + dOF.nonBending) < m_fitter.config().parsToUse.size()){
153 return result;
154 }
155 fitOpts.measurements = std::move(calibHits);
156 fitOpts.localToGlobal = localToGlobal;
157 fitOpts.startParameters = startPars;
159 constexpr auto t0idx = toUnderlying(ParamDefs::t0);
160 fitOpts.startParameters[t0idx] = ActsTrk::timeToActs(fitOpts.startParameters[t0idx]);
162 result = m_fitter.fit(std::move(fitOpts));
164 if (m_fitter.config().fitT0){
165 result.parameters[t0idx] = ActsTrk::timeToAthena(result.parameters[t0idx]);
166 result.covariance(t0idx, t0idx) = Acts::square(ActsTrk::timeToAthena(1.)) * result.covariance(t0idx, t0idx);
167 for (ParamDefs p : {ParamDefs::x0, ParamDefs::y0, ParamDefs::phi, ParamDefs::theta}) {
168 auto pidx = toUnderlying(p);
169 result.covariance(t0idx, pidx) = ActsTrk::timeToAthena(result.covariance(t0idx, pidx));
170 result.covariance(pidx, t0idx) = ActsTrk::timeToAthena(result.covariance(pidx, t0idx));
171 }
172 }
174 return result;
175 }
176 std::unique_ptr<Segment>
177 SegmentLineFitter::fitSegment(const EventContext& ctx,
178 const SegmentSeed* parent,
179 const Parameters& startPars,
180 const Amg::Transform3D& localToGlobal,
181 HitVec_t&& calibHits) const {
182
183 const Acts::CalibrationContext cctx = ActsTrk::getCalibrationContext(ctx);
184 if (m_cfg.visionTool) {
185 Result_t preFit{};
186 preFit.parameters = startPars;
187 preFit.measurements = copy(calibHits);
188 auto seedCopy = convertToSegment(localToGlobal, parent, std::move(preFit));
189 m_cfg.visionTool->visualizeSegment(ctx, *seedCopy, "Prefit");
190 }
191 Result_t segFit = callLineFit(cctx, startPars, localToGlobal, std::move(calibHits));
192 if (m_cfg.visionTool && segFit.converged) {
193 auto seedCopy = convertToSegment(localToGlobal, parent, copy(segFit));
194 m_cfg.visionTool->visualizeSegment(ctx, *seedCopy, "Intermediate fit");
195 }
196 if (!removeOutliers(cctx, *parent, localToGlobal,
197 segFit.converged? segFit.parameters : startPars,
198 segFit)) {
199 return nullptr;
200 }
201 if (!plugHoles(cctx, *parent, localToGlobal, segFit)) {
202 return nullptr;
203 }
204 auto finalSeg = convertToSegment(localToGlobal, parent, std::move(segFit));
205 if (m_cfg.visionTool) {
206 m_cfg.visionTool->visualizeSegment(ctx, *finalSeg, "Final fit");
207 }
208 return finalSeg;
209 }
210 std::unique_ptr<Segment>
212 const SegmentSeed* patternSeed,
213 Result_t&& data) const {
214 const auto [locPos, locDir] = makeLine(data.parameters);
215 Amg::Vector3D globPos = locToGlob * locPos;
216 Amg::Vector3D globDir = locToGlob.linear()* locDir;
217
218 std::ranges::sort(data.measurements, [](const Hit_t& a, const Hit_t& b){
219 return a->localPosition().z() < b->localPosition().z();
220 });
221 if (msgLvl(MSG::VERBOSE)) {
222 std::stringstream sstr{};
223 for (const Hit_t& h : data.measurements) {
224 sstr<<" **** "<<(*h)<<", pull: "
225 <<std::sqrt(SeedingAux::chi2Term(locPos, locDir, *h))<<std::endl;
226 }
227 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": Create new segment "
228 <<toString(data.parameters)<<" in "<<patternSeed->msSector()->identString()<<"built from:\n"<<sstr.str());
229 }
230
231 auto finalSeg = std::make_unique<Segment>(std::move(globPos), std::move(globDir),
232 patternSeed, std::move(data.measurements),
233 data.chi2, data.nDoF);
234 finalSeg->setCallsToConverge(data.nIter);
235 finalSeg->setParUncertainties(std::move(data.covariance));
236 if (m_fitter.config().fitT0) {
237 finalSeg->setSegmentT0(data.parameters[toUnderlying(ParamDefs::t0)]);
238 }
239 return finalSeg;
240 }
241
242 bool SegmentLineFitter::removeOutliers(const Acts::CalibrationContext& cctx,
243 const SegmentSeed& seed,
244 const Amg::Transform3D& localToGlobal,
245 const LinePar_t& startPars,
246 Result_t& fitResult) const {
247
248
249 if (countPrecHits(fitResult.measurements) < m_cfg.nPrecHitCut || fitResult.nDoF == 0
250 || fitResult.nIter > m_fitter.config().maxIter) {
251 for(const auto& meas : fitResult.measurements){
252 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__<<": Measurement from fitresult is" << (*meas));
253 }
254 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__
255 <<": No degree of freedom available. What shall be removed?!. nDoF: "
256 <<fitResult.nDoF<<", n-meas: "<<countPrecHits(fitResult.measurements));
257 return false;
258 }
259 if (fitResult.converged && fitResult.chi2 / std::max(fitResult.nDoF, 1ul) < m_cfg.outlierRemovalCut) {
260 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": The segment "<<toString(fitResult.parameters)
261 <<" is already of good quality "<<fitResult.chi2 / std::max(fitResult.nDoF, 1ul)
262 <<". Don't remove outliers");
263 return true;
264 }
265 if (fitResult.nDoF == 0u){
266 return false;
267 }
268 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": Segment "
269 <<toString(fitResult.parameters)<<" is of badish quality. Remove worst hit");
270
273 if (m_cfg.doBeamSpot) {
274 removeBeamSpot(fitResult.measurements);
275 }
276 const auto [segPos, segDir] = makeLine(fitResult.parameters);
278 std::ranges::sort(fitResult.measurements,
279 [&segPos, &segDir](const HitVec_t::value_type& a, const HitVec_t::value_type& b){
280 const double chiSqA = isGoodHit(*a) ? SeedingAux::chi2Term(segPos, segDir, *a) : 0.;
281 const double chiSqB = isGoodHit(*b) ? SeedingAux::chi2Term(segPos, segDir, *b) : 0.;
282 return chiSqA < chiSqB;
283 });
284 fitResult.measurements.back()->setFitState(HitState::Outlier);
285
287 Result_t newAttempt = callLineFit(cctx, startPars, localToGlobal,
288 std::move(fitResult.measurements));
289 if (newAttempt.converged) {
290 newAttempt.nIter+=fitResult.nIter;
291 fitResult = std::move(newAttempt);
292 if (m_cfg.visionTool) {
293 const EventContext& ctx{*cctx.get<const EventContext*>()};
294 auto seedCopy = convertToSegment(localToGlobal, &seed, copy(fitResult));
295 m_cfg.visionTool->visualizeSegment(ctx, *seedCopy, "Bad fit recovery");
296 }
297 } else {
298 fitResult.nIter+=newAttempt.nIter;
299 fitResult.measurements = std::move(newAttempt.measurements);
300 }
301 return removeOutliers(cctx, seed, localToGlobal,
302 fitResult.converged ? fitResult.parameters : startPars,
303 fitResult);
304 }
305
307 auto [segPos, segDir] = makeLine(candidate.parameters);
308 cleanStripLayers(segPos, segDir, candidate.measurements);
309 candidate.measurements.erase(std::remove_if(candidate.measurements.begin(),
310 candidate.measurements.end(),
311 [&](const HitVec_t::value_type& hit){
312 if (hit->fitState() == HitState::Valid) {
313 return false;
314 } else if (hit->fitState() == HitState::Duplicate) {
315 return true;
316 }
319 const double dist = Amg::lineDistance(segPos, segDir,
320 hit->localPosition(),
321 hit->sensorDirection());
322 const auto* dc = static_cast<const xAOD::MdtDriftCircle*>(hit->spacePoint()->primaryMeasurement());
323 return dist >= dc->readoutElement()->innerTubeRadius();
324 }
325 return false;
326 }), candidate.measurements.end());
327 }
329 const Amg::Vector3D& lineDir,
330 HitVec_t& hits) const {
331 const SpacePointPerLayerSorter sorter{};
333 std::ranges::sort(hits, [&sorter, &linePos, &lineDir](const Hit_t&a ,const Hit_t& b){
334 if (a->isStraw() || b->isStraw()) {
335 return !a->isStraw();
336 }
337 if (a->type() == xAOD::UncalibMeasType::Other ||
338 b->type() == xAOD::UncalibMeasType::Other) {
339 return a->type() != xAOD::UncalibMeasType::Other;
340 }
341 const unsigned lay_a = sorter.sectorLayerNum(*a->spacePoint());
342 const unsigned lay_b = sorter.sectorLayerNum(*b->spacePoint());
343 if (lay_a != lay_b) {
344 return lay_a < lay_b;
345 }
346 return SeedingAux::chi2Term(linePos, lineDir, *a) <
347 SeedingAux::chi2Term(linePos, lineDir, *b);
348 });
350 for (HitVec_t::iterator itr = hits.begin(); itr != hits.end(); ++itr) {
351 const Hit_t& hit_a{*itr};
352 if (hit_a->isStraw()){
353 break;
354 }
355 if(hit_a->fitState() == HitState::Duplicate ||
356 hit_a->type() == xAOD::UncalibMeasType::Other) {
357 continue;
358 }
359 const unsigned lay_a = sorter.sectorLayerNum(*hit_a->spacePoint());
361 for (HitVec_t::iterator itr2 = itr + 1; itr2 != hits.end(); ++itr2) {
362 const Hit_t& hit_b{*itr2};
363 if (hit_b->type() == xAOD::UncalibMeasType::Other) {
364 continue;
365 }
366 if (lay_a != sorter.sectorLayerNum(*hit_b->spacePoint())) {
367 break;
368 }
370 if ( (hit_a->measuresEta() && hit_b->measuresEta()) ||
371 (hit_a->measuresPhi() && hit_b->measuresPhi())) {
372 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": Reject "
373 <<m_cfg.idHelperSvc->toString(hit_b->spacePoint()->identify()) <<" in favour of "
374 <<m_cfg.idHelperSvc->toString(hit_a->spacePoint()->identify()));
375 hit_b->setFitState(HitState::Duplicate);
376 }
377 }
378 }
379 }
380 bool SegmentLineFitter::plugHoles(const Acts::CalibrationContext& cctx,
381 const SegmentSeed& seed,
382 const Amg::Transform3D& localToGlobal,
383 Result_t& toRecover) const {
385 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": segment "<<toString(toRecover.parameters)
386 <<", chi2: "<<toRecover.chi2 /std::max(toRecover.nDoF, 1ul)
387 <<", nDoF: "<<toRecover.nDoF);
389
390
391 std::unordered_set<const SpacePoint*> usedSpacePoints{};
392 for (auto& hit : toRecover.measurements) {
393 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__ <<": "<<(*hit)<<" is known");
394 usedSpacePoints.insert(hit->spacePoint());
395 }
397 const EventContext& ctx{*cctx.get<const EventContext*>()};
398
399 const double timeOff = toRecover.parameters[toUnderlying(ParamDefs::t0)];
400 HitVec_t candidateHits{};
401 bool hasCandidate{false};
402 const auto [locPos, locDir] = makeLine(toRecover.parameters);
403
405 for (const auto& hit : *seed.parentBucket()){
407 if (usedSpacePoints.count(hit.get())){
408 continue;
409 }
410 Hit_t calibHit{};
411 double pull{-1.};
412 if (hit->isStraw()) {
413 using namespace Acts::detail::LineHelper;
414 const double dist = signedDistance(locPos, locDir, hit->localPosition(), hit->sensorDirection());
415 const auto* dc = static_cast<const xAOD::MdtDriftCircle*>(hit->primaryMeasurement());
416 // Check whether the tube is crossed by the hit
417 if (Acts::abs(dist) >= dc->readoutElement()->innerTubeRadius()) {
418 continue;
419 }
420 } else {
422 if (!hit->measuresEta() &&
423 std::abs(hit->sensorDirection().dot(hit->localPosition() -
424 SeedingAux::extrapolateToPlane(locPos,locDir, *hit))) >
425 std::sqrt(hit->covariance()[toUnderlying(AxisDefs::etaCov)])){
426 continue;
427 }
430 pull = std::sqrt(SeedingAux::chi2Term(locPos, locDir, *hit));
431 if (pull > 1.1 * m_cfg.recoveryPull) {
432 continue;
433 }
434 }
435 calibHit = m_cfg.calibrator->calibrate(ctx, hit.get(), locPos, locDir, ActsTrk::timeToActs(timeOff));
436 pull = std::sqrt(SeedingAux::chi2Term(locPos, locDir, *calibHit));
437 if (pull <= m_cfg.recoveryPull) {
438 hasCandidate |= calibHit->fitState() == HitState::Valid;
439 } else {
440 calibHit->setFitState(HitState::Outlier);
441 }
442 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__<<": Candidate hit for recovery "
443 <<seed.msSector()->idHelperSvc()->toString(hit->identify())<<", pull: "<<pull);
444 candidateHits.push_back(std::move(calibHit));
445 }
447 if (!hasCandidate) {
448 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__<<": No space point candidates for recovery were found");
449 toRecover.measurements.insert(toRecover.measurements.end(),
450 std::make_move_iterator(candidateHits.begin()),
451 std::make_move_iterator(candidateHits.end()));
452 eraseWrongHits(toRecover);
453 return toRecover.nDoF > 0;
454 }
455
456
457 HitVec_t copied = copy(toRecover.measurements);
458 HitVec_t copiedCandidates = copy(candidateHits);
460 if (m_cfg.doBeamSpot) {
461 removeBeamSpot(copied);
462 }
463
464 candidateHits.insert(candidateHits.end(),
465 std::make_move_iterator(copied.begin()),
466 std::make_move_iterator(copied.end()));
467
468 cleanStripLayers(locPos, locDir, candidateHits);
469 Result_t recovered = callLineFit(cctx, toRecover.parameters, localToGlobal,
470 std::move(candidateHits));
471 if (!recovered.converged) {
472 return false;
473 }
475 if (recovered.nDoF <= toRecover.nDoF) {
476 for (HitVec_t::value_type& hit : copiedCandidates) {
477 hit->setFitState(HitState::Outlier);
478 toRecover.measurements.push_back(std::move(hit));
479 }
480 eraseWrongHits(toRecover);
481 return true;
482 }
483 std::vector<const CalibratedSpacePoint*> stripOutliers{};
484 stripOutliers.reserve(toRecover.measurements.size());
485 ATH_MSG_VERBOSE(__func__<<"() - "<<__LINE__<<": Before - chi2: "<<toRecover.chi2
486 <<", nDoF "<<toRecover.nDoF<<" <=> after recovery - chi2: "
487 <<recovered.chi2<<", nDoF: "<<recovered.nDoF);
488
489 double redChi2 = recovered.chi2 / std::max(recovered.nDoF, 1ul);
492 if (redChi2 < m_cfg.outlierRemovalCut ||
493 toRecover.nDoF == 0 || redChi2 < toRecover.chi2 / toRecover.nDoF) {
494 ATH_MSG_VERBOSE("Accept segment with recovered "<<(recovered.nDoF - toRecover.nDoF)<<" hits.");
495 recovered.nIter += toRecover.nIter;
496 toRecover = std::move(recovered);
499 unsigned recovLoop{0u};
500 while (++recovLoop <= m_cfg.nRecoveryLoops) {
501 copied = copy(toRecover.measurements);
502 if (m_cfg.doBeamSpot) {
503 removeBeamSpot(copied);
504 }
505 const auto [beforePos, beforeDir] = makeLine(toRecover.parameters);
506
507 for (HitVec_t::value_type& copyHit : copied) {
508 if (copyHit->fitState() != HitState::Outlier) {
509 continue;
510 }
511 if (std::sqrt(SeedingAux::chi2Term(beforePos, beforeDir, *copyHit)) < m_cfg.recoveryPull) {
512 copyHit->setFitState(HitState::Valid);
513 stripOutliers.push_back(copyHit.get());
514 }
515 }
516 // Nothing to recover
517 if (stripOutliers.empty()) {
518 break;
519 }
520 cleanStripLayers(beforePos, beforeDir, copied);
521 // Recovery turned out to be duplicates on the same layer
522 if (std::ranges::none_of(stripOutliers,[](const CalibratedSpacePoint* sp){
523 return sp->fitState() == HitState::Valid;
524 })) {
525 break;
526 }
527 stripOutliers.clear();
528 recovered = callLineFit(cctx, toRecover.parameters, localToGlobal, std::move(copied));
529 if (!recovered.converged) {
530 break;
531 }
532 if (recovered.nDoF <= toRecover.nDoF) {
533 break;
534 }
535 redChi2 = recovered.chi2 / std::max(recovered.nDoF, 1ul);
536 if (redChi2 < m_cfg.outlierRemovalCut || redChi2 < toRecover.chi2 / toRecover.nDoF) {
537 recovered.nIter += toRecover.nIter;
538 toRecover = std::move(recovered);
539 } else {
540 break;
541 }
542 }
543 } else{
544 for (HitVec_t::value_type& hit : copiedCandidates) {
545 hit->setFitState(HitState::Outlier);
546 toRecover.measurements.push_back(std::move(hit));
547 }
548 }
549 eraseWrongHits(toRecover);
550 return true;
551 }
553 if (std::ranges::any_of(result.measurements,[](const Hit_t& h){
554 return h->measuresPhi();
555 })) {
556 ATH_MSG_VERBOSE("The segment has phi measurements. No centering needed");
557 return;
558 }
559 double avgX{0.};
560 const double nHits = result.measurements.size();
561 std::ranges::for_each(result.measurements, [&avgX, nHits](const Hit_t& hit) {
562 return avgX += hit->localPosition().x() / nHits;
563 });
564 result.parameters[toUnderlying(ParamDefs::x0)] = avgX;
565 }
566
567
568}
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_MSG_VERBOSE(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
static Double_t sp
static Double_t a
static Double_t t0
if(febId1==febId2)
static const uint32_t nHits
std::unique_ptr< const Acts::Logger > makeActsAthenaLogger(IMessageSvc *svc, const std::string &name, int level, std::optional< std::string > parent_name)
Header file for AthHistogramAlgorithm.
bool msgLvl(const MSG::Level lvl) const
Test the output level.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
std::string identString() const
Returns a string encoding the chamber index & the sector of the MS sector.
The calibrated Space point is created during the calibration process.
const SpacePoint * spacePoint() const
The pointer to the space point out of which this space point has been built.
xAOD::UncalibMeasType type() const
Returns the space point type.
State
State flag to distinguish different space point states.
State fitState() const
Returns the state of the calibrated space point.
void centerAlongWire(Result_t &fitResult) const
Moves the segment to the average x0 position, if the segment does not contain any measurement.
bool plugHoles(const Acts::CalibrationContext &cctx, const SegmentSeed &seed, const Amg::Transform3D &localToGlobal, Result_t &toRecover) const
Recovery of missed hits.
std::unique_ptr< CalibratedSpacePoint > Hit_t
Abrivation of the space point type to use.
void cleanStripLayers(const Amg::Vector3D &linePos, const Amg::Vector3D &lineDir, HitVec_t &hits) const
Marks duplicate hits on a strip layer as outliers to avoid competing contributions from the same laye...
Selector_t m_goodHitSel
Selector to identify the valid hits.
Result_t callLineFit(const Acts::CalibrationContext &cctx, const Parameters &startPars, const Amg::Transform3D &localToGlobal, HitVec_t &&calibHits) const
Calls the underlying line fitter to determine the segment parameters.
ConfigSwitches m_cfg
Configuration switches of the ATLAS fitter implementation.
Fitter_t::FitResult< HitVec_t > Result_t
Abrivation of the fit result.
void eraseWrongHits(Result_t &candidate) const
Removes all hits from the segment which are obvious outliers.
std::unique_ptr< Segment > convertToSegment(const Amg::Transform3D &locToGlobTrf, const SegmentSeed *parentSeed, Result_t &&toConvert) const
Converts the fit result into a segment object.
Fitter_t::FitOptions< HitVec_t, ISpacePointCalibrator > FitOpts_t
Abrivation of the fit options.
Fitter_t::FitParameters FitPars_t
Abrivation of the fit parameters.
std::unique_ptr< Segment > fitSegment(const EventContext &ctx, const SegmentSeed *parent, const LinePar_t &startPars, const Amg::Transform3D &localToGlobal, HitVec_t &&calibHits) const
Fit a set of measurements to a straight segment line.
Fitter_t m_fitter
Actual implementation of the straight line fit.
Fitter_t::ParamVec_t LinePar_t
Abrivation of the fitted line parameters.
std::vector< Hit_t > HitVec_t
Collection of space points.
SegmentLineFitter(const std::string &name, Config &&config)
Standard constructor.
bool removeOutliers(const Acts::CalibrationContext &cctx, const SegmentSeed &seed, const Amg::Transform3D &localToGlobal, const LinePar_t &startPars, Result_t &fitResult) const
Cleans the fitted segment from the most outlier hit and then attempts to refit the segment.
Representation of a segment seed (a fully processed hough maximum) produced by the hough transform.
Definition SegmentSeed.h:14
const MuonGMR4::SpectrometerSector * msSector() const
Returns the associated chamber.
The SpacePointPerLayerSorter sort two given space points by their layer Identifier.
std::array< double, 3 > Cov_t
Abrivation of the covariance type.
const xAOD::UncalibratedMeasurement * primaryMeasurement() const
constexpr double timeToAthena(const double actsT)
Converts a time unit from Acts to Athena units.
Acts::CalibrationContext getCalibrationContext(const EventContext &ctx)
The Acts::Calibration context is piped through the Acts fitters to (re)calibrate the Acts::SourceLink...
constexpr double timeToActs(const double athenaT)
Converts a time unit from Athena to Acts units.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
SegmentLineFitter::Result_t Result_t
SegmentLineFitter::HitVec_t HitVec_t
SeedingAux::FitParIndex ParamDefs
Use the same parameter indices as used by the CompSpacePointAuxiliaries.
SegmentLineFitter::Hit_t Hit_t
std::pair< Amg::Vector3D, Amg::Vector3D > makeLine(const Parameters &pars)
Returns the parsed parameters into an Eigen line parametrization.
Acts::Experimental::CompositeSpacePointLineFitter::ParamVec_t Parameters
std::string toString(const Parameters &pars)
Dumps the parameters into a string with labels in front of each number.
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.
MdtDriftCircle_v1 MdtDriftCircle
UncalibMeasType
Define the type of the uncalibrated measurement.
sTgcMeasurement_v1 sTgcMeasurement
static RangeArray defaultRanges()
Function that returns a set of predefined ranges for testing.
Tell the compiler to optimize assuming that FP may trap.
#define CXXUTILS_TRAPPING_FP
Definition trapping_fp.h:24