ATLAS Offline Software
Loading...
Searching...
No Matches
JiveXML::TrackRetrieverHelpers Namespace Reference

Namespace for all helper functions. More...

Functions

const Trk::PerigeegetPerigeeParameters (const Trk::Track *track, DataVect &pt, DataVect &d0, DataVect &z0, DataVect &phi0, DataVect &cotTheta, DataVect &covMatrix)
 Obtain the perigee parameters for a given track, if available, and fill them in the corresponding data vectors.
std::vector< const Trk::TrackStateOnSurface * > getTrackStateOnSurfaces (const Trk::Track *track, const Trk::Perigee *perigee, bool doHitsSorting)
 Get a list of track-State on Surfaces for measurement and outlier hits, sorted using the perigee comparison functions.
void getPolylineFromHits (const std::vector< const Trk::TrackStateOnSurface * > &TSoSVec, DataVect &polylineX, DataVect &polylineY, DataVect &polylineZ, DataVect &numPolyline)
 Get polyline hits if available.
const Trk::RIO_OnTrackgetBaseInfoFromHit (const Trk::TrackStateOnSurface *tsos, const AtlasDetectorID *idHelper, DataVect &isOutlier, DataVect &hits, DataVect &driftSign, DataVect &tsosDetType)
 Retrieve all the basic hit information from the Trk::TrackStateOnSurface.
void getResidualPullFromHit (const Trk::TrackStateOnSurface *tsos, const Trk::RIO_OnTrack *rot, const ToolHandle< Trk::IResidualPullCalculator > &residualPullCalculator, DataVect &tsosResLoc1, DataVect &tsosResLoc2, DataVect &tsosPullLoc1, DataVect &tsosPullLoc2)
 Get the residual pull information from the Trk::TrackStateOnSurface hit.
void getTruthFromTrack (const Trk::Track *track, const TrackCollection *trackCollection, SG::ReadHandle< TrackTruthCollection > &truthCollection, DataVect &barcode)
 Get the barcode of the associated truth track.

Detailed Description

Namespace for all helper functions.

Function Documentation

◆ getBaseInfoFromHit()

const Trk::RIO_OnTrack * JiveXML::TrackRetrieverHelpers::getBaseInfoFromHit ( const Trk::TrackStateOnSurface * tsos,
const AtlasDetectorID * idHelper,
DataVect & isOutlier,
DataVect & hits,
DataVect & driftSign,
DataVect & tsosDetType )

Retrieve all the basic hit information from the Trk::TrackStateOnSurface.

Returns
the reconstruction input object (RIO) for further use

Definition at line 181 of file TrackRetriever.cxx.

