ATLAS Offline Software
Loading...
Searching...
No Matches
SegmentFitterEventData.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
5
10#include <GaudiKernel/PhysicalConstants.h>
11
12#include "Acts/Surfaces/PlaneSurface.hpp"
13#include "Acts/Definitions/Units.hpp"
14
15
16#include <sstream>
17#include <format>
18using namespace Acts;
19using namespace Acts::UnitLiterals;
20
21namespace {
22 constexpr double straightQoverP = 1. / 20._TeV;
23}
24
25namespace MuonR4{
26 double houghTanBeta(const Amg::Vector3D& v){
27 constexpr double eps = std::numeric_limits<float>::epsilon();
28 return v.y() / ( std::abs(v.z()) > eps ? v.z() : eps);
29 }
30 double houghTanAlpha(const Amg::Vector3D& v){
31 constexpr double eps = std::numeric_limits<float>::epsilon();
32 return v.x() / ( std::abs(v.z()) > eps ? v.z() : eps);
33 }
34 namespace SegmentFit {
35 std::pair<Amg::Vector3D, Amg::Vector3D> makeLine(const Parameters& pars) {
36 using enum ParamDefs;
37 return std::make_pair(Amg::Vector3D(pars[toUnderlying(x0)],
38 pars[toUnderlying(y0)],0.),
39 Amg::dirFromAngles(pars[toUnderlying(phi)],
40 pars[toUnderlying(theta)]));
41 }
43 static const SG::Accessor<xAOD::MeasVector<toUnderlying(ParamDefs::nPars)>> acc{"localSegPars"};
44 Parameters segPars{};
45 for (std::size_t p =0 ; p < segPars.size(); ++p) {
46 segPars[p] = acc(seg)[p];
47 }
48 return segPars;
49 }
51 const Segment& segment) {
52 Parameters pars{};
53 const Amg::Transform3D globToLoc = segment.msSector()->globalToLocalTrans(gctx);
54 const Amg::Vector3D locPos = globToLoc * segment.position();
55 const Amg::Vector3D locDir = globToLoc.linear() * segment.direction();
56 pars[toUnderlying(ParamDefs::x0)] = locPos.x();
57 pars[toUnderlying(ParamDefs::y0)] = locPos.y();
58 pars[toUnderlying(ParamDefs::theta)] = locDir.theta();
59 pars[toUnderlying(ParamDefs::phi)] = locDir.phi();
60 pars[toUnderlying(ParamDefs::t0)] = segment.segementT0();
61 return pars;
62 }
63
64 std::string makeLabel(const Parameters&pars) {
65 std::stringstream sstr{};
66 sstr<<std::format("x_{{0}}={:.2f}", pars[toUnderlying(ParamDefs::x0)])<<", ";
67 sstr<<std::format("y_{{0}}={:.2f}", pars[toUnderlying(ParamDefs::y0)])<<", ";
68 sstr<<std::format("#theta={:.2f}^{{#circ}}", pars[toUnderlying(ParamDefs::theta)] / Gaudi::Units::deg )<<", ";
69 sstr<<std::format("#phi={:.2f}^{{#circ}}", pars[toUnderlying(ParamDefs::phi)] / Gaudi::Units::deg)<<", ";
70 sstr<<std::format("t_{{0}}={:.1f}", pars[toUnderlying(ParamDefs::t0)]);
71 return sstr.str();
72 }
73 std::string toString(const Parameters& pars) {
74 std::stringstream sstr{};
75 sstr<< std::format("{}={:.2f}, ",toString(ParamDefs::x0), pars[toUnderlying(ParamDefs::x0)]);
76 sstr<< std::format("{}={:.2f}, ",toString(ParamDefs::y0), pars[toUnderlying(ParamDefs::y0)]);
77 sstr<< std::format("{}={:.2f}, ",toString(ParamDefs::theta), pars[toUnderlying(ParamDefs::theta)]/Gaudi::Units::deg);
78 sstr<< std::format("{}={:.2f}, ",toString(ParamDefs::phi), pars[toUnderlying(ParamDefs::phi)]/Gaudi::Units::deg);
79 sstr<< std::format("{}={:.2f}",toString(ParamDefs::t0), pars[toUnderlying(ParamDefs::t0)]);
80 return sstr.str();
81 }
82 std::string toString(const ParamDefs a) {
83 return SeedingAux::parName(a);
84 }
85 Acts::BoundTrackParameters boundSegmentPars(const MuonGMR4::MuonDetectorManager& detMgr,
86 const xAOD::MuonSegment& segment,
87 std::optional<Acts::BoundMatrix> cov,
88 Acts::ParticleHypothesis hypot) {
89
90 const auto* msSector = detMgr.getSectorEnvelope(segment.chamberIndex(),
91 segment.sector(),
92 segment.etaIndex());
93 const Acts::Surface& surface = msSector->surface();
94
95 const auto locSegPars = localSegmentPars(segment);
96
97 const Amg::Vector3D globDir = segment.direction();
98
99
100 Acts::BoundVector boundPars{};
101 boundPars[Acts::eBoundLoc0] = locSegPars[toUnderlying(ParamDefs::x0)];
102 boundPars[Acts::eBoundLoc1] = locSegPars[toUnderlying(ParamDefs::y0)];
103 boundPars[Acts::eBoundPhi] = globDir.phi();
104 boundPars[Acts::eBoundTheta] = globDir.theta();
105 boundPars[Acts::eBoundQOverP] = straightQoverP;
106 boundPars[Acts::eBoundTime] = ActsTrk::timeToActs(segment.position().mag() / Gaudi::Units::c_light +
107 segment.t0());
108
109 return Acts::BoundTrackParameters{surface.getSharedPtr(), std::move(boundPars),
110 cov, hypot};
111 }
112 Acts::BoundTrackParameters boundSegmentPars(const ActsTrk::GeometryContext& gctx,
113 const Segment& segment,
114 const Acts::ParticleHypothesis hypot) {
115 const auto& surface = segment.msSector()->surface();
116
117 const Amg::Vector3D locPos = surface.transform(gctx.context()).inverse() *
118 segment.position();
119 Acts::BoundVector boundPars{};
120 boundPars[Acts::eBoundLoc0] = locPos.x();
121 boundPars[Acts::eBoundLoc1] = locPos.y();
122 boundPars[Acts::eBoundPhi] = segment.direction().phi();
123 boundPars[Acts::eBoundTheta] = segment.direction().theta();
124 boundPars[Acts::eBoundQOverP] = straightQoverP;
125 boundPars[Acts::eBoundTime] = ActsTrk::timeToActs(segment.position().mag() / Gaudi::Units::c_light +
126 segment.segementT0());
127 Acts::BoundMatrix cov{Acts::BoundMatrix::Identity()};
128
129 return Acts::BoundTrackParameters{surface.getSharedPtr(), std::move(boundPars),
130 cov, hypot};
131 }
132
133 }
134}
Scalar phi() const
phi method
Scalar theta() const
theta method
static Double_t a
Acts::GeometryContext context() const
const SpectrometerSector * getSectorEnvelope(const Identifier &channelId) const
Retrieves the spectrometer envelope enclosing the channel's readout element.
Amg::Transform3D globalToLocalTrans(const ActsTrk::GeometryContext &gctx) const
Returns the global -> local transformation from the ATLAS global.
const Acts::PlaneSurface & surface() const
Returns the associated surface.
Placeholder for what will later be the muon segment EDM representation.
double segementT0() const
Returns the fitted segment time, if there's any.
const MuonGMR4::SpectrometerSector * msSector() const
Returns the associated MS sector.
const Amg::Vector3D & position() const
Returns the global segment position.
const Amg::Vector3D & direction() const
Returns the global segment direction.
Helper class to provide type-safe access to aux data.
float t0() const
Amg::Vector3D direction() const
Returns the direction as Amg::Vector.
::Muon::MuonStationIndex::ChIndex chamberIndex() const
Returns the chamber index.
Amg::Vector3D position() const
Returns the position as Amg::Vector.
int etaIndex() const
Returns the eta index, which corresponds to stationEta in the offline identifiers (and the ).
constexpr double timeToActs(const double athenaT)
Converts a time unit from Athena to Acts units.
Amg::Vector3D dirFromAngles(const double phi, const double theta)
Constructs a direction vector from the azimuthal & polar angles.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
SeedingAux::FitParIndex ParamDefs
Use the same parameter indices as used by the CompSpacePointAuxiliaries.
Parameters localSegmentPars(const xAOD::MuonSegment &seg)
Returns the localSegPars decoration from a xAODMuon::Segment.
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 makeLabel(const Parameters &pars)
Dumps the parameters into a string in the form of TLatex.
std::string toString(const Parameters &pars)
Dumps the parameters into a string with labels in front of each number.
Acts::BoundTrackParameters boundSegmentPars(const MuonGMR4::MuonDetectorManager &detMgr, const xAOD::MuonSegment &segment, std::optional< Acts::BoundMatrix > cov=std::nullopt, Acts::ParticleHypothesis hypot=Acts::ParticleHypothesis::muon())
Returns the segment parameters as boundTrackParameters.
This header ties the generic definitions in this package.
double houghTanBeta(const Amg::Vector3D &v)
Returns the hough tanBeta [y] / [z].
double houghTanAlpha(const Amg::Vector3D &v)
: Returns the hough tanAlpha [x] / [z]
Eigen::Matrix< float, N, 1 > MeasVector
Abrivation of the Matrix & Covariance definitions.
MuonSegment_v1 MuonSegment
Reference the current persistent version: