ATLAS Offline Software
Loading...
Searching...
No Matches
FPGAActsTrkConverter Class Reference

#include <FPGAActsTrkConverter.h>

Inheritance diagram for FPGAActsTrkConverter:
Collaboration diagram for FPGAActsTrkConverter:

Public Member Functions

 FPGAActsTrkConverter (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~FPGAActsTrkConverter ()=default
virtual StatusCode initialize () override final
virtual StatusCode findProtoTracks (const EventContext &ctx, const xAOD::PixelClusterContainer &pixelContainer, const xAOD::StripClusterContainer &stripContainer, std::vector< ActsTrk::ProtoTrack > &foundProtoTracks, const FPGATrackSimHitContainer &hitsInRoads, const std::vector< FPGATrackSimRoad > &roads) const override final
virtual StatusCode findProtoTracks (const EventContext &ctx, const xAOD::PixelClusterContainer &pixelContainer, const xAOD::StripClusterContainer &stripContainer, std::vector< ActsTrk::ProtoTrack > &foundProtoTracks, const std::vector< FPGATrackSimTrack > &tracks) const override final

Protected Member Functions

std::unique_ptr< Acts::BoundTrackParameters > makeParams (const FPGATrackSimRoad &road) const
std::unique_ptr< Acts::BoundTrackParameters > makeParams (const FPGATrackSimTrack &track) const
template<typename XAOD_CLUSTER>
StatusCode matchTrackMeasurements (const EventContext &ctx, const XAOD_CLUSTER &cluster, const FPGATrackSimHit &trackHit, std::vector< ActsTrk::ATLASUncalibSourceLink > &measurements, const DataVector< XAOD_CLUSTER > &clusterContainer) const
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 FPGATrackSimHitCollection &hits) const

Private Attributes

const PixelIDm_pixelId {nullptr}
const SCT_IDm_SCTId {nullptr}

Detailed Description

Definition at line 17 of file FPGAActsTrkConverter.h.

Constructor & Destructor Documentation

◆ FPGAActsTrkConverter()

FPGAActsTrkConverter::FPGAActsTrkConverter ( const std::string & type,
const std::string & name,
const IInterface * parent )

Definition at line 13 of file FPGAActsTrkConverter.cxx.

15 : base_class(type,name,parent) { }

◆ ~FPGAActsTrkConverter()

virtual FPGAActsTrkConverter::~FPGAActsTrkConverter ( )
virtualdefault

Member Function Documentation

◆ findPrototrackMeasurements()

StatusCode FPGAActsTrkConverter::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 FPGATrackSimHitCollection & hits ) const
protected

Definition at line 100 of file FPGAActsTrkConverter.cxx.

106 {
107 if (hits.empty()) {
108 ATH_MSG_ERROR("Found FPGATrack without hits");
109 return StatusCode::FAILURE;
110 }
111
112 for (const FPGATrackSimHit* h : hits) {
113 if (h->isReal()) {
114 if (h->isPixel()) {
115 ATH_MSG_DEBUG("Looking for Pixel cluster to match");
116 auto range = pixelClusterMap.equal_range(h->getRdoIdentifier());
117 for (auto it = range.first; it != range.second; ++it) {
118 ATH_CHECK(matchTrackMeasurements<xAOD::PixelCluster>(ctx, *(it->second), *h, measurements, pixelContainer));
119 }
120 }
121 else if (h->isStrip()) {
122 ATH_MSG_DEBUG("Looking for Strip cluster to match");
123 auto range = stripClusterMap.equal_range(h->getRdoIdentifier());
124 for (auto it = range.first; it != range.second; ++it) {
125 ATH_CHECK(matchTrackMeasurements<xAOD::StripCluster>(ctx, *(it->second), *h, measurements, stripContainer));
126 }
127 }
128 else {
129 ATH_MSG_ERROR("FPGA hit not classified as pixel or strip");
130 return StatusCode::FAILURE;
131 }
132 }
133 else {
134 ATH_MSG_DEBUG("Skipping hit as non-Real");
135 }
136 }
137 return StatusCode::SUCCESS;
138}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
StatusCode matchTrackMeasurements(const EventContext &ctx, const XAOD_CLUSTER &cluster, const FPGATrackSimHit &trackHit, std::vector< ActsTrk::ATLASUncalibSourceLink > &measurements, const DataVector< XAOD_CLUSTER > &clusterContainer) const