182 {
183
184 //Get corresponding measurement
185 const Trk::MeasurementBase *measurement = tsos->measurementOnTrack();
186
187 //If measurement is invalid, return imediately
188 if (!measurement) return nullptr;
189
190 // get RIO_OnTrack
191 const Trk::RIO_OnTrack* rot = dynamic_cast<const Trk::RIO_OnTrack*>(measurement);
192
193 // check for competing RIO_OnTracks
194 if (!rot) {
195 // try to get identifier by CompetingROT:
196 const Trk::CompetingRIOsOnTrack* comprot = dynamic_cast<const Trk::CompetingRIOsOnTrack*>(measurement);
197 //Get the input object with highest probability
198 if (comprot) rot = &(comprot->rioOnTrack(comprot->indexOfMaxAssignProb()));
199 }
200
201 //If there is still no RIO_onTrack, return Null
202 if (!rot) return nullptr;
203
204 // Now start writing out values:
205 // Check if this is an outlier hit
206 isOutlier.emplace_back(tsos->type(Trk::TrackStateOnSurface::Outlier));
207
208 //Now try to get the identifier, create an empty invalid one if no rot
209 Identifier hitId (rot->identify());
210 //Check if it is valid, othwerise store 0
211 hits.push_back( hitId.is_valid()?DataType( hitId.get_compact() ):DataType(0) );
212
213 // get sign of drift radius for TRT measurements
214 int theDriftSign = 0;
215 if (idHelper->is_trt(hitId)) {
216 // get local parameters
217 theDriftSign = measurement->localParameters()[Trk::driftRadius] > 0. ? 1 : -1;
218 }
219 driftSign.emplace_back(theDriftSign);
220
221 //Now get the detector type of the hit
222 if ( !hitId.is_valid() ){
223 tsosDetType.emplace_back("unident");
224 } else if (idHelper->is_pixel(hitId) ) {
225 tsosDetType.emplace_back("PIX"); // is PIX in Atlantis
226 } else if (idHelper->is_sct(hitId)) {
227 tsosDetType.emplace_back("SIL"); // is SIL in Atlantis
228 } else if (idHelper->is_trt(hitId)) {
229 tsosDetType.emplace_back("TRT");
230 } else if (idHelper->is_mdt(hitId)) {
231 tsosDetType.emplace_back("MDT");
232 } else if (idHelper->is_csc(hitId)) {
233 tsosDetType.emplace_back("CSC");
234 } else if (idHelper->is_rpc(hitId)) {
235 tsosDetType.emplace_back("RPC");
236 } else if (idHelper->is_tgc(hitId)) {
237 tsosDetType.emplace_back("TGC");
238 } else {
239 tsosDetType.emplace_back("unident");
240 }
241
242 //Return the reco input object
243 return rot;
244 }
bool is_mdt(Identifier id) const
bool is_rpc(Identifier id) const
bool is_sct(Identifier id) const
bool is_tgc(Identifier id) const
bool is_pixel(Identifier id) const
bool is_csc(Identifier id) const
bool is_trt(Identifier id) const
Templated class to convert any object that is streamable in a ostringstream in a string.
Definition DataType.h:21
Base class for all CompetingRIOsOnTack implementations, extends the common MeasurementBase.
unsigned int indexOfMaxAssignProb() const
Index of the ROT with the highest assignment probability.
virtual const RIO_OnTrack & rioOnTrack(unsigned int) const =0
returns the RIO_OnTrack (also known as ROT) objects depending on the integer.
This class is the pure abstract base class for all fittable tracking measurements.
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Identifier identify() const
return the identifier -extends MeasurementBase
const MeasurementBase * measurementOnTrack() const
returns MeasurementBase const overload
bool type(const TrackStateOnSurfaceType type) const
Use this method to find out if the TSoS is of a certain type: i.e.
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
@ driftRadius
trt, straws
Definition ParamDefs.h:53

◆ getPerigeeParameters()

const Trk::Perigee * JiveXML::TrackRetrieverHelpers::getPerigeeParameters ( const Trk::Track * track,
DataVect & pt,
DataVect & d0,
DataVect & z0,
DataVect & phi0,
DataVect & cotTheta,
DataVect & covMatrix )

Obtain the perigee parameters for a given track, if available, and fill them in the corresponding data vectors.

Perigee parameters are written out in the old format using \( \cot\theta \) and \( q/p_T \)

Returns
the perigee parameter object for further use

Get perigee parameters in old format ( \( d_0 \), \( z_0 \), \( \phi \), \( \cot\theta \), \( q/p_T \)), whereas tracking uses ( \( d_0 \), \( z_0 \), \( \phi \), \( \theta \), q/p), therefore a transformation of the covariance matrix is needed

get transformed covariance matrix

Definition at line 48 of file TrackRetriever.cxx.

