ATLAS Offline Software
Loading...
Searching...
No Matches
MuonTrackToSegmentTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <set>
8
19#include "TrkTrack/Track.h"
21
22namespace Muon {
23 MuonTrackToSegmentTool::MuonTrackToSegmentTool(const std::string& t, const std::string& n, const IInterface* p) : AthAlgTool(t, n, p) {
24 declareInterface<IMuonTrackToSegmentTool>(this);
25 }
26
28 ATH_CHECK(m_propagator.retrieve());
29 ATH_CHECK(m_idHelperSvc.retrieve());
30 ATH_CHECK(m_edmHelperSvc.retrieve());
31 ATH_CHECK(m_printer.retrieve());
32 ATH_CHECK(m_chamberGeoKey.initialize());
33 ATH_CHECK(m_DetectorManagerKey.initialize());
34 return StatusCode::SUCCESS;
35 }
36
37 MuonSegment* MuonTrackToSegmentTool::convert(const EventContext& ctx, const Trk::Track& track) const {
38 // Avoids FPE with clang.
41
42 ATH_MSG_DEBUG(" creating MuonSegment from track ");
43
44 const Trk::Perigee* perigee = track.perigeeParameters();
45 if (!perigee) {
46 ATH_MSG_WARNING(" was expecting a perigee here... ");
47 return nullptr;
48 }
49
50 const Trk::FitQuality* fq = track.fitQuality();
51 if (!fq) {
52 ATH_MSG_WARNING(" was expecting a FitQuality here... ");
53 return nullptr;
54 }
55
56 std::set<Identifier> chIds;
57
58 // copy rots, get surface
60 rots.reserve(track.measurementsOnTrack()->size());
61
62 // loop over TSOS
63 const Trk::TrackStates* states = track.trackStateOnSurfaces();
64 if (!states) {
65 ATH_MSG_WARNING(" track without states, discarding track ");
66 return nullptr;
67 }
68 // track direction vector
69 const Amg::Vector3D dir = perigee->momentum().unit();
70
71 const Amg::Transform3D* surfaceTransform = nullptr;
72 const Amg::Transform3D* backupTransform = nullptr;
73 std::unique_ptr<Amg::Transform3D> surfaceTransformToBeDeleted;
74 double weightedDistanceSquared{0}, weightSquared{0};
75 for (const Trk::TrackStateOnSurface* tsos : *states) {
76 // require TrackParameters
77 const Trk::TrackParameters* pars = tsos->trackParameters();
78 if (!pars) continue;
79
80 // check whether state is a measurement
81 const Trk::MeasurementBase* meas = tsos->measurementOnTrack();
82 if (!meas || tsos->type(Trk::TrackStateOnSurface::Outlier)) continue;
83 rots.push_back(meas->clone());
84
85 // only consider eta hits
86 Identifier id = m_edmHelperSvc->getIdentifier(*meas);
87 if (!id.is_valid() || m_idHelperSvc->measuresPhi(id)) continue;
88
89 double distance = (pars->position() - perigee->position()).dot(dir);
90 double weight = 1. / meas->localCovariance()(Trk::locX, Trk::locX);
91 ATH_MSG_VERBOSE(" distance " << distance << " error " << Amg::error(meas->localCovariance(), Trk::locX) << " weight " << weight
92 << " " << m_idHelperSvc->toString(id));
93 weightedDistanceSquared += distance * weight;
94 weightSquared += weight;
95 if (m_idHelperSvc->isMdt(id)) {
96 chIds.insert(m_idHelperSvc->chamberId(id));
97 if (!surfaceTransform) {
98 const MdtDriftCircleOnTrack* mdt = dynamic_cast<const MdtDriftCircleOnTrack*>(meas);
99 if (mdt) {
100 // create new surface using AMDB reference frame
101 surfaceTransformToBeDeleted = std::make_unique<Amg::Transform3D>(mdt->detectorElement()->AmdbLRSToGlobalTransform().rotation());
102 surfaceTransformToBeDeleted->pretranslate(mdt->detectorElement()->center());
103 surfaceTransform = surfaceTransformToBeDeleted.get();
104 }
105 }
106 } else if ((m_idHelperSvc->isMM(id) || m_idHelperSvc->isCsc(id)) && !surfaceTransform) {
107 surfaceTransform = &(meas)->associatedSurface().transform();
108 } else if (!surfaceTransform && !backupTransform) {
109 backupTransform = &(meas)->associatedSurface().transform();
110 }
111 }
112 if (!surfaceTransform) surfaceTransform = backupTransform;
113 // calculate distance new reference point, shift it 100 mm towards the start of the segment
114 double refDistance = (weightSquared > 0 ? weightedDistanceSquared / weightSquared : 1) - 100;
115 ATH_MSG_DEBUG(" weighted distance " << refDistance);
116
117 const Amg::Vector3D refPos = perigee->position() + refDistance * dir;
118
119 // find closest measured parameters
120 double minDist = -1e6;
121 const Trk::TrackParameters* closestPars = nullptr;
122 for (const Trk::TrackStateOnSurface* tsos : *states) {
123 // require TrackParameters
124 const Trk::TrackParameters* pars = tsos->trackParameters();
125 if (!pars || !pars->covariance()) continue;
126
127 // look for the closest measured parameters to the reference point
128 double distance = (pars->position() - refPos).dot(dir);
129 if (distance < 0 && std::abs(distance) < std::abs(minDist)) {
130 minDist = distance;
131 closestPars = pars;
132 }
133 }
134
135 if (!surfaceTransform) {
136 ATH_MSG_DEBUG(" failed to create a PlaneSurface for the track, cannot make segment!!! " << std::endl
137 << m_printer->print(track)
138 << std::endl
139 << m_printer->printStations(track));
140 return nullptr;
141 }
142
143 if (!closestPars) {
144 closestPars = perigee;
145 minDist = (perigee->position() - refPos).dot(dir);
146 }
147
148 Amg::Transform3D transform(surfaceTransform->rotation());
149 transform.pretranslate(refPos);
150 constexpr double surfDim = 500.;
151 std::unique_ptr<Trk::PlaneSurface> surf = std::make_unique<Trk::PlaneSurface>(transform, surfDim, surfDim);
152 std::unique_ptr<Trk::TrackParameters> exPars = m_propagator->propagate(ctx, *closestPars, *surf, minDist > 0 ? Trk::oppositeMomentum : Trk::alongMomentum, false,
154 if (!exPars || !exPars->covariance()) {
155 ATH_MSG_VERBOSE("First trial reaching the surface failed. This is presumably due to a too large momentum. Let's try with a dummy 1 GeV momentum");
156 std::unique_ptr<Trk::TrackParameters> cloned_pars {closestPars->clone()};
157 constexpr double OneOverGeV = 1./ Gaudi::Units::GeV;
158 cloned_pars->parameters()[Trk::qOverP] = cloned_pars->charge() *OneOverGeV;
159 exPars = m_propagator->propagate(ctx, *cloned_pars, *surf, Trk::anyDirection, false,
161
162 if (!exPars){
163 ATH_MSG_DEBUG(" propagation failed!!! "<<*cloned_pars<<std::endl<<std::endl<<*surf);
164 return nullptr;
165 }
166 exPars->parameters()[Trk::qOverP] = closestPars->parameters()[Trk::qOverP];
167 }
168 Amg::Vector2D locPos{Amg::Vector2D::Zero()};
169 if (!surf->globalToLocal(exPars->position(), exPars->momentum(), locPos)) {
170 ATH_MSG_WARNING(" localToGlobal failed!!! ");
171 return nullptr;
172 }
173 Trk::LocalDirection locDir;
174 surf->globalToLocalDirection(exPars->momentum(), locDir);
175
176 // convert errors on global angles theta/phi to errors on local angles angleYZ/angleXZ
177 Trk::JacobianPhiThetaLocalAngles globalToLocalMeasAnglesJacobian(exPars->parameters()[Trk::phi], exPars->parameters()[Trk::theta],
178 exPars->associatedSurface().transform().rotation().inverse());
179
180 // make the Jacobian to convert all in one go from global to local
181 // so that the correlations are calculated correctly
182 AmgSymMatrix(5) globalToLocalMeasJacobian{AmgSymMatrix(5)::Zero()};
183 globalToLocalMeasJacobian(Trk::locX, Trk::locX) = 1.0;
184 globalToLocalMeasJacobian(Trk::locY, Trk::locY) = 1.0;
185 globalToLocalMeasJacobian(Trk::phi, Trk::phi) = globalToLocalMeasAnglesJacobian(0, 0);
186 globalToLocalMeasJacobian(Trk::theta, Trk::theta) = globalToLocalMeasAnglesJacobian(1, 1);
187 globalToLocalMeasJacobian(Trk::theta, Trk::phi) = globalToLocalMeasAnglesJacobian(0, 1); // also fills (Trk::phi,Trk::theta)
188 globalToLocalMeasJacobian(Trk::phi, Trk::theta) = globalToLocalMeasJacobian(Trk::theta, Trk::phi); // also fills (Trk::theta,Trk::phi)
189 globalToLocalMeasJacobian(Trk::qOverP, Trk::qOverP) = 1.0;
190
191 AmgSymMatrix(5) cov = exPars->covariance()->similarity(globalToLocalMeasJacobian);
192
193 Trk::FitQuality* quality = nullptr;
194 if (!chIds.empty()) {
195 // calculate holes
196 std::vector<Identifier> holes;
197 for (const Identifier& chid : chIds) {
198 std::vector<Identifier> holesChamber = calculateHoles(ctx, chid, *exPars, rots.stdcont());
199 holes.insert(holes.end(), holesChamber.begin(), holesChamber.end());
200 }
201 quality = new MuonSegmentQuality(fq->chiSquared(), fq->numberDoF(), holes);
202
203 } else {
204 quality = new Trk::FitQuality(fq->chiSquared(), fq->numberDoF());
205 }
206 MuonSegment* seg = new MuonSegment(
207 locPos, locDir, cov, surf.release(), std::move(rots), quality);
208 return seg;
209 }
210
211 std::vector<Identifier> MuonTrackToSegmentTool::calculateHoles(const EventContext& ctx, const Identifier& chid,
212 const Trk::TrackParameters& pars,
213 const MuonTrackToSegmentTool::MeasVec& measurements) const {
215 if (!InterSectSvc.isValid()) {
216 ATH_MSG_ERROR("Failed to retrieve chamber intersection service");
217 }
218
220 const MuonGM::MuonDetectorManager* MuonDetMgr = detMgr.cptr();
221
222 const MuonStationIntersect intersect = InterSectSvc->tubesCrossedByTrack(MuonDetMgr, chid, pars.position(), pars.momentum().unit());
223
224 // set to identify the hit on the segment
225 std::set<Identifier> hitsOnSegment;
226 for (const Trk::MeasurementBase* mit : measurements) {
227 const MdtDriftCircleOnTrack* mdt = dynamic_cast<const MdtDriftCircleOnTrack*>(mit);
228 if (mdt) hitsOnSegment.insert(mdt->identify());
229 }
230
231 // clear hole vector
232 std::vector<Identifier> holes;
233 for (unsigned int ii = 0; ii < intersect.tubeIntersects().size(); ++ii) {
234 const MuonTubeIntersect& tint = intersect.tubeIntersects()[ii];
235
236 // skip hole check if there is a hit in this tube
237 if (hitsOnSegment.count(tint.tubeId)) continue;
238
239 // if track goes through a tube which did not have a hit count as hole
240 if (std::abs(tint.rIntersect) < MuonDetMgr->getMdtReadoutElement(tint.tubeId)->innerTubeRadius() && tint.xIntersect < -200.) {
241 holes.push_back(tint.tubeId);
242 }
243 }
244 return holes;
245 }
246
247} // namespace Muon
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
MuonSegment_v1 MuonSegment
Reference the current persistent version:
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Derived DataVector<T>.
Definition DataVector.h:795
void reserve(size_type n)
Attempt to preallocate enough memory for a specified number of elements.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
const PtrVector & stdcont() const
Return the underlying std::vector of the container.
double innerTubeRadius() const
Returns the inner tube radius excluding the aluminium walls.
virtual const Amg::Vector3D & center(const Identifier &) const override final
Return the center of the surface associated with this identifier In the case of silicon it returns th...
The MuonDetectorManager stores the transient representation of the Muon Spectrometer geometry and pro...
const MdtReadoutElement * getMdtReadoutElement(const Identifier &id) const
access via extended identifier (requires unpacking)
This class represents the corrected MDT measurements, where the corrections include the effects of wi...
virtual const MuonGM::MdtReadoutElement * detectorElement() const override final
Returns the detector element, assoicated with the PRD of this class.
This is the common muon segment quality object.
This is the common class for 3D segments used in the muon spectrometer.
StatusCode initialize()
initialize method, method taken from bass-class AlgTool
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
SG::ReadCondHandleKey< MuonGM::MuonDetectorManager > m_DetectorManagerKey
std::vector< const Trk::MeasurementBase * > MeasVec
MuonTrackToSegmentTool(const std::string &, const std::string &, const IInterface *)
default AlgTool constructor
std::vector< Identifier > calculateHoles(const EventContext &ctx, const Identifier &chid, const Trk::TrackParameters &pars, const MeasVec &measurements) const
calculate holes
MuonSegment * convert(const EventContext &ctx, const Trk::Track &track) const
convert track to segment
SG::ReadCondHandleKey< Muon::MuonIntersectGeoData > m_chamberGeoKey
PublicToolHandle< MuonEDMPrinterTool > m_printer
ServiceHandle< IMuonEDMHelperSvc > m_edmHelperSvc
ToolHandle< Trk::IPropagator > m_propagator
const_pointer_type cptr()
Class to represent and store fit qualities from track reconstruction in terms of and number of degre...
Definition FitQuality.h:97
int numberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as integer
Definition FitQuality.h:60
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
This jacobian transforms the polar and azimutal angle of a global momentum representation to a local...
represents the three-dimensional global direction with respect to a planar surface frame.
magnetic field properties to steer the behavior of the extrapolation
This class is the pure abstract base class for all fittable tracking measurements.
virtual MeasurementBase * clone() const =0
Pseudo-Constructor.
virtual bool type(MeasurementBaseType::Type type) const =0
Interface method checking the type.
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
const Amg::Vector3D & momentum() const
Access method for the momentum.
virtual ParametersBase< DIM, T > * clone() const override=0
clone method for polymorphic deep copy
const Amg::Vector3D & position() const
Access method for the position.
Identifier identify() const
return the identifier -extends MeasurementBase
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
@ oppositeMomentum
@ alongMomentum
@ anyDirection
DataVector< const Trk::TrackStateOnSurface > TrackStates
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ NoField
Field is set to 0., 0., 0.,.
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
ParametersBase< TrackParametersDim, Charged > TrackParameters
Tell the compiler to optimize assuming that FP may trap.
#define CXXUTILS_TRAPPING_FP
Definition trapping_fp.h:24