◆ findProtoTracks() [1/2]

StatusCode FPGAActsTrkConverter::findProtoTracks ( const EventContext & ctx,
const xAOD::PixelClusterContainer & pixelContainer,
const xAOD::StripClusterContainer & stripContainer,
std::vector< ActsTrk::ProtoTrack > & foundProtoTracks,
const FPGATrackSimHitContainer & hitsInRoads,
const std::vector< FPGATrackSimRoad > & roads ) const
finaloverridevirtual

Definition at line 30 of file FPGAActsTrkConverter.cxx.

35 {
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}
#define ATH_MSG_INFO(x)
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 FPGATrackSimHitCollection &hits) const
std::unique_ptr< Acts::BoundTrackParameters > makeParams(const FPGATrackSimRoad &road) const
StripCluster_v1 StripCluster
Define the version of the strip cluster class.
PixelCluster_v1 PixelCluster
Define the version of the pixel cluster class.

◆ findProtoTracks() [2/2]

StatusCode FPGAActsTrkConverter::findProtoTracks ( const EventContext & ctx,
const xAOD::PixelClusterContainer & pixelContainer,
const xAOD::StripClusterContainer & stripContainer,
std::vector< ActsTrk::ProtoTrack > & foundProtoTracks,
const std::vector< FPGATrackSimTrack > & tracks ) const
finaloverridevirtual

Definition at line 64 of file FPGAActsTrkConverter.cxx.

68 {
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 auto hitCollection = std::make_unique<FPGATrackSimHitCollection>();
86 hitCollection->reserve(hits.size());
87 for (const auto& hit : hits) {
88 hitCollection->push_back(new FPGATrackSimHit(hit));
89 }
90 ATH_CHECK(findPrototrackMeasurements(ctx, pixelContainer, stripContainer, pixelClusterMap, stripClusterMap, points, *hitCollection));
91 if (points.size()) {
92 ATH_MSG_DEBUG("\tMaking a proto-track with " << points.size() << " clusters");
93 std::unique_ptr<Acts::BoundTrackParameters> inputPerigee = makeParams(track);
94 foundProtoTracks.emplace_back(points, std::move(inputPerigee));
95 }
96 }
97 return StatusCode::SUCCESS;
98}

◆ initialize()

StatusCode FPGAActsTrkConverter::initialize ( )
finaloverridevirtual

Definition at line 18 of file FPGAActsTrkConverter.cxx.

18 {
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}
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ makeParams() [1/2]

std::unique_ptr< Acts::BoundTrackParameters > FPGAActsTrkConverter::makeParams ( const FPGATrackSimRoad & road) const
protected

Definition at line 160 of file FPGAActsTrkConverter.cxx.

160 {
161 using namespace Acts::UnitLiterals;
162
163 std::shared_ptr<const Acts::Surface> actsSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(0., 0., 0.));
164 Acts::BoundVector params;
165
166 constexpr double GeVToMeV = 1000;
167 double d0=0.; //?
168 double z0=0.; //?
169 double phi=road.getX();
170 double eta=0.2;
171 double theta=2*std::atan(std::exp(-eta));
172 double qop = (std::abs(road.getY()) > 1E-9) ? road.getY()/GeVToMeV : 1E-12;
173 double t=0.; //?
174 ATH_MSG_DEBUG("\tphi=" <<phi << " eta=" << eta << " qop=" << qop);
175
176 params << d0, z0, phi, theta, qop, t;
177
178 // Covariance - TODO
179 Acts::BoundSquareMatrix cov = Acts::BoundSquareMatrix::Identity();
180 cov *= (GeVToMeV*GeVToMeV);
181
182 // some ACTS paperwork
184 float mass = Trk::ParticleMasses::mass[hypothesis] * Acts::UnitConstants::MeV;
185 Acts::PdgParticle absPdg = Acts::makeAbsolutePdgParticle(Acts::ePionPlus);
186 Acts::ParticleHypothesis actsHypothesis{
187 absPdg, mass, Acts::AnyCharge{1.0f}};
188
189 return std::make_unique<Acts::BoundTrackParameters>(actsSurface, params,
190 cov, actsHypothesis);
191
192}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
float getY() const
float getX() const
constexpr double mass[PARTICLEHYPOTHESES]
the array of masses
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.