50 {
51
57 const Trk::Perigee *perigee = track->perigeeParameters();
58
59 //return immediately if there is no perigee information
60 if (!perigee) return nullptr;
61
62 //write out p_T
63 if ((perigee->parameters())[Trk::qOverP]==0) pt.emplace_back(9999.);
64 else pt.push_back( (perigee->charge() > 0) ? DataType(perigee->pT()/Gaudi::Units::GeV) : DataType((-perigee->pT())/Gaudi::Units::GeV));
65
66 d0.emplace_back((perigee->parameters())[Trk::d0]/Gaudi::Units::cm);
67 z0.emplace_back(perigee->parameters()[Trk::z0]/Gaudi::Units::cm);
68 phi0.emplace_back(perigee->parameters()[Trk::phi0]);
69
70 if (perigee->parameters()[Trk::theta] == 0.) cotTheta.emplace_back(9999.);
71 else cotTheta.emplace_back(1./tan(perigee->parameters()[Trk::theta]));
72
73 // CLHEP->Eigen migration. jpt Dec'13
74 // https://twiki.cern.ch/twiki/bin/viewauth/Atlas/MigrationCLHEPtoEigen
75 // https://twiki.cern.ch/twiki/bin/viewauth/Atlas/MigrationToUpdatedEDM#Changes_to_TrkParameters
76
78 AmgSymMatrix(5) covVert;
79
80 const AmgSymMatrix(5)* covariance = perigee->covariance(); //perigee cannot be null here
81 if (perigee && covariance) {
82 // do trafo to old format
83 double measuredTheta = perigee->parameters()[Trk::theta];
84 double measuredQoverp = perigee->parameters()[Trk::qOverP];
85 const Trk::JacobianThetaPToCotThetaPt theJac( measuredTheta, measuredQoverp );
86 covVert = covariance->similarity(theJac);
87 }else{
88 for ( int ii=0; ii<20; ii++){ // placeholder. Do this nicer.
89 covVert(ii) = 0.;
90 }
91 }
92 //Scale covariance matrix values to get good resolution with fixed
93 //precision in JiveXML data
94
95 const long scale = 10000;
96 const double thisScale(scale/100.);
97 // Migration: Now only has diagonal elements from covariance matrix ?
98 covMatrix.emplace_back(covVert(0)*thisScale); // 5 elements
99 covMatrix.emplace_back(covVert(1)*thisScale);
100 covMatrix.emplace_back(covVert(2)*thisScale);
101 covMatrix.emplace_back(covVert(3)*thisScale);
102 covMatrix.emplace_back(covVert(4)*thisScale);
103
104 // Used to be 15 elements before migration, so need to put 10 placeholders
105 for ( int i=0; i<10; i++){
106 covMatrix.emplace_back( 0. );
107 }
108
109 //All for perigee, return object for use by other functions
110 return perigee;
111 }
#define AmgSymMatrix(dim)
if(febId1==febId2)
This is the 5x5 jacobian for the transformation of track parameters and errors having the new standar...
double charge() const
Returns the charge.
double pT() const
Access method for transverse momentum.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64

◆ getPolylineFromHits()

void JiveXML::TrackRetrieverHelpers::getPolylineFromHits ( const std::vector< const Trk::TrackStateOnSurface * > & TSoSVec,
DataVect & polylineX,
DataVect & polylineY,
DataVect & polylineZ,
DataVect & numPolyline )

Get polyline hits if available.

Polyline tracks that have less than 2 points are not useful - skip

Definition at line 155 of file TrackRetriever.cxx.

156 {
157 int numPoly = 0 ;
158 if (TSoSVec.size() > 1) {
159 //Loop over track state on surfaces
160 std::vector<const Trk::TrackStateOnSurface*>::const_iterator tsosIter;
161 const double onetenth(0.1);
162 for (tsosIter=TSoSVec.begin(); tsosIter!=TSoSVec.end(); ++tsosIter) {
163 // get global track position
164 if (!(*tsosIter)->trackParameters()) continue ;
165 const Amg::Vector3D& pos = (*tsosIter)->trackParameters()->position();
166 polylineX.emplace_back(pos.x()*onetenth);
167 polylineY.emplace_back(pos.y()*onetenth);
168 polylineZ.emplace_back(pos.z()*onetenth);
169 ++numPoly;
170 }
171 }
172 //Store counter as well
173 numPolyline.emplace_back(numPoly);
174 }
Eigen::Matrix< double, 3, 1 > Vector3D

◆ getResidualPullFromHit()

void JiveXML::TrackRetrieverHelpers::getResidualPullFromHit ( const Trk::TrackStateOnSurface * tsos,
const Trk::RIO_OnTrack * rot,
const ToolHandle< Trk::IResidualPullCalculator > & residualPullCalculator,
DataVect & tsosResLoc1,
DataVect & tsosResLoc2,
DataVect & tsosPullLoc1,
DataVect & tsosPullLoc2 )

Get the residual pull information from the Trk::TrackStateOnSurface hit.

Using track residual tool: ResidualPullCalculator -excerpt from Tracking/TrkValidation/TrkValTools/src/BasicValidationNtupleTool.cxx

Definition at line 250 of file TrackRetriever.cxx.

