ATLAS Offline Software
Loading...
Searching...
No Matches
TruthTrackBuilder.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
6// TruthTrackBuilder.cxx, (c) ATLAS Detector software
8
9// package include
10#include "TruthTrackBuilder.h"
11// Trk includes
12#include "TrkTrack/Track.h"
13
17// Gaudi
18#include "GaudiKernel/IPartPropSvc.h"
19#include "GaudiKernel/SystemOfUnits.h"
20// HepMC
24#include "HepPDT/ParticleDataTable.hh"
25
27
29Trk::TruthTrackBuilder::TruthTrackBuilder(const std::string& t, const std::string& n, const IInterface* p) :
30 AthAlgTool(t,n,p),
31 m_DetID(nullptr),
32 m_particlePropSvc("PartPropSvc",n),
33 m_particleDataTable(nullptr)
34{
35 declareInterface<Trk::ITruthTrackBuilder>(this);
36}
37
38
39// Athena algtool's Hooks
41{
42 ATH_MSG_VERBOSE("Initializing ...");
43 if (m_rotcreator.retrieve().isFailure()) {
44 msg(MSG::FATAL) << "Could not get " << m_rotcreator.type() << endmsg;
45 return StatusCode::FAILURE;
46 }
47
48 if (m_rotcreatorbroad.retrieve().isFailure()) {
49 msg(MSG::FATAL) << "Could not get " << m_rotcreatorbroad.type() << endmsg;
50 return StatusCode::FAILURE;
51 }
52
53 // track fitter
54 if (m_trackFitter.retrieve().isFailure()){
55 ATH_MSG_ERROR("Could not retrieve " << m_trackFitter << ". Aborting ...");
56 return StatusCode::FAILURE;
57 }
58
59 // extrapolator
60 if (m_extrapolator.retrieve().isFailure()){
61 ATH_MSG_ERROR("Could not retrieve " << m_extrapolator << ". Aborting ...");
62 return StatusCode::FAILURE;
63 }
64
65 // particle property service
66 if (m_particlePropSvc.retrieve().isFailure())
67 {
68 ATH_MSG_ERROR("Can not retrieve " << m_particlePropSvc << " . Aborting ... " );
69 return StatusCode::FAILURE;
70 }
71
72 // and the particle data table
74 if (m_particleDataTable==nullptr)
75 {
76 ATH_MSG_ERROR( "Could not get ParticleDataTable! Cannot associate pdg code with charge. Aborting. " );
77 return StatusCode::FAILURE;
78 }
79
80 // need an Atlas id-helper to identify sub-detectors, take the one from detStore
81 if (detStore()->retrieve(m_DetID, "AtlasID").isFailure()) {
82 ATH_MSG_ERROR ("Could not get AtlasDetectorID helper" );
83 return StatusCode::FAILURE;
84 }
85 return StatusCode::SUCCESS;
86}
87
88
90{
91 ATH_MSG_VERBOSE("Finalizing ...");
92 return StatusCode::SUCCESS;
93}
94
95
97{
98 const EventContext& ctx = Gaudi::Hive::currentContext();
99 if( segs != nullptr ){
100 ATH_MSG_WARNING("Requested to fill segment collection but mode not implemented");
101 }
102 ATH_MSG_VERBOSE("The PRD Truth trajectory contains " << prdTraj.prds.size() << " PRDs.");
103
104 // get the associated GenParticle
105 auto genPart = prdTraj.genParticle;
106 if (!genPart) {
107 ATH_MSG_WARNING("No GenParticle associated to this PRD_TruthTrajectory. Ignoring track creation.");
108 return nullptr;
109 }
110 ATH_MSG_VERBOSE("Got PRD Truth trajectory with " << prdTraj.nDoF << " degrees of freedom.");
111 // check min degrees of freedom
112 if ( m_minNdof > 0 && prdTraj.nDoF < m_minNdof) return nullptr;
113 // get the startPosition : math library madness as usual
114 Amg::Vector3D startPos = genPart->production_vertex() ?
115 Amg::Vector3D(genPart->production_vertex()->position().x(),
116 genPart->production_vertex()->position().y(),
117 genPart->production_vertex()->position().z()) : Amg::Vector3D(0.,0.,0.);
118 Amg::Vector3D startMom(genPart->momentum().x(),
119 genPart->momentum().y(),
120 genPart->momentum().z());
122
123 int pdgCode = genPart->pdg_id();
124 // get the charge: ap->charge() is used later, DOES NOT WORK RIGHT NOW
125 const HepPDT::ParticleData* ap = m_particleDataTable->particle(std::abs(pdgCode));
126 double charge = 1.;
127 if (ap) charge = ap->charge();
128 // since the PDT table only has abs(PID) values for the charge
129 charge *= (pdgCode > 0.) ? 1. : -1.;
130
131 // ----------------------- get teh PRDS and start
132 const std::vector<const Trk::PrepRawData*> & clusters = prdTraj.prds;
133 // nominal 0,0,0 position for track fit seeding
134 Trk::PerigeeSurface persurf;
135 Trk::CurvilinearParameters startParams(startPos,startMom,charge);
136 //minimal conversion; ideally the extrapolator would return a unique_ptr
137 auto per = std::unique_ptr<Trk::TrackParameters>(
138 m_extrapolator->extrapolate(ctx,
139 startParams,
140 persurf,
142 false,
144 if (!per) {
145 ATH_MSG_DEBUG("Perigee creation for genParticle start position failed. "
146 "Skipping track creation.");
147 return nullptr;
148 }
149 // first TrackStateOnSurface is the Perigee
150 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> typePattern;
151 typePattern.set(Trk::TrackStateOnSurface::Perigee);
152
153 const Trk::TrackStateOnSurface *pertsos=new Trk::TrackStateOnSurface(nullptr,std::move(per),nullptr,typePattern);
154 auto traj = std::make_unique<Trk::TrackStates>();
155 traj->push_back(pertsos);
156
157
158
159 std::unique_ptr<const Trk::TrackParameters> prevpar(startParams.uniqueClone());
160 // First create a Trk::Track object 'traj' that will go into the fitter for refitting
161 int i=0;
162 for ( ;i<(int)clusters.size();i++){
163 if (m_DetID->is_trt(clusters[i]->identify())) break;
164 const Trk::Surface &surf=clusters[i]->detectorElement()->surface(clusters[i]->identify());
165 if (surf==prevpar->associatedSurface()) continue;
166 bool ispixel=false;
167 if (m_DetID->is_pixel(clusters[i]->identify())) ispixel=true;
168
169 auto thispar = std::unique_ptr<const Trk::TrackParameters>(
170 m_extrapolator->extrapolate(ctx,
171 *prevpar,
172 surf,
174 false,
176 if (!thispar)
177 break;
178 if (!surf.insideBounds(thispar->localPosition(),20*Gaudi::Units::mm,50*Gaudi::Units::mm)) {
179 continue;
180 }
181 AmgVector(5) params=thispar->parameters();
182 params[Trk::loc1]=clusters[i]->localPosition().x();
183 if (ispixel) params[Trk::loc2]=clusters[i]->localPosition().y();
184 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> typePattern;
186 std::unique_ptr<Trk::RIO_OnTrack> rot{m_rotcreator->correct(*clusters[i],*thispar,ctx)};
187 if (!rot) {
188 continue;
189 }
190 // create the ROTs for the reference trajectory
191 const Trk::TrackStateOnSurface *tsos=new Trk::TrackStateOnSurface(std::move(rot),thispar->uniqueClone(),nullptr,typePattern);
192 traj->push_back(tsos);
193 prevpar=std::move(thispar);
194 }
195 // this is the reference trajectory to be refitted
196 Trk::TrackInfo info;
197 Trk::Track track(info, std::move(traj), nullptr);
198 if (/* ndof<0 */ (track.measurementsOnTrack()->size() < m_minSiHits &&
199 std::abs(genPart->momentum().eta()) <= m_forwardBoundary) ||
200 (track.measurementsOnTrack()->size() < m_minSiHitsForward &&
201 std::abs(genPart->momentum().eta()) > m_forwardBoundary) ||
204 "Track does not fulfill requirements for refitting. Skipping it.");
205 return nullptr;
206 }
207 // choose the material effects
209
210 // create the refitted track
211 //
212 //
213 //
215
216 Trk::Track *refittedtrack=(m_trackFitter->fit(Gaudi::Hive::currentContext(),track,false,materialInteractions)).release();
217
219 Trk::Track *refittedtrack2=nullptr;
220 if (refittedtrack && (int)clusters.size()-i>=9){
221
222 //owner of measurements so they get cleaned up automatically if needs be
223 std::vector<std::unique_ptr<MeasurementBase>> meassetOwn;
224 //vector of plain ptr to the owned ones acts as a "view"
225 Trk::MeasurementSet measset;
226
227 std::unique_ptr<const Trk::TrackParameters> prevpar(refittedtrack->trackParameters()->back()->uniqueClone());
228 for (;i<(int)clusters.size();i++) {
229 const Trk::Surface *surf=&clusters[i]->detectorElement()->surface(clusters[i]->identify());
230 std::unique_ptr<const Trk::TrackParameters> thispar(m_extrapolator->extrapolate(ctx,*prevpar,*surf,Trk::alongMomentum,false,Trk::nonInteracting));
231 if (!thispar) break;
232 Trk::RIO_OnTrack *rot=m_rotcreatorbroad->correct(*clusters[i],*thispar,ctx);
233
234 if (rot) {
235 meassetOwn.emplace_back(rot);
236 measset.push_back(meassetOwn.back().get());
237 }
238 prevpar=std::move(thispar);
239 }
240
241 refittedtrack2=(m_trackFitter->fit(Gaudi::Hive::currentContext(),*refittedtrack,measset,false,materialInteractions)).release();
242
243 if (!refittedtrack2){
244 auto traj2 = std::make_unique<Trk::TrackStates>();
245 for (const Trk::TrackStateOnSurface *j : *refittedtrack->trackStateOnSurfaces()) traj2->push_back(new Trk::TrackStateOnSurface(*j));
246 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> typePattern2;
247 typePattern2.set(Trk::TrackStateOnSurface::Outlier);
248 //measset needs to be unique_ptr before progress further
249 for (auto & j : meassetOwn) {
250 traj2->push_back(new Trk::TrackStateOnSurface(
251 std::move(j),
252 nullptr,
253 nullptr,
254 typePattern2));
255 }
256 refittedtrack2 = new Trk::Track(refittedtrack->info(),
257 std::move(traj2),
258 refittedtrack->fitQuality()->uniqueClone());
259 }
260 } else if(!refittedtrack){
261 ATH_MSG_VERBOSE("Track fit of truth trajectory NOT successful, NO track created. ");
262 return nullptr;
263 }
264
265 if (refittedtrack2) delete refittedtrack;
266 if (!refittedtrack2 && refittedtrack) refittedtrack2=refittedtrack;
267
268 //
269 ATH_MSG_DEBUG("Track fit of truth trajectory successful, track created. ");
270 // return what you have
271 // Before returning, fix the creator
272 if (refittedtrack2){
274 }
275 return refittedtrack2;
276}
277
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
double charge(const T &p)
Definition AtlasPID.h:997
#define AmgVector(rows)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
std::unique_ptr< FitQuality > uniqueClone() const
NVI uniqueClone.
std::unique_ptr< ParametersBase< DIM, T > > uniqueClone() const
clone method for polymorphic deep copy returning unique_ptr; it is not overriden, but uses the existi...
Class describing the Line to which the Perigee refers to.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Abstract Base Class for tracking surfaces.
virtual bool insideBounds(const Amg::Vector2D &locpos, double tol1=0., double tol2=0.) const =0
virtual methods to be overwritten by the inherited surfaces
Contains information about the 'fitter' of this track.
void setPatternRecognitionInfo(const TrackPatternRecoInfo &patternReco)
Method setting the pattern recognition algorithm.
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
const DataVector< const TrackParameters > * trackParameters() const
Return a pointer to a vector of TrackParameters.
const TrackInfo & info() const
Returns a const ref to info of a const tracks.
const FitQuality * fitQuality() const
return a pointer to the fit quality const-overload
ToolHandle< IRIO_OnTrackCreator > m_rotcreator
Gaudi::Property< unsigned int > m_minSiHits
min number of Si hits for refit
Gaudi::Property< float > m_forwardBoundary
Boundary eta value defining the forward region.
Gaudi::Property< unsigned int > m_minSiHitsForward
min number of Si hits for refit in forward region (ITk specific)
const AtlasDetectorID * m_DetID
TruthTrackBuilder(const std::string &t, const std::string &n, const IInterface *p)
Constructor.
Track * createTrack(const PRD_TruthTrajectory &prdTraj, SegmentCollection *segs=0) const
return a map of GenParticles to PRDs for further processing
ServiceHandle< IPartPropSvc > m_particlePropSvc
Pointer to the particle properties svc *‍/.
Gaudi::Property< size_t > m_minNdof
checks min degrees of freedom if bigger -1
Gaudi::Property< bool > m_onlyPrimaries
restrict track creation to primaries
const HepPDT::ParticleDataTable * m_particleDataTable
ParticleDataTable needed to get connection pdg_code <-> charge *‍/.
ToolHandle< ITrackFitter > m_trackFitter
fits the PRDs
Gaudi::Property< int > m_matEffects
ToolHandle< IRIO_OnTrackCreator > m_rotcreatorbroad
ToolHandle< IExtrapolator > m_extrapolator
extrapolator
Eigen::Matrix< double, 3, 1 > Vector3D
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
@ alongMomentum
@ anyDirection
std::vector< const MeasurementBase * > MeasurementSet
vector of fittable measurements
Definition FitterTypes.h:30
DataVector< Trk::Segment > SegmentCollection
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ loc1
Definition ParamDefs.h:34
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
simple definitiion of a PRD_TruhtTrajectory
std::vector< const Trk::PrepRawData * > prds
public members
HepMC::ConstGenParticlePtr genParticle
MsgStream & msg
Definition testRead.cxx:32