◆ makeParams() [2/2]

std::unique_ptr< Acts::BoundTrackParameters > FPGAActsTrkConverter::makeParams ( const FPGATrackSimTrack & track) const
protected

Definition at line 195 of file FPGAActsTrkConverter.cxx.

195 {
196
197 using namespace Acts::UnitLiterals;
198 std::shared_ptr<const Acts::Surface> actsSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(0., 0., 0.));
199 Acts::BoundVector params;
200
201 constexpr double GeVToMeV = 1000.;
202 double d0=track.getD0();
203 double z0=track.getZ0();
204 double phi=track.getPhi();
205 double eta=track.getEta();
206 double theta=track.getTheta();
207 double qopt=track.getQOverPt()*GeVToMeV;
208 double pt=track.getPt()/GeVToMeV;
209 double px=pt*std::cos(phi);
210 double py=pt*std::sin(phi);
211 double pz=pt*std::sinh(eta);
212 double p=std::sqrt(px*px+py*py+pz*pz);
213 double qop=((p > 1e-10) ? (1/p) : 1e10);
214 if (qopt < 0) qop *= -1;
215 double t=0.;
216
217 params << d0, z0, phi, theta, qop, t;
218 ATH_MSG_DEBUG("\td0= " << d0 << " z0=" <<z0 << " phi=" <<phi << " theta=" << theta<< " qoverp=" << qop);
219
220 // Covariance - let's be honest and say we have no clue ;-)
221 Acts::BoundSquareMatrix cov = Acts::BoundSquareMatrix::Identity();
222
223 (cov)(0,0) *= 0.16; // d0: 0.4 **2 (conservative)
224 (cov)(1,1) *= 25; // z0: 5**2 = 25 (conservative)
225 (cov)(2,2) *= 0.0008; // phi: 0.02**2 = 0.0004, increase a bit = double
226 (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
227 (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
228
229
230
231 // some ACTS paperwork
233 float mass = Trk::ParticleMasses::mass[hypothesis] * Acts::UnitConstants::MeV;
234 Acts::PdgParticle absPdg = Acts::makeAbsolutePdgParticle(Acts::ePionPlus);
235 Acts::ParticleHypothesis actsHypothesis{
236 absPdg, mass, Acts::AnyCharge{1.0f}};
237
238 return std::make_unique<Acts::BoundTrackParameters>(actsSurface, params,
239 cov, actsHypothesis);
240
241}

◆ matchTrackMeasurements()

template<typename XAOD_CLUSTER>
StatusCode FPGAActsTrkConverter::matchTrackMeasurements ( const EventContext & ctx,
const XAOD_CLUSTER & cluster,
const FPGATrackSimHit & trackHit,
std::vector< ActsTrk::ATLASUncalibSourceLink > & measurements,
const DataVector< XAOD_CLUSTER > & clusterContainer ) const
protected

Definition at line 142 of file FPGAActsTrkConverter.cxx.

147{
148 const Identifier::value_type trackId = (trackHit.getHitType() == HitType::spacepoint)
149 ? trackHit.getOriginalHit().getRdoIdentifier()
150 : trackHit.getRdoIdentifier();
151
152 if(cluster.identifier() == trackId) {
153 measurements.emplace_back(ActsTrk::makeATLASUncalibSourceLink(&clusterContainer, cluster.index(), ctx));
154 ATH_MSG_DEBUG("Matched FPGATrackSimHit to xAOD cluster");
155 }
156
157 return StatusCode::SUCCESS;
158}
Identifier::value_type getRdoIdentifier() const
const FPGATrackSimHit getOriginalHit() const
HitType getHitType() const
ATLASUncalibSourceLink makeATLASUncalibSourceLink(const xAOD::UncalibratedMeasurementContainer *container, std::size_t index, const EventContext &ctx)

Member Data Documentation

◆ m_pixelId

const PixelID* FPGAActsTrkConverter::m_pixelId {nullptr}
private

Definition at line 55 of file FPGAActsTrkConverter.h.

55{nullptr};

◆ m_SCTId

const SCT_ID* FPGAActsTrkConverter::m_SCTId {nullptr}
private

Definition at line 56 of file FPGAActsTrkConverter.h.

56{nullptr};

The documentation for this class was generated from the following files: