ATLAS Offline Software
Loading...
Searching...
No Matches
MsTrackSeeder.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
8
9#include "Acts/Utilities/Helpers.hpp"
10#include "Acts/Surfaces/detail/LineHelper.hpp"
11#include "Acts/Definitions/Units.hpp"
12#include "GaudiKernel/PhysicalConstants.h"
13namespace {
14 using namespace Acts::UnitLiterals;
15 constexpr double inDeg(const double a) {
16 return a / 1._degree;
17 }
19 void average(const Amg::Vector3D& segPos, std::optional<Amg::Vector3D>& posSlot) {
20 if (!posSlot) {
21 posSlot = segPos;
22 } else {
23 *posSlot = 0.5 * (*posSlot + segPos);
24 }
25 }
27 inline Amg::Vector3D dirForBField(const xAOD::MuonSegment& seg,
28 const double circPhi) {
29 //if (seg.nPhiLayers() > 0) {
30 // return seg.direction();
31 //}
32 const double theta = std::atan2(Acts::fastHypot(seg.px(), seg.py()), seg.pz());
33 return Acts::makeDirectionFromPhiTheta(circPhi, theta);
34 }
35 std::string print(const Amg::Vector3D& v){
36 return std::format("r: {:.2f}, z: {:.2f}, phi: {:.2f}, theta: {:.2f}",
37 v.perp(), v.z(), inDeg((v.phi())), inDeg(v.theta()));
38 }
39
40 float reducedChi2(const xAOD::MuonSegment& seg) {
41 return seg.chiSquared() / std::max(1.f, seg.numberDoF());
42 }
43 std::string print(const xAOD::MuonSegment& seg) {
44 return std::format("{:}, nPrecHits: {:}, nPhiHits: {:}", MuonR4::printID(seg),
45 seg.nPrecisionHits(), seg.nPhiLayers());
46 }
47 static const Muon::MuonSectorMapping sectorMap{};
48}
49
50namespace MuonR4{
54 using enum SectorProjector;
55 switch (proj) {
56 case leftOverlap:
57 return "overlap with left sector";
58 case center:
59 return "sector center";
60 case rightOverlap:
61 return "overlap with right sector";
62 default:
63 return "";
64 }
65 }
66 int MsTrackSeeder::ringSector(const int sector) {
67 constexpr int nSectors = Muon::MuonStationIndex::numberOfSectors();
68 return sector == 0 ? nSectors : (sector > nSectors ? 1 : sector);
69 }
70 int MsTrackSeeder::ringOverlap(const int sector) {
71 constexpr int nSectors = 2*Muon::MuonStationIndex::numberOfSectors();
72 return sector > nSectors ? 1 : sector;
73 }
74
75 MsTrackSeeder::MsTrackSeeder(const std::string& msgName, Config&& cfg):
76 AthMessaging{msgName},
77 m_cfg{std::move(cfg)}{
78 auto& v{m_cfg.fieldExtpSteps};
79 v.insert(0.); v.insert(1.);
80 if (std::ranges::any_of(v, [](const double x){
81 return x < 0. || x > 1.;
82 })) {
83 THROW_EXCEPTION("Found invalid extrapolation steps.");
84 }
85 }
86
87
88 double MsTrackSeeder::projectedPhi(const int sector,
89 const SectorProjector proj) {
90 return sectorMap.sectorOverlapPhi(sector, ringSector(sector + Acts::toUnderlying(proj)));
91 }
94 return m_cfg.detMgr->getSectorEnvelope(segment.chamberIndex(),
95 segment.sector(),
96 segment.etaIndex());
97 }
100 const MsTrackSeed& refSeed) {
101
102 int doubSector = 2 * seg.sector();
103 constexpr int nSectors = 2*Muon::MuonStationIndex::numberOfSectors();
104 if (refSeed.sector() == nSectors && doubSector ==2) {
105 return SectorProjector::leftOverlap;
106 } else if (refSeed.sector() == 1 && doubSector == nSectors) {
107 return SectorProjector::rightOverlap;
108 }
109 return static_cast<SectorProjector>(refSeed.sector() - doubSector );
110 }
112 const xAOD::MuonSegment& segment,
113 const MsTrackSeed& seed) const {
114 return projectOntoSector(gctx, segment, projectorFromSeed(segment, seed));
115 }
117 const xAOD::MuonSegment& segment,
118 const SectorProjector proj) const {
120 const double sectorPhi{projectedPhi(segment.sector(), proj)};
121 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Project onto "<<to_string(proj)<<" in "
122 <<envelope(segment)->identString());
123 return projectOntoPhiPlane(gctx, segment, sectorPhi);
124 }
126 const xAOD::MuonSegment& segment,
127 const double projectPhi) const {
128 using enum SectorProjector;
129 const Amg::Vector3D segPos3D{segment.position()};
132 const Amg::Vector3D dirAlongTube{envelope(segment)->localToGlobalTransform(gctx).linear().col(0)};
133 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Project onto phi: "<<inDeg(projectPhi));
134 const Amg::Vector3D radialDir = Amg::getRotateZ3D(projectPhi) * Amg::Vector3D::UnitX();
135 using namespace Acts::detail::LineHelper;
137 const auto isect = lineIntersect<3>(segPos3D.z()*Amg::Vector3D::UnitZ(), radialDir,
138 segPos3D, dirAlongTube);
139 using namespace Muon::MuonStationIndex;
140 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Projected "<<printID(segment)
141 <<" segment in @"<<Amg::toString(segPos3D)<<" + "
142 <<Amg::toString(segment.direction())<<", chi2: "<<(segment.chiSquared() / std::max(1.f, segment.numberDoF()))
143 <<", nDoF: "<<segment.numberDoF()<<" --> "<<Amg::toString(isect.position()));
144 return isect.position();
145 }
147 const xAOD::MuonSegment& segment,
148 const Location loc,
149 const SectorProjector proj) const {
151 const Amg::Vector3D pos{projectOntoSector(gctx, segment, proj)};
152 const Amg::Vector3D dir{segment.direction()};
153
154 const Amg::Vector2D projPos{pos.perp(), pos.z()};
155 const Amg::Vector2D projDir{dir.perp(), dir.z()};
156
157 ATH_MSG_VERBOSE( "segment position:" << segment.position() << ", direction: " << segment.direction() );
158 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Express segment in @"<<Amg::toString(pos)
159 <<", direction: "<<Amg::toString(dir)<< " sector projector: " << Acts::toUnderlying(proj) << " location: " << Acts::toUnderlying(loc));
160 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Projected position onto sector: "<<Amg::toString(projPos)
161 <<", projected direction: "<<Amg::toString(projDir));
162
163 double lambda{0.};
164 if (Location::Barrel == loc) {
165 lambda = Amg::intersect<2>(projPos, projDir, Amg::Vector2D::UnitX(),
166 m_cfg.barrelRadius).value_or(10. * Gaudi::Units::km);
167 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Intersect with barrel at radius: "<<m_cfg.barrelRadius<<" --> "<<Amg::toString(projPos + lambda * projDir));
168 } else {
169 lambda = Amg::intersect<2>(projPos, projDir, Amg::Vector2D::UnitY(),
170 Acts::copySign(m_cfg.endcapDiscZ, projPos[1])).value_or(10. * Gaudi::Units::km);
171 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Intersect with endcap at z: "<<Acts::copySign(m_cfg.endcapDiscZ, projPos[1])<<" --> "<<Amg::toString(projPos + lambda * projDir));
172 }
173 return projPos + lambda * projDir;
174 }
176 const Location loc) const {
177 using enum Location;
178 if (loc == Barrel && std::abs(projPos[1]) > std::min(m_cfg.endcapDiscZ, m_cfg.barrelLength)) {
179 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Position "<<Amg::toString(projPos)<<
180 " exceeds cylinder boundaries ("<<m_cfg.barrelRadius<<", "
181 <<std::min(m_cfg.endcapDiscZ, m_cfg.barrelLength)<<")");
182 return false;
183 } else if (loc == Endcap && (0 > projPos[0] || projPos[0] > m_cfg.endcapDiscRadius)) {
184 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Position "<<Amg::toString(projPos)<<
185 " exceeds endcap boundaries ("<<m_cfg.endcapDiscRadius<<", "<<Acts::copySign(m_cfg.endcapDiscZ, projPos[1])<<")");
186 return false;
187 }
188 return true;
189 }
190 std::optional<double> MsTrackSeeder::calculateRadius(VecOpt_t&& pI, VecOpt_t&& pM, VecOpt_t&& pO,
191 const Amg::Vector3D& planeNorm) const {
192 if (!pI || !pM || !pO) {
193 return std::nullopt;
194 }
195 const Amg::Vector3D leverL = (*pO) - (*pI);
196
197 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Construct sagitta from "
198 <<"\n --- Inner: "<<print(*pI) <<"\n --- Middle: "<<print(*pM)
199 <<"\n --- Outer: "<<print(*pO));
200 if (std::abs(planeNorm.dot(leverL)) > Acts::s_onSurfaceTolerance ||
201 std::abs(planeNorm.dot( (*pM) - (*pI))) > Acts::s_onSurfaceTolerance) {
202 THROW_EXCEPTION("The lever arm is in the bending plane: "<<Amg::toString(planeNorm)
203 <<", "<<Amg::toString(leverL.unit())<<", "<<Amg::toString( ((*pM) - (*pI)).unit()));
204 }
205 const Amg::Vector3D sagittaDir = leverL.cross(planeNorm).unit();
206 std::optional<double> sagitta = Amg::intersect<3>(*pI, leverL.unit(), *pM, sagittaDir);
207 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Estimated sagitta: "<<(sagitta ?
208 std::to_string(sagitta.value_or(0.)) : "---")<<", lever arm: "
209 <<(leverL.mag() / Gaudi::Units::m)<<" [m]" );
210 if (!sagitta) {
211 return std::nullopt;
212 }
221 return leverL.dot(leverL) / (8. * (*sagitta));
222 }
224 const AtlasFieldCacheCondObj& magField,
225 const MsTrackSeed& seed) const {
226
227 MagField::AtlasFieldCache fieldCache{};
228 magField.getInitializedCache(fieldCache);
229
231 double circPhi{0.};
232 unsigned nSegsWithPhi{0};
233 for (const xAOD::MuonSegment* segment : seed.segments()) {
234 if (segment->nPhiLayers() > 0) {
235 circPhi += std::atan2(segment->y(), segment->x());
236 ++nSegsWithPhi;
237 }
238 }
239 if (!nSegsWithPhi){
240 circPhi = projectedPhi(seed.segments()[0]->sector(),
241 projectorFromSeed(*seed.segments()[0], seed));
242 } else {
243 circPhi /=nSegsWithPhi;
244 }
245 auto layerPos{Acts::filledArray<VecOpt_t, 6>(std::nullopt)};
246 double avgBField{0.}, avgTheta{0.};
247 const xAOD::TruthParticle* truthMuon{nullptr};
248 const Amg::Vector3D planeNorm = Acts::makeDirectionFromPhiTheta(circPhi+ 90._degree, 90._degree);
249
250 {
251 Amg::Vector3D projPos = projectOntoPhiPlane(gctx, *seed.segments()[0], circPhi);
252 Amg::Vector3D projDir = dirForBField(*seed.segments()[0], circPhi);
253 Amg::Vector3D locField{Amg::Vector3D::Zero()};
254 unsigned nBFieldPoints{0};
255 const unsigned nSeedSeg = seed.segments().size();
256 using namespace Muon::MuonStationIndex;
257 for (unsigned int s = 0; s < nSeedSeg; ++s) {
258 avgTheta += projDir.theta();
260 const xAOD::MuonSegment& segment{*seed.segments()[s]};
261 if (!truthMuon) {
262 truthMuon = getTruthMatchedParticle(segment);
263 }
264 switch (toStationIndex(segment.chamberIndex())) {
265 using enum StIndex;
266 case BI:
267 case BE:{
268 average(projPos, layerPos[0]);
269 break;
270 } case BM: {
271 average(projPos, layerPos[1]);
272 break;
273 } case BO: {
274 average(projPos, layerPos[2]);
275 break;
276 } case EI:
277 case EE: {
278 average(projPos, layerPos[3]);
279 break;
280 } case EM : {
281 average(projPos, layerPos[4]);
282 break;
283 } case EO : {
284 average(projPos, layerPos[5]);
285 break;
286 } default: {
287 break;
288 }
289 }
290 if (s + 1 == nSeedSeg) {
291 break;
292 }
294 Amg::Vector3D nextDir = dirForBField(*seed.segments()[s+1], circPhi);
295 Amg::Vector3D nextPos = projectOntoPhiPlane(gctx, *seed.segments()[s+1], circPhi);
296
297 for (double fieldStep : m_cfg.fieldExtpSteps) {
298 /* The */
299 if (fieldStep == 0. && s > 0) {
300 continue;
301 }
302 const Amg::Vector3D extPos = (1. -fieldStep) * projPos + fieldStep * nextPos;
303 const Amg::Vector3D extDir = ((1. -fieldStep) * projDir + fieldStep * nextDir).unit();
304 // const Amg::Vector3D loc
305 fieldCache.getField(extPos.data(), locField.data());
306 const double localB = extDir.cross(locField).mag();
307
308 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Segment "<<s<<", step: "<<fieldStep
309 <<", position: "<<Amg::toString(extPos)<<", dir: "<<Amg::toString(extDir)
310 <<" --> localB: "<<(localB * 1000.)<<" T."
311 <<", B: "<<Amg::toString(locField* 1000.)
312 <<", PxB:"<<Amg::toString(extDir.cross(locField).unit())
313 <<", planeNorm: "<<Amg::toString(planeNorm));
314 avgBField += localB;
315 ++nBFieldPoints;
316 }
317 projPos = std::move(nextPos);
318 projDir = std::move(nextDir);
319 }
320 avgBField /= std::max(nBFieldPoints, 1u);
321 ATH_MSG_DEBUG(__func__<<"() "<<__LINE__<<" - averaged field: "<<(avgBField * 1000.)
322 <<" T, number of points: "<<nBFieldPoints<<".");
323 }
324 avgTheta /= seed.segments().size();
326 std::optional<double> barrelR = calculateRadius(std::move(layerPos[0]),
327 std::move(layerPos[1]),
328 std::move(layerPos[2]), planeNorm);
330 std::optional<double> endcapR = calculateRadius(std::move(layerPos[3]),
331 std::move(layerPos[4]),
332 std::move(layerPos[5]), planeNorm);
334 if (!barrelR && !endcapR) {
335 return 5.*Gaudi::Units::TeV;
336 }
337 const double r = 0.5* ((barrelR ? *barrelR : *endcapR) +
338 (endcapR ? *endcapR : *barrelR));
340 const double P = 0.3* Gaudi::Units::GeV* avgBField * r / std::abs(std::sin(avgTheta));
341
342 ATH_MSG_DEBUG(__func__<<"() "<<__LINE__<<" - Estimated radius "<<r / Gaudi::Units::m<<" [m] --> P: "<<
343 (P / Gaudi::Units::GeV)<<" [GeV]");
344
345 if (truthMuon && Acts::copySign(1.f, P) != truthMuon->charge() && truthMuon->abseta() < 2.5 &&
346 (truthMuon->abseta() < 1.3 || truthMuon->abseta() > 1.4) ) {
347 ATH_MSG_WARNING("Invalid charge, pT: "<<(truthMuon->pt() / Gaudi::Units::GeV)<<" [GeV], eta: "
348 <<truthMuon->eta()<<", phi: "<<(truthMuon->phi() / 1._degree)<<", q: "<<truthMuon->charge());
349 }
350 return P;
351 }
353 const xAOD::MuonSegment* segment,
354 const Location loc,
355 TreeRawVec_t& outContainer) const {
356
357 const int segSector = segment->sector();
358 for (const auto proj : {SectorProjector::leftOverlap, SectorProjector::center, SectorProjector::rightOverlap}) {
360 const int projSector = ringSector(segSector + Acts::toUnderlying(proj));
361 if (segment->nPhiLayers() > 0 && proj != SectorProjector::center &&
362 !sectorMap.insideSector(projSector, segment->position().phi())) {
363 ATH_MSG_VERBOSE(__func__<<"() "<<__LINE__<<" - Segment @"<<Amg::toString(segment->position())
364 <<" is not in sector "<<projSector<<" which is "<<to_string(proj) <<" to "<<segment->sector());
365 continue;
366 }
367 const Amg::Vector2D refPoint{expressOnCylinder(gctx, *segment, loc, proj)};
368 if (!withinBounds(refPoint, loc)) {
369 continue;
370 }
371 using enum SeedCoords;
372 std::array<double, 3> coords{};
376 const int treeSector = 2*segSector + Acts::toUnderlying(proj);
377 coords[Acts::toUnderlying(eSector)] = ringOverlap(treeSector);
380 coords[Acts::toUnderlying(eDetSection)] = Acts::copySign(Acts::toUnderlying(loc), refPoint[1]);
382 coords[Acts::toUnderlying(ePosOnCylinder)] = refPoint[Location::Barrel == loc];
383 ATH_MSG_VERBOSE("Add segment "<<print(*segment)<<", seed quality: "
384 <<m_cfg.selector->passSeedingQuality(Gaudi::Hive::currentContext(), *detailedSegment(*segment))
385 <<" with "<<coords<<" to the search tree");
386 outContainer.emplace_back(std::move(coords), segment);
387 }
388 }
390 const xAOD::MuonSegmentContainer& segments) const{
391 TreeRawVec_t rawData{};
392 rawData.reserve(3*segments.size());
393 for (const xAOD::MuonSegment* segment : segments){
394 appendSegment(gctx, segment, Location::Barrel, rawData);
395 appendSegment(gctx, segment, Location::Endcap, rawData);
396 }
397 ATH_MSG_VERBOSE("Create a new tree with "<<rawData.size()<<" entries. ");
398 return SearchTree_t{std::move(rawData)};
399 }
400 std::unique_ptr<MsTrackSeedContainer> MsTrackSeeder::findTrackSeeds(const EventContext& ctx,
401 const ActsTrk::GeometryContext& gctx,
402 const xAOD::MuonSegmentContainer& segments) const {
403 SearchTree_t orderedSegs{constructTree(gctx, segments)};
404 MsTrackSeedContainer trackSeeds{};
405 using enum SeedCoords;
406 for (const auto& [coords, seedCandidate] : orderedSegs) {
409 const Segment* recoSeedCandidate = detailedSegment(*seedCandidate);
410 if (!m_cfg.selector->passSeedingQuality(ctx, *recoSeedCandidate)){
411 ATH_MSG_VERBOSE("Segment "<<print(*seedCandidate)<<" does not pass the seeding quality.");
412 continue;
413 }
415 SearchTree_t::range_t selectRange{};
418 selectRange[Acts::toUnderlying(eDetSection)].shrink(coords[Acts::toUnderlying(eDetSection)] - 0.1,
419 coords[Acts::toUnderlying(eDetSection)] + 0.1);
421 selectRange[Acts::toUnderlying(ePosOnCylinder)].shrink(coords[Acts::toUnderlying(ePosOnCylinder)] - m_cfg.seedHalfLength,
422 coords[Acts::toUnderlying(ePosOnCylinder)] + m_cfg.seedHalfLength);
424 selectRange[Acts::toUnderlying(eSector)].shrink(coords[Acts::toUnderlying(eSector)] -0.25,
425 coords[Acts::toUnderlying(eSector)] +0.25);
426
427 MsTrackSeed newSeed{static_cast<Location>(std::abs(coords[Acts::toUnderlying(eDetSection)])),
428 static_cast<int>(coords[Acts::toUnderlying(eSector)])};
430 ATH_MSG_VERBOSE("Search for compatible segments to "<<print(*seedCandidate)<<".");
431 orderedSegs.rangeSearchMapDiscard(selectRange, [&](
432 const SearchTree_t::coordinate_t& /*coords*/,
433 const xAOD::MuonSegment* extendWithMe) {
435 const Segment* extendCandidate = detailedSegment(*extendWithMe);
436 if (!m_cfg.selector->compatibleForTrack(ctx, *recoSeedCandidate, *extendCandidate)) {
437 ATH_MSG_VERBOSE("Segment "<<print(*extendWithMe)<<" is not compatible.");
438 return;
439 }
440 auto itr = std::ranges::find_if(newSeed.segments(), [extendWithMe](const xAOD::MuonSegment* onSeed){
441 return extendWithMe->chamberIndex() == onSeed->chamberIndex();
442 });
443 if (itr == newSeed.segments().end()){
444 ATH_MSG_VERBOSE("Add segment "<<print(*extendWithMe)<<" to seed.");
445 newSeed.addSegment(extendWithMe);
446 }
447 else if (reducedChi2(**itr) > reducedChi2(*extendWithMe) &&
448 (*itr)->nPhiLayers() <= extendWithMe->nPhiLayers()) {
449
450 ATH_MSG_VERBOSE("Replace segment "<<print(**itr)<<" with "<<print(*extendWithMe)
451 <<" on seed due to better chi2.");
452 newSeed.replaceSegment(*itr, extendWithMe);
453 }
454 });
456 if (newSeed.segments().empty()) {
457 continue;
458 }
459 newSeed.addSegment(seedCandidate);
460
461 //Check if we have multiple segments from the same station, if so split the seed and create duplicate seeds
462
464 const double r = newSeed.location() == Location::Barrel ? m_cfg.barrelRadius
465 : coords[Acts::toUnderlying(ePosOnCylinder)];
466 const double z = newSeed.location() == Location::Barrel ? coords[Acts::toUnderlying(ePosOnCylinder)]
467 : coords[Acts::toUnderlying(eDetSection)]* m_cfg.endcapDiscZ;
468 const int secCoord = coords[Acts::toUnderlying(eSector)];
469
470 const double phi = sectorMap.sectorOverlapPhi( (secCoord - secCoord % 2) / 2, (secCoord + secCoord % 2) / 2 );
471 Amg::Vector3D pos = r * Acts::makeDirectionFromPhiTheta(phi, 90._degree)
472 + z * Amg::Vector3D::UnitZ();
473
474 newSeed.setPosition(std::move(pos));
475 ATH_MSG_VERBOSE("Add new seed "<<newSeed);
476 trackSeeds.emplace_back(std::move(newSeed));
477 }
478 ATH_MSG_VERBOSE("Found in total "<<trackSeeds.size()<<" before overlap removal");
479 return resolveOverlaps(std::move(trackSeeds));
480 }
481 std::unique_ptr<MsTrackSeedContainer>
483
485 std::ranges::sort(unresolved, [](const MsTrackSeed& a, const MsTrackSeed&b) {
486 return a.segments().size() > b.segments().size();
487 });
488 auto outputSeeds = std::make_unique<MsTrackSeedContainer>();
489 outputSeeds->reserve(unresolved.size());
490 std::ranges::copy_if(std::move(unresolved), std::back_inserter(*outputSeeds),
491 [this, &outputSeeds](const MsTrackSeed& testMe) {
492 MsTrackSeedContainer::iterator test_itr =
493 std::ranges::find_if(*outputSeeds, [&testMe](const MsTrackSeed& good) {
494 if (ringOverlap(good.sector() - testMe.sector()) > 1) {
495 return false;
496 }
497 return std::ranges::find_if(testMe.segments(),
498 [&good](const xAOD::MuonSegment* segInTest) {
499 return std::ranges::find(good.segments(), segInTest) != good.segments().end();
500 }) != testMe.segments().end();
501 });
503 if (test_itr == outputSeeds->end()) {
504 ATH_MSG_VERBOSE("Add new seed "<<testMe);
505 return true;
506 }
508 if ( (*test_itr).segments().size() < testMe.segments().size()){
509 (*test_itr) = testMe;
510 }
511 return false;
512 });
513
514 ATH_MSG_VERBOSE("Found in total "<<outputSeeds->size()<<" after overlap removal");
515 return outputSeeds;
516 }
517}
Scalar phi() const
phi method
Scalar theta() const
theta method
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t a
static Double_t P(Double_t *tt, Double_t *par)
void print(char *figname, TCanvas *c1)
#define x
#define z
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
size_type size() const noexcept
Returns the number of elements in the collection.
Local cache for magnetic field (based on MagFieldServices/AtlasFieldSvcTLS.h)
void getField(const double *ATH_RESTRICT xyz, double *ATH_RESTRICT bxyz, double *ATH_RESTRICT deriv=nullptr)
get B field value at given position xyz[3] is in mm, bxyz[3] is in kT if deriv[9] is given,...
A spectrometer sector forms the envelope of all chambers that are placed in the same MS sector & laye...
const Amg::Transform3D & localToGlobalTransform(const ActsTrk::GeometryContext &gctx) const
Returns the local -> global tarnsformation from the sector.
int sector() const
Returns the seed's sector.
Definition MsTrackSeed.h:55
void addSegment(const xAOD::MuonSegment *seg)
Append a segment to the seed.
Location location() const
Returns the location of the seed.
const std::vector< const xAOD::MuonSegment * > & segments() const
Returns the vector of associated segments.
void replaceSegment(const xAOD::MuonSegment *exist, const xAOD::MuonSegment *updated)
Replaces an already added segment in the seed with a better suited one.
void setPosition(Amg::Vector3D &&pos)
set the seed's position
std::optional< double > calculateRadius(VecOpt_t &&pI, VecOpt_t &&pM, VecOpt_t &&pO, const Amg::Vector3D &planeNorm) const
Calculates the radius of the bending circle from three points using the sagitta.
bool withinBounds(const Amg::Vector2D &projPos, const Location loc) const
Returns whether the expression on the cylinder is within the surface bounds.
static std::string to_string(const SectorProjector proj)
Print the sector projector.
SearchTree_t constructTree(const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegmentContainer &segments) const
Construct a complete search tree from a MuonSegment container.
static int ringSector(const int sector)
Ensure that the parsed sector number is following the MS sector schema 0 is mapped to 16 and 17 is ma...
Amg::Vector3D projectOntoSector(const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegment &segment, const SectorProjector proj) const
Projects the segment's position onto the sector centre or onto the overlap point with one of the neig...
SearchTree_t::vector_t TreeRawVec_t
Abbrivation of the KDTree raw data vector.
Amg::Vector2D expressOnCylinder(const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegment &segment, const Location loc, const SectorProjector proj) const
Expresses the segment on the cylinder surface.
std::unique_ptr< MsTrackSeedContainer > resolveOverlaps(MsTrackSeedContainer &&unresolved) const
Removes exact duplciates or partial subsets of the MsTrackSeeds.
static int ringOverlap(const int sector)
Maps the sector 33 -> 0 to close the extended MS symmetry ring.
void appendSegment(const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegment *segment, const Location loc, TreeRawVec_t &outContainer) const
Append the to the raw data container.
Acts::KDTree< 3, const xAOD::MuonSegment *, double, std::array, 6 > SearchTree_t
Definition of the search tree class.
const MuonGMR4::SpectrometerSector * envelope(const xAOD::MuonSegment &segment) const
Returns the spectrometer envelope associated to the segment (Coord system where the parameter are exp...
static double projectedPhi(const int sector, const SectorProjector proj)
Returns the projected phi for a given sector and projector.
MsTrackSeeder(const std::string &msgName, Config &&cfg)
Standard constructor.
MsTrackSeed::Location Location
Enum toggling whether the segment is in the endcap or barrel.
SeedCoords
Abrivation of the seed coordinates.
@ eSector
Sector of the associated spectrometer sector.
@ ePosOnCylinder
Extrapolation position along the cylinder surface.
@ eDetSection
Encode the seed location (-1,1 -> endcaps, 0 -> barrel.
double estimateQtimesP(const ActsTrk::GeometryContext &gctx, const AtlasFieldCacheCondObj &magField, const MsTrackSeed &seed) const
Estimate the q /p of the seed candidate from the contained segments.
std::unique_ptr< MsTrackSeedContainer > findTrackSeeds(const EventContext &ctx, const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegmentContainer &segments) const
Constructs the MS track seeds from the segment container.
Amg::Vector3D projectOntoPhiPlane(const ActsTrk::GeometryContext &gctx, const xAOD::MuonSegment &segment, const double projectPhi) const
static SectorProjector projectorFromSeed(const xAOD::MuonSegment &seg, const MsTrackSeed &refSeed)
Returns the Sector projector within the context of a MsTrackSeed.
std::optional< Amg::Vector3D > VecOpt_t
SectorProjector
Enumeration to select the sector projection.
@ rightOverlap
Project the segment onto the sector centre.
@ center
Project the segment onto the overlap with the previous sector.
Placeholder for what will later be the muon segment EDM representation.
float px() const
float pz() const
Returns the pz.
float numberDoF() const
Returns the numberDoF.
int nPrecisionHits() const
Amg::Vector3D direction() const
Returns the direction as Amg::Vector.
float chiSquared() const
float py() const
Returns the py.
::Muon::MuonStationIndex::ChIndex chamberIndex() const
Returns the chamber index.
Amg::Vector3D position() const
Returns the position as Amg::Vector.
int nPhiLayers() const
Returns the number of phi layers.
int etaIndex() const
Returns the eta index, which corresponds to stationEta in the offline identifiers (and the ).
int r
Definition globals.cxx:22
std::optional< double > intersect(const AmgVector(N)&posA, const AmgVector(N)&dirA, const AmgVector(N)&posB, const AmgVector(N)&dirB)
Calculates the point B' along the line B that's closest to a second line A.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Amg::Transform3D getRotateZ3D(double angle)
get a rotation transformation around Z-axis
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
This header ties the generic definitions in this package.
const xAOD::TruthParticle * getTruthMatchedParticle(const xAOD::MuonSegment &segment)
Returns the particle truth-matched to the segment.
std::vector< MsTrackSeed > MsTrackSeedContainer
Definition MsTrackSeed.h:63
std::string printID(const xAOD::MuonSegment &seg)
Print the chamber ID of a segment, e.g.
MsTrackSeeder::SearchTree_t SearchTree_t
MsTrackSeeder::SectorProjector SectorProjector
const Segment * detailedSegment(const xAOD::MuonSegment &seg)
Helper function to navigate from the xAOD::MuonSegment to the MuonR4::Segment.
StIndex
enum to classify the different station layers in the muon spectrometer
StIndex toStationIndex(ChIndex index)
convert ChIndex into StIndex
constexpr unsigned numberOfSectors()
return total number of sectors
STL namespace.
MuonSegmentContainer_v1 MuonSegmentContainer
Definition of the current "MuonSegment container version".
TruthParticle_v1 TruthParticle
Typedef to implementation.
MuonSegment_v1 MuonSegment
Reference the current persistent version:
Configuration object.
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10