252 {
253
254 //Define default return values for invalid states
255 double ResLoc1 = -99.;
256 double ResLoc2 = -99.;
257 double PullLoc1 = -99.;
258 double PullLoc2 = -99.;
259
260 // get TrackParameters on the surface to calculate residual
261 const Trk::TrackParameters* tsosParameters = tsos->trackParameters();
262
263 //Check we got the parameters
264 if (tsosParameters){
265
270
271 //Get the residualPull object
272 std::optional<Trk::ResidualPull> residualPull = residualPullCalculator->residualPull(rot, tsosParameters,Trk::ResidualPull::Biased);
273
274 if (residualPull) {
275 //Get the first residual
276 ResLoc1 = residualPull->residual()[Trk::loc1];
277 //Get the second residual for more than 1 dimension
278 if (residualPull->dimension() >= 2) ResLoc2 = residualPull->residual()[Trk::loc2];
279
280 if ((residualPull->isPullValid()) ) {
281 //Get the first residual
282 PullLoc1 = residualPull->pull()[Trk::loc1];
283 //Get the second residual for more than 1 dimension
284 if (residualPull->dimension() >= 2) PullLoc2 = residualPull->pull()[Trk::loc2];
285 }
286 } // end if residualPull
287 } // end if tsosParameters
288
289 //Finally store the values
290 tsosResLoc1.emplace_back(ResLoc1 );
291 tsosResLoc2.emplace_back(ResLoc2 );
292 tsosPullLoc1.emplace_back(PullLoc1 );
293 tsosPullLoc2.emplace_back(PullLoc2 );
294 }
@ Biased
RP with track state including the hit.
const TrackParameters * trackParameters() const
return ptr to trackparameters const overload
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ loc1
Definition ParamDefs.h:34
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ getTrackStateOnSurfaces()

std::vector< const Trk::TrackStateOnSurface * > JiveXML::TrackRetrieverHelpers::getTrackStateOnSurfaces ( const Trk::Track * track,
const Trk::Perigee * perigee,
bool doHitsSorting )

Get a list of track-State on Surfaces for measurement and outlier hits, sorted using the perigee comparison functions.

Returns
a std::vector of Trk::TrackStateOnSurface*

Definition at line 118 of file TrackRetriever.cxx.

118 {
119 // vector for the return object
120 std::vector<const Trk::TrackStateOnSurface*> TSoSVec;
121
122 // loop over TrackStateOnSurfaces to extract interesting ones
123 Trk::TrackStates::const_iterator tsos = track->trackStateOnSurfaces()->begin();
124 for (; tsos!=track->trackStateOnSurfaces()->end(); ++tsos) {
125 // include measurements AND outliers:
126 if ((*tsos)->type(Trk::TrackStateOnSurface::Measurement) || (*tsos)->type(Trk::TrackStateOnSurface::Outlier) ) {
127 // add to temp vector
128 TSoSVec.push_back(*tsos);
129 } // end if TSoS is measurement or outlier
130 } // end loop over TSoS
131
132 // sort the selected TSoS, if not already sorted using a comparison functor
133
134 if (perigee) {
135 //Get hold of the comparison functor
138 if (compFunc){
139 if (doHitsSorting) {
140 //Sort track state on surface if needed
141 if (TSoSVec.size() > 2 && !is_sorted(TSoSVec.begin(), TSoSVec.end(), *compFunc))
142 std::sort(TSoSVec.begin(), TSoSVec.end(), *compFunc);
143 }
144 }
145 delete compFunc;
146 } // end if compFunc
147
148 //Now return that vector
149 return TSoSVec;
150 }
bool(*)(const xAOD::Jet *, const xAOD::Jet *) compFunc
Definition JetSorter.cxx:12
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
Class providing comparison function, or relational definition, for sorting MeasurementBase objects.
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.

◆ getTruthFromTrack()

void JiveXML::TrackRetrieverHelpers::getTruthFromTrack ( const Trk::Track * track,
const TrackCollection * trackCollection,
SG::ReadHandle< TrackTruthCollection > & truthCollection,
DataVect & barcode )

Get the barcode of the associated truth track.

Definition at line 299 of file TrackRetriever.cxx.

300 {
301
302 if (!truthCollection.isValid()) {
303 //Fill with zero if none found
304 barcode.emplace_back(0);
305 return;
306 }
307
308 //Get the element link to this track collection
310 tracklink.setElement(track);
311 tracklink.setStorableObject(*trackCollection);
312
313 //try to find it in the truth collection
314 std::map<Trk::TrackTruthKey,TrackTruth>::const_iterator tempTrackTruthItr = truthCollection->find(tracklink);
315
316 //Fill with zero if none found
317 if (tempTrackTruthItr == truthCollection->end()){
318 //Fill with zero if none found
319 barcode.emplace_back(0);
320 return;
321 }
322
323 //if found, Store the barcode
324 barcode.emplace_back((*tempTrackTruthItr).second.particleLink().barcode());
325 }
virtual bool isValid() override final
Can the handle be successfully dereferenced?