ATLAS Offline Software
Loading...
Searching...
No Matches
TrackPairsSelector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5/***************************************************************************
6 TrackPairsSelector.cxx - Description
7 ------------------
8 begin : 01-01-2008
9 authors : Tatjana Lenz, Thomas Koffas
10 email : tatjana.lenz@cern.ch, Thomas.Koffas@cern.ch
11 changes : M.Elsing
12***************************************************************************/
16
18
19
20// Would be nice to use something like Amg::distance instead.
21// But that rounds slightly differently.
22// Do it like this so that results are identical with the pre-MT version.
23namespace {
24 inline double square(const double tosquare) {
25 return std::pow(tosquare,2);
26 }
27 double dist(const std::pair<Amg::Vector3D,Amg::Vector3D>& pairofpos) {
28 Amg::Vector3D diff(pairofpos.first-pairofpos.second);
29 return std::sqrt(square(diff.x())+square(diff.y())+square(diff.z()));
30 }
31}
32
33
34namespace InDet {
35
36 // -------------------------------------------------------------
37 static const InterfaceID IID_ITrackPairsSelector("InDet::TrackPairsSelector", 1, 0);
38
39 // -------------------------------------------------------------
41 const std::string& name,
42 const IInterface* parent)
43 : AthAlgTool(type, name, parent)
44 {
45 declareInterface<TrackPairsSelector>(this);
46 }
47
48 // -------------------------------------------------------------
49 const InterfaceID& TrackPairsSelector::interfaceID() {
51 }
52
53 // -------------------------------------------------------------
55
56 /* Get the track selector tool from ToolSvc */
57 if ( m_helpertool.retrieve().isFailure() ) {
58 ATH_MSG_ERROR("Failed to retrieve tool " << m_helpertool);
59 return StatusCode::FAILURE;
60 }
61 ATH_MSG_DEBUG("Retrieved tool " << m_helpertool);
62
63 /* Get the distance tool from ToolsSvc */
64 if(m_distanceTool.retrieve().isFailure()) {
65 ATH_MSG_ERROR("Could not get " << m_distanceTool);
66 return StatusCode::FAILURE;
67 }
68 ATH_MSG_DEBUG("Got the distance tool " << m_distanceTool);
69
70 ATH_MSG_DEBUG("Initialization successful");
71 return StatusCode::SUCCESS;
72 }
73
74 // -------------------------------------------------------------
76 return StatusCode::SUCCESS;
77 }
78
79 // -------------------------------------------------------------
80 bool
82 const xAOD::TrackParticle* trkPpos,
83 const xAOD::TrackParticle* trkPneg,
84 TrackPairsSelector::Cache& cache) const
85 {
86
87 bool pass = true;
88 //Getting the track perigee parameters
89 const Trk::TrackParameters* perPos = &(trkPpos->perigeeParameters());
90 const Trk::TrackParameters* perNeg = &(trkPneg->perigeeParameters());
91 if (!(InDet::ConversionFinderUtils::momFraction(perPos, perNeg))) pass = false;
92
93 //Track summary information
94
95
96 uint8_t nclusPos(0);
97 uint8_t dummy(0);
98 if(trkPpos->summaryValue(dummy,xAOD::numberOfSCTHits)){
99 nclusPos += dummy;
100 }
101 if(trkPpos->summaryValue(dummy,xAOD::numberOfPixelHits)){
102 nclusPos += dummy;
103 }
104
105 uint8_t nclusNeg(0);
106 if(trkPneg->summaryValue(dummy,xAOD::numberOfSCTHits)){
107 nclusNeg += dummy;
108 }
109 if(trkPneg->summaryValue(dummy,xAOD::numberOfPixelHits)){
110 nclusNeg += dummy;
111 }
112
113
114 int sCase = 100;
115 if(nclusNeg>0 && nclusPos>0) sCase = 0;
116 if((nclusNeg>0 && nclusPos==0) || (nclusNeg==0 && nclusPos>0)) sCase = 1;
117 if(nclusNeg==0 && nclusPos==0) sCase = 2;
118
119 //Position of first hit in track particle
122
123 int index(-1);
124 for(unsigned int i(0); i< trkPpos->numberOfParameters() ; ++i ){
125 if( xAOD::FirstMeasurement == trkPpos->parameterPosition(i) ){
126 index = i;
127 break;
128 }
129 }
130 if(index!=-1){
131 parPos = trkPpos->curvilinearParameters(index);
132 } else {
133 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
134 return false;
135 }
136
137 index = -1;
138 for(unsigned int i(0); i< trkPneg->numberOfParameters() ; ++i ){
139 if( xAOD::FirstMeasurement == trkPneg->parameterPosition(i) ){
140 index = i;
141 break;
142 }
143 }
144 if(index!=-1){
145 parNeg = trkPneg->curvilinearParameters(index);
146 } else {
147 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
148 return false;
149 }
150
151 double firstRpos = parPos.position().perp();
152 double firstRneg = parNeg.position().perp();
153
154 //Cut on Deta
155 double detaCut = 0.0;
156 if(sCase == 0) {
157 detaCut = m_etaCut[0];
158 } else if(sCase == 1) {
159 detaCut = m_etaCut[1];
160 } else if(sCase == 2) {
161 if(fabs(perPos->eta())<0.6 && fabs(perNeg->eta())<0.6) detaCut = 10000.; // No eta cut for barrel TRT tracks
162 else detaCut = m_etaCut[2];
163 }
164
165 cache.m_deltaCotTheta = fabs(1. / tan(perPos->parameters()[Trk::theta]) -
166 1. / tan(perNeg->parameters()[Trk::theta]));
167 if (cache.m_deltaCotTheta > detaCut)
168 return false;
169
170 //Cut on distance between the initial hit position of the two tracks.
171 double dinit = 1000.;
172 if(sCase == 0) {
173 dinit = m_initCut[0];
174 } else if(sCase == 1) {
175 dinit = m_initCut[1];
176 } else if(sCase == 2) {
177 dinit = m_initCut[2];
178 }
179
180 cache.m_deltaInit = fabs(firstRpos - firstRneg);
181 if (cache.m_deltaInit > dinit) return false;
182
183 //Cut on distance of minimum approach between the two tracks.
184 double maxDist = 1000.;
185 if(sCase == 0) {
186 maxDist = m_maxDist[0];
187 } else if(sCase == 1) {
188 maxDist = m_maxDist[1];
189 } else if(sCase == 2) {
190 maxDist = m_maxDist[2];
191 }
192
193 cache.m_distance = 1000000.;
194 std::optional<Trk::ITrkDistanceFinder::TwoPoints> result
195 = m_distanceTool->CalculateMinimumDistance(trkPneg->perigeeParameters(),
196 trkPpos->perigeeParameters() );
197 if (!result) return false;
198 cache.m_distance = dist (result.value());
199 if (cache.m_distance>maxDist) return false;
200
201 //3D angle cut in the case of V0s, not used in the case of conversions
202 double d_beta = (perPos->momentum().dot(perNeg->momentum())) /
203 (perPos->momentum().mag() * perNeg->momentum().mag());
204 if(d_beta <m_MinTrkAngle) pass = false;
205
206 return pass;
207 }
208
209 // -------------------------------------------------------------
210 bool
212 const Trk::Track* trkneg) const
213 {
214
215 bool pass = true;
217 const Trk::Perigee* perPos = trkpos->perigeeParameters();
218 const Trk::Perigee* perNeg = trkneg->perigeeParameters();
219 if (!(InDet::ConversionFinderUtils::momFraction(perPos, perNeg))) pass = false;
220
222 double init_pos = 0.; double init_neg = 0.;
225 if(!mb_pos || !mb_neg) {pass = false; return pass;}
227 init_pos = (*itp_pos)->globalPosition().perp();
229 init_neg = (*itp_neg)->globalPosition().perp();
230 int sCase = 100;
231 if (init_neg<=m_maxR && init_pos<=m_maxR) sCase = 0;
232 if ((init_neg<=m_maxR && init_pos>m_maxR) || (init_neg>m_maxR && init_pos<=m_maxR)) sCase = 1;
233 if (init_neg>m_maxR && init_pos>m_maxR) sCase = 2;
234
235 //Cut on Deta
236 double detaCut = 0.0;
237 if(sCase == 0) {
238 detaCut = m_etaCut[0];
239 } else if(sCase == 1) {
240 detaCut = m_etaCut[1];
241 } else if(sCase == 2) {
242 detaCut = m_etaCut[2];
243 }
244
245 if (fabs(1. / tan(perPos->parameters()[Trk::theta]) -
246 1. / tan(perNeg->parameters()[Trk::theta])) > detaCut)
247 pass = false;
248
249 //Cut on distance between the initial hit position of the two tracks.
250 double dinit = 1000.;
251 if(sCase == 0) {
252 dinit = m_initCut[0];
253 } else if(sCase == 1) {
254 dinit = m_initCut[1];
255 } else if(sCase == 2) {
256 dinit = m_initCut[2];
257 }
258 if(fabs(init_pos - init_neg) > dinit) pass = false;
259
260 //Cut on distance of minimum approach between the two tracks.
261 double maxDist = 1000.;
262 if(sCase == 0) {
263 maxDist = m_maxDist[0];
264 } else if(sCase == 1) {
265 maxDist = m_maxDist[1];
266 } else if(sCase == 2) {
267 maxDist = m_maxDist[2];
268 }
269
270 double newDistance = 1000000.;
271 std::optional<Trk::ITrkDistanceFinder::TwoPoints> result
272 = m_distanceTool->CalculateMinimumDistance(*trkpos, *trkneg);
273 if (!result) {
274 pass = false;
275 }
276 else {
277 newDistance = dist (result.value());
278 if (newDistance>maxDist) pass = false;
279 }
280
281 //3D angle cut in the case of V0s, not used in the case of conversions
282 double d_beta = (perPos->momentum().dot(perNeg->momentum())) /
283 (perPos->momentum().mag() * perNeg->momentum().mag());
284 if(d_beta <m_MinTrkAngle) pass = false;
285
286 return pass;
287 }
288
289 // -------------------------------------------------------------
290 std::map<std::string, float>
292 const TrackPairsSelector::Cache& cache)
293 {
294 return {{"minimumDistanceTrk", cache.m_distance},
295 {"deltaCotThetaTrk", cache.m_deltaCotTheta},
296 {"deltaInitRadius", cache.m_deltaInit} };
297 }
298
299} // namespace InDet
300
301
302
303
Scalar mag() const
mag method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
Definition Jet.cxx:631
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Derived DataVector<T>.
Definition DataVector.h:795
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
static double momFraction(const Trk::TrackParameters *per1, const Trk::TrackParameters *per2)
helper functions
ToolHandle< Trk::ITrkDistanceFinder > m_distanceTool
Distance of minimum approach tool.
static const InterfaceID & interfaceID()
ToolHandle< InDet::ConversionFinderUtils > m_helpertool
Conversion helper tool.
DoubleArrayProperty m_etaCut
bool selectTrackParticlePair(const xAOD::TrackParticle *trkPpos, const xAOD::TrackParticle *trkPneg, Cache &cache) const
Track pair selectors.Return true if the argument track fulfills the selection.
virtual StatusCode finalize() override
bool selectTrackPair(const Trk::Track *trkpos, const Trk::Track *trkneg) const
DoubleProperty m_maxR
Properties for track selection: all cuts are ANDed.
DoubleArrayProperty m_initCut
static std::map< std::string, float > getLastValues(const Cache &cache)
Return a map with the values calculated for the last pair to decorate the vertex once it is created.
TrackPairsSelector(const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode initialize() override
DoubleArrayProperty m_maxDist
double eta() const
Access method for pseudorapidity - from momentum.
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
const DataVector< const MeasurementBase > * measurementsOnTrack() const
return a pointer to a vector of MeasurementBase (NOT including any that come from outliers).
const Perigee * perigeeParameters() const
return Perigee.
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
const Trk::CurvilinearParameters curvilinearParameters(unsigned int index) const
Returns a curvilinear representation of the parameters at 'index'.
xAOD::ParameterPosition parameterPosition(unsigned int index) const
Return the ParameterPosition of the parameters at 'index'.
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
size_t numberOfParameters() const
Returns the number of additional parameters stored in the TrackParticle.
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
static const InterfaceID IID_ITrackPairsSelector("InDet::TrackPairsSelector", 1, 0)
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
@ theta
Definition ParamDefs.h:66
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition dot.py:1
Definition index.py:1
TrackParticle_v1 TrackParticle
Reference the current persistent version:
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.
float m_deltaCotTheta
Delta cot theta between the tracks.
float m_distance
Distance of closest approach between the tracks.
float m_deltaInit
Distance difference between initial hits of tracks.