ATLAS Offline Software
Loading...
Searching...
No Matches
FPGAActsTrkConverter.cxx
Go to the documentation of this file.
1// Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
2
5#include "Acts/Surfaces/PerigeeSurface.hpp"
8
9#include <format>
10#include <sstream>
11//Heavily inspired by https://gitlab.cern.ch/atlas/athena/-/blob/main/Tracking/Acts/ActsTrackReconstruction/src/RandomProtoTrackCreator.cxx
12
14 const std::string& name,
15 const IInterface* parent): base_class(type,name,parent) { }
16
17
19
20 ATH_MSG_DEBUG("Initializing FPGAActsTrkConverter...");
21
22 // Get SCT & pixel Identifier helpers
23 ATH_CHECK(detStore()->retrieve(m_pixelId, "PixelID"));
24 ATH_CHECK(detStore()->retrieve(m_SCTId, "SCT_ID"));
25
26 return StatusCode::SUCCESS;
27
28}
29
30StatusCode FPGAActsTrkConverter::findProtoTracks(const EventContext& ctx,
31 const xAOD::PixelClusterContainer & pixelContainer,
32 const xAOD::StripClusterContainer & stripContainer,
33 std::vector<ActsTrk::ProtoTrack> & foundProtoTracks,
34 const std::vector<std::vector<FPGATrackSimHit>>& hitsInRoads,
35 const std::vector<FPGATrackSimRoad>& roads) const {
36
37 ATH_MSG_INFO("Creating Acts proto-tracks from FPGA roads...");
38
39 if (hitsInRoads.size() > 0) {
40 std::multimap<xAOD::DetectorIdentType, const xAOD::PixelCluster*> pixelClusterMap;
41 for (const xAOD::PixelCluster* cluster : pixelContainer) {
42 pixelClusterMap.emplace(cluster->identifier(), cluster);
43 }
44
45 std::multimap<xAOD::DetectorIdentType, const xAOD::StripCluster*> stripClusterMap;
46 for (const xAOD::StripCluster* cluster : stripContainer) {
47 stripClusterMap.emplace(cluster->identifier(), cluster);
48
49 }
50 for(size_t roadIndex=0; roadIndex<=hitsInRoads.size()-1;roadIndex++) {
51 std::vector<ActsTrk::ATLASUncalibSourceLink> points;
52 ATH_CHECK(findPrototrackMeasurements(ctx, pixelContainer, stripContainer, pixelClusterMap, stripClusterMap, points, hitsInRoads.at(roadIndex)));
53 if (points.size()) {
54 std::unique_ptr<Acts::BoundTrackParameters> inputPerigee = makeParams(roads.at(roadIndex));
55 foundProtoTracks.emplace_back(points, std::move(inputPerigee));
56 ATH_MSG_INFO("Made a prototrack with " << points.size() << " measurements");
57 }
58 }
59 }
60
61 return StatusCode::SUCCESS;
62}
63
64StatusCode FPGAActsTrkConverter::findProtoTracks(const EventContext& ctx,
65 const xAOD::PixelClusterContainer& pixelContainer,
66 const xAOD::StripClusterContainer& stripContainer,
67 std::vector<ActsTrk::ProtoTrack>& foundProtoTracks,
68 const std::vector<FPGATrackSimTrack>& tracks) const {
69
70 ATH_MSG_INFO("Creating Acts proto-tracks from FPGA tracks...");
71 // Initialize multimaps for pixel and strip clusters
72 std::multimap<xAOD::DetectorIdentType, const xAOD::PixelCluster*> pixelClusterMap;
73 for (const xAOD::PixelCluster* cluster : pixelContainer) {
74 pixelClusterMap.emplace(cluster->identifier(), cluster);
75 }
76
77 std::multimap<xAOD::DetectorIdentType, const xAOD::StripCluster*> stripClusterMap;
78 for (const xAOD::StripCluster* cluster : stripContainer) {
79 stripClusterMap.emplace(cluster->identifier(), cluster);
80 }
81 for (const FPGATrackSimTrack& track : tracks) {
82 if (not track.passedOR()) continue;
83 std::vector<ActsTrk::ATLASUncalibSourceLink> points;
84 const std::vector <FPGATrackSimHit>& hits = track.getFPGATrackSimHits();
85 ATH_CHECK(findPrototrackMeasurements(ctx, pixelContainer, stripContainer, pixelClusterMap, stripClusterMap, points, hits));
86 if (points.size()) {
87 ATH_MSG_DEBUG("\tMaking a proto-track with " << points.size() << " clusters");
88 std::unique_ptr<Acts::BoundTrackParameters> inputPerigee = makeParams(track);
89 foundProtoTracks.emplace_back(points, std::move(inputPerigee));
90 }
91 }
92 return StatusCode::SUCCESS;
93}
94
95StatusCode FPGAActsTrkConverter::findPrototrackMeasurements( const EventContext& ctx,
96 const xAOD::PixelClusterContainer& pixelContainer,
97 const xAOD::StripClusterContainer& stripContainer,
98 const std::multimap<xAOD::DetectorIdentType, const xAOD::PixelCluster*> & pixelClusterMap,
99 const std::multimap<xAOD::DetectorIdentType, const xAOD::StripCluster*> & stripClusterMap,
100 std::vector<ActsTrk::ATLASUncalibSourceLink>& measurements,
101 const std::vector <FPGATrackSimHit>& hits) const {
102 if (hits.empty()) {
103 ATH_MSG_ERROR("Found FPGATrack without hits");
104 return StatusCode::FAILURE;
105 }
106
107 for (const FPGATrackSimHit& h : hits) {
108 if (h.isReal()) {
109 if (h.isPixel()) {
110 ATH_MSG_DEBUG("Looking for Pixel cluster to match");
111 auto range = pixelClusterMap.equal_range(h.getRdoIdentifier());
112 for (auto it = range.first; it != range.second; ++it) {
113 ATH_CHECK(matchTrackMeasurements<xAOD::PixelCluster>(ctx, *(it->second), h, measurements, pixelContainer));
114 }
115 }
116 else if (h.isStrip()) {
117 ATH_MSG_DEBUG("Looking for Strip cluster to match");
118 auto range = stripClusterMap.equal_range(h.getRdoIdentifier());
119 for (auto it = range.first; it != range.second; ++it) {
120 ATH_CHECK(matchTrackMeasurements<xAOD::StripCluster>(ctx, *(it->second), h, measurements, stripContainer));
121 }
122 }
123 else {
124 ATH_MSG_ERROR("FPGA hit not classified as pixel or strip");
125 return StatusCode::FAILURE;
126 }
127 }
128 else {
129 ATH_MSG_DEBUG("Skipping hit as non-Real");
130 }
131 }
132 return StatusCode::SUCCESS;
133}
134
135
136template <typename XAOD_CLUSTER>
137StatusCode FPGAActsTrkConverter::matchTrackMeasurements(const EventContext& ctx,
138 const XAOD_CLUSTER& cluster,
139 const FPGATrackSimHit & trackHit,
140 std::vector<ActsTrk::ATLASUncalibSourceLink>& measurements,
141 const DataVector<XAOD_CLUSTER>& clusterContainer) const
142{
143 const Identifier::value_type trackId = (trackHit.getHitType() == HitType::spacepoint)
144 ? trackHit.getOriginalHit().getRdoIdentifier()
145 : trackHit.getRdoIdentifier();
146
147 if(cluster.identifier() == trackId) {
148 measurements.emplace_back(ActsTrk::makeATLASUncalibSourceLink(&clusterContainer, cluster.index(), ctx));
149 ATH_MSG_DEBUG("Matched FPGATrackSimHit to xAOD cluster");
150 }
151
152 return StatusCode::SUCCESS;
153}
154
155std::unique_ptr<Acts::BoundTrackParameters> FPGAActsTrkConverter::makeParams (const FPGATrackSimRoad &road) const{
156 using namespace Acts::UnitLiterals;
157
158 std::shared_ptr<const Acts::Surface> actsSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(0., 0., 0.));
159 Acts::BoundVector params;
160
161 constexpr double GeVToMeV = 1000;
162 double d0=0.; //?
163 double z0=0.; //?
164 double phi=road.getX();
165 double eta=0.2;
166 double theta=2*std::atan(std::exp(-eta));
167 double qop = (std::abs(road.getY()) > 1E-9) ? road.getY()/GeVToMeV : 1E-12;
168 double t=0.; //?
169 ATH_MSG_DEBUG("\tphi=" <<phi << " eta=" << eta << " qop=" << qop);
170
171 params << d0, z0, phi, theta, qop, t;
172
173 // Covariance - TODO
174 Acts::BoundSquareMatrix cov = Acts::BoundSquareMatrix::Identity();
175 cov *= (GeVToMeV*GeVToMeV);
176
177 // some ACTS paperwork
179 float mass = Trk::ParticleMasses::mass[hypothesis] * Acts::UnitConstants::MeV;
180 Acts::PdgParticle absPdg = Acts::makeAbsolutePdgParticle(Acts::ePionPlus);
181 Acts::ParticleHypothesis actsHypothesis{
182 absPdg, mass, Acts::AnyCharge{1.0f}};
183
184 return std::make_unique<Acts::BoundTrackParameters>(actsSurface, params,
185 cov, actsHypothesis);
186
187}
188
189
190std::unique_ptr<Acts::BoundTrackParameters> FPGAActsTrkConverter::makeParams (const FPGATrackSimTrack &track) const{
191
192 using namespace Acts::UnitLiterals;
193 std::shared_ptr<const Acts::Surface> actsSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(0., 0., 0.));
194 Acts::BoundVector params;
195
196 constexpr double GeVToMeV = 1000.;
197 double d0=track.getD0();
198 double z0=track.getZ0();
199 double phi=track.getPhi();
200 double eta=track.getEta();
201 double theta=track.getTheta();
202 double qopt=track.getQOverPt()*GeVToMeV;
203 double pt=track.getPt()/GeVToMeV;
204 double px=pt*std::cos(phi);
205 double py=pt*std::sin(phi);
206 double pz=pt*std::sinh(eta);
207 double p=std::sqrt(px*px+py*py+pz*pz);
208 double qop=((p > 1e-10) ? (1/p) : 1e10);
209 if (qopt < 0) qop *= -1;
210 double t=0.;
211
212 params << d0, z0, phi, theta, qop, t;
213 ATH_MSG_DEBUG("\td0= " << d0 << " z0=" <<z0 << " phi=" <<phi << " theta=" << theta<< " qoverp=" << qop);
214
215 // Covariance - let's be honest and say we have no clue ;-)
216 Acts::BoundSquareMatrix cov = Acts::BoundSquareMatrix::Identity();
217
218 (cov)(0,0) *= 0.16; // d0: 0.4 **2 (conservative)
219 (cov)(1,1) *= 25; // z0: 5**2 = 25 (conservative)
220 (cov)(2,2) *= 0.0008; // phi: 0.02**2 = 0.0004, increase a bit = double
221 (cov)(3,3) *= 0.0008; // width in eta is nearly 0.2, but width in theta = 2*atan(e^-eta) will vary. Take biggest one, which is at eta of 0 when width is 0.02, so get 0.02**2 = 0.0004, increase a bit = double
222 (cov)(4,4) *= 0.36; // qop also varies with eta. Error on q/pt conservatively = 0.0003 in mev ^-1, or 0.3 in gev ^-1, double to start giving us 0.6. then square that to get 0.36
223
224
225
226 // some ACTS paperwork
228 float mass = Trk::ParticleMasses::mass[hypothesis] * Acts::UnitConstants::MeV;
229 Acts::PdgParticle absPdg = Acts::makeAbsolutePdgParticle(Acts::ePionPlus);
230 Acts::ParticleHypothesis actsHypothesis{
231 absPdg, mass, Acts::AnyCharge{1.0f}};
232
233 return std::make_unique<Acts::BoundTrackParameters>(actsSurface, params,
234 cov, actsHypothesis);
235
236}
237
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Header file for AthHistogramAlgorithm.
Derived DataVector<T>.
Definition DataVector.h:795
FPGAActsTrkConverter(const std::string &type, const std::string &name, const IInterface *parent)
StatusCode findPrototrackMeasurements(const EventContext &ctx, const xAOD::PixelClusterContainer &pixelClusterContainer, const xAOD::StripClusterContainer &stripClusterContainer, const std::multimap< xAOD::DetectorIdentType, const xAOD::PixelCluster * > &pixelClusterMap, const std::multimap< xAOD::DetectorIdentType, const xAOD::StripCluster * > &stripClusterMap, std::vector< ActsTrk::ATLASUncalibSourceLink > &measurements, const std::vector< FPGATrackSimHit > &hits) const
virtual StatusCode findProtoTracks(const EventContext &ctx, const xAOD::PixelClusterContainer &pixelContainer, const xAOD::StripClusterContainer &stripContainer, std::vector< ActsTrk::ProtoTrack > &foundProtoTracks, const std::vector< std::vector< FPGATrackSimHit > > &hitsInRoads, const std::vector< FPGATrackSimRoad > &roads) const override final
std::unique_ptr< Acts::BoundTrackParameters > makeParams(const FPGATrackSimRoad &road) const
virtual StatusCode initialize() override final
StatusCode matchTrackMeasurements(const EventContext &ctx, const XAOD_CLUSTER &cluster, const FPGATrackSimHit &trackHit, std::vector< ActsTrk::ATLASUncalibSourceLink > &measurements, const DataVector< XAOD_CLUSTER > &clusterContainer) const
Identifier::value_type getRdoIdentifier() const
const FPGATrackSimHit getOriginalHit() const
HitType getHitType() const
float getY() const
float getX() const
ATLASUncalibSourceLink makeATLASUncalibSourceLink(const xAOD::UncalibratedMeasurementContainer *container, std::size_t index, const EventContext &ctx)
constexpr double mass[PARTICLEHYPOTHESES]
the array of masses
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
PixelClusterContainer_v1 PixelClusterContainer
Define the version of the pixel cluster container.
StripCluster_v1 StripCluster
Define the version of the strip cluster class.
StripClusterContainer_v1 StripClusterContainer
Define the version of the strip cluster container.
PixelCluster_v1 PixelCluster
Define the version of the pixel cluster class.