ATLAS Offline Software
Loading...
Searching...
No Matches
Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Header include
9
11
16
17#include "TH1D.h"
18#include "TNtuple.h"
19#include "TTree.h"
20#include "TROOT.h"
21
22#include <iostream>
23#include <map>
24#include <vector>
25#include <deque>
26
27//-------------------------------------------------
28
29using namespace std;
30
31namespace VKalVrtAthena {
32
33 //____________________________________________________________________________________________________
35
36 bool is_pv_associated = false;
37 const xAOD::TrackParticle *trk_from_gsf;
38 // get ID track matched to GSF track
40 for( const auto* vtx : *vertices ) {
41 for( size_t iv = 0; iv < vtx->nTrackParticles(); iv++ ) {
42 const auto* pvtrk = vtx->trackParticle( iv );
43 if (pvtrk == nullptr) continue;
44 // when using lepton-only selection, also need to check if the ID track matched to the GSF track is associated to the PV
45 if ( (trk_from_gsf == pvtrk) or (trk == pvtrk) ) {
46 is_pv_associated = true;
47 break;
48 }
49 }
50 }
51 return is_pv_associated;
52 }
53
54 //____________________________________________________________________________________________________
55 double vtxVtxDistance( const Amg::Vector3D& v1, const Amg::Vector3D& v2 ) {
56 return (v1-v2).norm();
57 }
58
59
60 //____________________________________________________________________________________________________
61 double VrtSecInclusive::significanceBetweenVertices( const WrkVrt& v1, const WrkVrt& v2 ) const {
62 try {
63 const auto distance = v2.vertex - v1.vertex;
64 AmgSymMatrix(3) sumCov;
65
66 sumCov.fillSymmetric(0, 0, v1.vertexCov.at(0) + v2.vertexCov.at(0));
67 sumCov.fillSymmetric(1, 0, v1.vertexCov.at(1) + v2.vertexCov.at(1));
68 sumCov.fillSymmetric(1, 1, v1.vertexCov.at(2) + v2.vertexCov.at(2));
69 sumCov.fillSymmetric(2, 0, v1.vertexCov.at(3) + v2.vertexCov.at(3));
70 sumCov.fillSymmetric(2, 1, v1.vertexCov.at(4) + v2.vertexCov.at(4));
71 sumCov.fillSymmetric(2, 2, v1.vertexCov.at(5) + v2.vertexCov.at(5));
72
73 const double s2 = distance.transpose() * sumCov.inverse() * distance;
74
75 return s2 > 0. ? sqrt( s2 ) : AlgConsts::maxValue;
76 } catch(...) {
77 ATH_MSG_WARNING( " >>> " << __FUNCTION__ << ": detected covariance matrix broken exception" );
79 }
80 }
81
82 //____________________________________________________________________________________________________
83 template<>
84 void genSequence( const xAOD::Muon*, std::vector<unsigned>& trackTypes ) {
85 trackTypes = { xAOD::Muon::Primary,
86 xAOD::Muon::InnerDetectorTrackParticle,
87 xAOD::Muon::MuonSpectrometerTrackParticle,
88 xAOD::Muon::CombinedTrackParticle,
89 xAOD::Muon::ExtrapolatedMuonSpectrometerTrackParticle
90 };
91 }
92
93 template<>
94 void genSequence( const xAOD::Electron* electron, std::vector<unsigned>& trackTypes ) {
95 for( size_t i=0; i<electron->nTrackParticles(); i++ ) trackTypes.emplace_back( i );
96 }
97
98
99 //____________________________________________________________________________________________________
100 template<>
101 const xAOD::TrackParticle* getLeptonTrackParticle( const xAOD::Muon* muon, const unsigned& trackType ) {
102 return muon->trackParticle( static_cast<xAOD::Muon::TrackParticleType>( trackType ) );
103 }
104
105 template<>
106 const xAOD::TrackParticle* getLeptonTrackParticle( const xAOD::Electron* electron, const unsigned& trackType ) {
107 return electron->trackParticle( trackType );
108 }
109
110
111 //____________________________________________________________________________________________________
112 double VrtSecInclusive::distanceBetweenVertices( const WrkVrt& v1, const WrkVrt& v2 ) const {
113 return (v2.vertex - v1.vertex).norm();
114 }
115
116
117 //____________________________________________________________________________________________________
118 StatusCode VrtSecInclusive::disassembleVertex(std::vector<WrkVrt> *workVerticesContainer, const unsigned& iv )
119 {
120
121 auto& wrkvrt = workVerticesContainer->at(iv);
122
123 ATH_MSG_VERBOSE(" >> disassembleVertex(): begin: disassembling vertex[" << iv << "], workVerticesContainer.size() = " << workVerticesContainer->size() );
124 ATH_MSG_VERBOSE(" >> disassembleVertex(): Vertex: r = " << wrkvrt.vertex.perp() << ", z = " << wrkvrt.vertex.z() );
125
126 // Loop over the tracks associated to the vertex and slect the maximum chi2 track
127 const auto& ntrk = wrkvrt.selectedTrackIndices.size();
128 size_t maxChi2TrackIndex = AlgConsts::invalidUnsigned;
129
130 // If the size of the tracks is less than 2, this algorithm is meaningless.
131 if( wrkvrt.selectedTrackIndices.size() <= 2 ) return StatusCode::SUCCESS;
132
133 for( auto& index : wrkvrt.selectedTrackIndices ) {
134 const xAOD::TrackParticle* trk = m_selectedTracks.at( index );
135
136 ATH_MSG_VERBOSE(" >> disassembleVertex(): > track at vertex[" << iv << "]: "
137 << "index = " << trk->index()
138 << ", pT = " << trk->pt()
139 << ", phi = " << trk->phi()
140 << ", d0 = " << trk->d0()
141 << ", z0 = " << trk->z0());
142 }
143
144 // find the track with the maximum chi2
145 const auto& max = std::max_element( wrkvrt.Chi2PerTrk.begin(), wrkvrt.Chi2PerTrk.end() );
146
147 if( max == wrkvrt.Chi2PerTrk.end() ) return StatusCode::SUCCESS;
148
149 maxChi2TrackIndex = max - wrkvrt.Chi2PerTrk.begin();
150
151 // return if no track is found.
152 if(maxChi2TrackIndex == AlgConsts::invalidUnsigned ) return StatusCode::SUCCESS;
153
154
155 // define work variables
156 vector<const xAOD::NeutralParticle*> dummyNeutrals;
157
158 vector<WrkVrt> new_vertices;
159
160 // Loop over the tracks associated to the vertex other than the selected tracks
161 ATH_MSG_VERBOSE(" >> disassembleVertex(): Loop over the tracks associated to the vertex other than the selected tracks.");
162 for(size_t itrk=0; itrk<ntrk; itrk++) {
163
164 ATH_MSG_VERBOSE(" >> disassembleVertex(): > Loop itrk = " << itrk << " / " << ntrk );
165
166 // reject the selected track
167 if( itrk == maxChi2TrackIndex ) {
168 ATH_MSG_VERBOSE(" >> disassembleVertex(): > skipped." );
169 continue;
170 }
171
172 const size_t this_trk_id = wrkvrt.selectedTrackIndices[itrk];
173 const size_t selected_trk_id = wrkvrt.selectedTrackIndices[maxChi2TrackIndex];
174
175 ATH_MSG_VERBOSE(" >> disassembleVertex(): > this_trk_id = " << this_trk_id << ", selected_trk_id = " << selected_trk_id << ", alltrks_size = " << m_selectedTracks.size() );
176 if( this_trk_id >= m_selectedTracks.size() ) {
177 ATH_MSG_VERBOSE(" >> disassembleVertex(): > this_trk_id is invalid. continue!" );
178 continue;
179 }
180 if( selected_trk_id >= m_selectedTracks.size() ) {
181 ATH_MSG_VERBOSE(" >> disassembleVertex(): > selected_trk_id is invalid. continue!" );
182 continue;
183 }
184
185 ATH_MSG_VERBOSE(" >> disassembleVertex(): > Storing tracks to ListBaseTracks" );
186 ATH_MSG_VERBOSE(" >> disassembleVertex(): > m_selectedTracks.at( this_trk_id ) = " << m_selectedTracks.at( this_trk_id )->index() );
187 ATH_MSG_VERBOSE(" >> disassembleVertex(): > m_selectedTracks.at( this_trk_id ) = " << m_selectedTracks.at( selected_trk_id )->index() );
188
189 vector<const xAOD::TrackParticle*> ListBaseTracks;
190 ListBaseTracks.emplace_back( m_selectedTracks.at( this_trk_id ) );
191 ListBaseTracks.emplace_back( m_selectedTracks.at( selected_trk_id ) );
192
193 ATH_MSG_VERBOSE(" >> disassembleVertex(): > ListBaseTracks was stored." );
194
195 WrkVrt newvrt;
196 newvrt.selectedTrackIndices.emplace_back( this_trk_id );
197 newvrt.selectedTrackIndices.emplace_back( selected_trk_id );
198
199 // Fit the new vertex
200 ATH_MSG_VERBOSE(" >> disassembleVertex(): > Fast Fit" );
201
202 std::unique_ptr<Trk::IVKalState> state = m_fitSvc->makeState();
203 ATH_CHECK( m_fitSvc->VKalVrtFitFast( ListBaseTracks, newvrt.vertex, *state ) );
204
205 ATH_MSG_VERBOSE( " >> disassembleVertex(): > ApproxVertex: r = " << newvrt.vertex.perp() << ", z = " << newvrt.vertex.z() );
206
207 if( vtxVtxDistance( wrkvrt.vertex, newvrt.vertex ) > 10. )
208 {
209 m_fitSvc->setApproximateVertex( wrkvrt.vertex[0], wrkvrt.vertex[1], wrkvrt.vertex[2], *state );
210 }
211 else
212 {
213 m_fitSvc->setApproximateVertex( newvrt.vertex[0], newvrt.vertex[1], newvrt.vertex[2], *state );
214 }
215
216 ATH_MSG_VERBOSE(" >> disassembleVertex(): > Fit the new vertex" );
217 StatusCode sc = m_fitSvc->VKalVrtFit(ListBaseTracks,
218 dummyNeutrals,
219 newvrt.vertex,
220 newvrt.vertexMom,
221 newvrt.Charge,
222 newvrt.vertexCov,
223 newvrt.Chi2PerTrk,
224 newvrt.TrkAtVrt,
225 newvrt.Chi2,
226 *state);
227
228 if( sc.isFailure() ) continue;
229
230 newvrt.closestWrkVrtIndex = 0;
232
233 // register the new vertex to the vertex list
234 ATH_MSG_VERBOSE(" >> disassembleVertex(): > register the new vertex to the vertex list" );
235 new_vertices.emplace_back( newvrt );
236 }
237
238 // remove the selected track from the original vertex
239 wrkvrt.selectedTrackIndices.erase( wrkvrt.selectedTrackIndices.begin() + maxChi2TrackIndex ); //remove track
240 ATH_MSG_VERBOSE(" >> disassembleVertex(): removed the selected track from the original vertex. wrkvrt.selectedTrackIndices.size = " << wrkvrt.selectedTrackIndices.size() );
241
242 // refit the original vertex
243 ATH_MSG_VERBOSE(" >> disassembleVertex(): refit the original vertex" );
244
245 StatusCode sc = refitVertex( wrkvrt );
246 if( sc.isFailure() ) {
247 // WARNING CODE ATLASRECTS-3145::001 refitVertex Failure, vertex lost
248 ATH_MSG_WARNING("ATLASRECTS-3145::001" );
249 return StatusCode::SUCCESS;
250 }
251 // end of workaround
252
253 for( const auto& vertex : new_vertices ) {
254 ATH_MSG_VERBOSE(" >> disassembleVertex(): > emplace_back new vertex" );
255 workVerticesContainer->emplace_back( vertex );
256 }
257
258 ATH_MSG_VERBOSE(" >> disassembleVertex(): end. workVerticesContainer.size() = " << workVerticesContainer->size() );
259 return StatusCode::SUCCESS;
260 }
261
262
263
264 //____________________________________________________________________________________________________
266 {
267 //
268 // Iterate track removal until vertex get good Chi2
269 //
270
271 auto fitQuality_begin = vertex.fitQuality();
272
273 auto removeCounter = 0;
274
275 if( vertex.nTracksTotal() <= 2 ) return 0.;
276
277 {
278 WrkVrt backup = vertex;
279 StatusCode sc = refitVertexWithSuggestion( vertex, vertex.vertex );
280 if( sc.isFailure() ) {
281 vertex = backup;
282 return 0;
283 }
284 }
285
286 double chi2Probability = TMath::Prob( vertex.Chi2, vertex.ndof() );
287
288 while (chi2Probability < m_jp.improveChi2ProbThreshold ) {
289 if( vertex.nTracksTotal() == 2 ) return chi2Probability;
290
291 WrkVrt vertex_backup = vertex;
292
293 auto maxChi2 = std::max_element( vertex.Chi2PerTrk.begin(), vertex.Chi2PerTrk.end() );
294 size_t index = maxChi2 - vertex.Chi2PerTrk.begin();
295
296
297 ATH_MSG_DEBUG( " >>> " << __FUNCTION__ << ": worst chi2 trk index = " << index << ", #trks = " << vertex.Chi2PerTrk.size() );
298
299 if( index < vertex.selectedTrackIndices.size() ) {
300 vertex.selectedTrackIndices.erase( vertex.selectedTrackIndices.begin() + index ); //remove track
301 removeCounter++;
302 } else {
303 index -= vertex.selectedTrackIndices.size();
304 if( index >= vertex.associatedTrackIndices.size() ) {
305 ATH_MSG_WARNING( " >>> " << __FUNCTION__ << ": invalid index" );
306 break;
307 }
308 vertex.associatedTrackIndices.erase( vertex.associatedTrackIndices.begin() + index ); //remove track
309 removeCounter++;
310 }
311
312 StatusCode sc = refitVertexWithSuggestion( vertex, vertex.vertex );
313
314 if( sc.isFailure() || vertex_backup.fitQuality() < vertex.fitQuality() ) {
315 vertex = vertex_backup;
316 chi2Probability = 0;
317 break;
318 }
319
320 chi2Probability = TMath::Prob( vertex.Chi2, vertex.ndof() );
321 }
322
323 auto fitQuality_end = vertex.fitQuality();
324
325 if( 0 == removeCounter ) {
326 ATH_MSG_DEBUG( " >>> " << __FUNCTION__ << ": no improvement was found." );
327 } else {
328 ATH_MSG_DEBUG( " >>> " << __FUNCTION__ << ": Removed " << removeCounter << " tracks; Fit quality improvement: " << fitQuality_begin << " ==> " << fitQuality_end );
329 }
330
331 return chi2Probability;
332 }
333
334
335
336 //____________________________________________________________________________________________________
337 void VrtSecInclusive::removeTrackFromVertex(std::vector<WrkVrt> *workVerticesContainer,
338 std::vector< std::deque<long int> > *TrkInVrt,
339 const long int & trackIndexToRemove,
340 const long int & SelectedVertex)
341 {
342
343 auto& wrkvrt = workVerticesContainer->at(SelectedVertex);
344 auto& tracks = wrkvrt.selectedTrackIndices;
345
346 {
347 auto end = std::remove_if( tracks.begin(), tracks.end(), [&](long int index) { return index == trackIndexToRemove; } );
348 tracks.erase( end, tracks.end() );
349 }
350
351 {
352 for( auto& trks : *TrkInVrt ) {
353 auto end = std::remove_if( trks.begin(), trks.end(), [&](long int index) { return index == trackIndexToRemove; } );
354 trks.erase( end, trks.end() );
355 }
356 }
357
358
359 // Check if track is removed from two-track vertex => then sharing of track left should also be decreased
360 if( wrkvrt.selectedTrackIndices.size() == 1 ) {
361
362 const auto& leftTrackIndex = *( tracks.begin() );
363 auto& list = TrkInVrt->at(leftTrackIndex);
364 auto end = std::remove_if( list.begin(), list.end(), [&](long int index) { return index == trackIndexToRemove; } );
365 list.erase( end, list.end() );
366
367 }
368
369 }
370
371
372
373 //____________________________________________________________________________________________________
374 double VrtSecInclusive::findMinVerticesPair( std::vector<WrkVrt> *workVerticesContainer, std::pair<unsigned, unsigned>& indexPair, const VrtSecInclusive::AlgForVerticesPair& algorithm )
375 {
376 //
377 // Minimal normalized vertex-vertex distance
378 //
379
380 for( auto& workVertex : *workVerticesContainer ) {
381 workVertex.closestWrkVrtValue = AlgConsts::maxValue;
382 workVertex.closestWrkVrtIndex = 0;
383 }
384
386
387 for( auto iv = workVerticesContainer->begin(); iv != workVerticesContainer->end(); ++iv ) {
388 if( (*iv).selectedTrackIndices.size()< 2) continue; /* Bad vertices */
389
390 auto i_index = iv - workVerticesContainer->begin();
391
392 for( auto jv = std::next(iv); jv != workVerticesContainer->end(); ++jv ) {
393 if( (*jv).selectedTrackIndices.size()< 2) continue; /* Bad vertices */
394
395 auto j_index = iv - workVerticesContainer->begin();
396
397 double value = (this->*algorithm)( (*iv), (*jv) );
398
399 if( value < minValue ){
400 minValue = value;
401 indexPair.first = i_index;
402 indexPair.second = j_index;
403 }
404 if( value < (*iv).closestWrkVrtValue ) {(*iv).closestWrkVrtValue = value; (*iv).closestWrkVrtIndex = j_index; }
405 if( value < (*jv).closestWrkVrtValue ) {(*jv).closestWrkVrtValue = value; (*jv).closestWrkVrtIndex = i_index; }
406 }
407 }
408
409 return minValue;
410 }
411
412
413 //____________________________________________________________________________________________________
414 double VrtSecInclusive::findMinVerticesNextPair( std::vector<WrkVrt> *workVerticesContainer, std::pair<unsigned, unsigned>& indexPair )
415 {
416 //
417 // Give minimal distance between nonmodifed yet vertices
418 //
419
420 indexPair.first = 0;
421 indexPair.second = 0;
422
424
425 for(unsigned iv=0; iv<workVerticesContainer->size()-1; iv++) {
426 auto& vertex = workVerticesContainer->at(iv);
427
428 if( vertex.selectedTrackIndices.size() < 2) continue; /* Bad vertex */
429 if( vertex.closestWrkVrtIndex == 0 ) continue; /* Used vertex */
430
431 if( vertex.closestWrkVrtValue < minValue ) {
432
433 unsigned jv = vertex.closestWrkVrtIndex;
434
435 // Close vertex to given [iv] one is modified already
436 if( workVerticesContainer->at(jv).closestWrkVrtIndex == 0 ) continue;
437
438 minValue = vertex.closestWrkVrtValue;
439
440 indexPair.first = iv;
441 indexPair.second = jv;
442
443 }
444 }
445
446 return minValue;
447 }
448
449
450
451 //____________________________________________________________________________________________________
453 {
454 //
455 // Merge two close vertices into one (first) and set NTr=0 for second vertex
456 //
457
458 // firstly, take a backup of the original vertices
459 auto v1_bak = v1;
460 auto v2_bak = v2;
461
462 for( auto& index : v2.selectedTrackIndices ) { v1.selectedTrackIndices.emplace_back( index ); }
463
464 // Cleaning
465 deque<long int>::iterator TransfEnd;
466 sort( v1.selectedTrackIndices.begin(), v1.selectedTrackIndices.end() );
467 TransfEnd = unique(v1.selectedTrackIndices.begin(), v1.selectedTrackIndices.end() );
468 v1.selectedTrackIndices.erase( TransfEnd, v1.selectedTrackIndices.end());
469 //
470 //----------------------------------------------------------
471 v2.selectedTrackIndices.clear(); //Clean dropped vertex
472 v2.closestWrkVrtValue = AlgConsts::maxValue; //Clean dropped vertex
473 v2.closestWrkVrtIndex = 0; //Clean dropped vertex
474 v2.isGood = false; //Clean dropped vertex
475
476 v1.closestWrkVrtValue = AlgConsts::maxValue; //Clean new vertex
477 v1.closestWrkVrtIndex = 0; //Clean new vertex
478 v1.isGood = true; //Clean new vertex
479
480 StatusCode sc = refitVertex( v1 );
481 if( sc.isFailure() ) {
482 v1 = v1_bak;
483 v2 = v2_bak;
484
485 ATH_MSG_DEBUG(" >>> " << __FUNCTION__ << ": failure in merging" );
486
487 return StatusCode::FAILURE;
488 }
489
490 return StatusCode::SUCCESS;
491 }
492
493
494
495 //____________________________________________________________________________________________________
496 StatusCode VrtSecInclusive::refitVertex( WrkVrt& workVertex )
497 {
498 std::unique_ptr<Trk::IVKalState> state = m_fitSvc->makeState();
499 return refitVertex (workVertex, *state);
500 }
501
502
503
504 //____________________________________________________________________________________________________
505 StatusCode VrtSecInclusive::refitVertex( WrkVrt& workVertex,
506 Trk::IVKalState& istate)
507 {
508
509 //
510 vector<const xAOD::NeutralParticle*> dummyNeutrals;
511
512 int nth = workVertex.selectedTrackIndices.size();
513
514 if(nth<2) {
515 workVertex.isGood = false;
516 return StatusCode::SUCCESS;
517 }
518
519 vector<const xAOD::TrackParticle*> ListBaseTracks;
520
521 workVertex.Chi2PerTrk.clear();
522
523 for( const auto& index : workVertex.selectedTrackIndices ) {
524 ListBaseTracks.emplace_back( m_selectedTracks.at( index ) );
525 workVertex.Chi2PerTrk.emplace_back( AlgConsts::chi2PerTrackInitValue );
526 }
527
528 for( const auto& index : workVertex.associatedTrackIndices ) {
529 ListBaseTracks.emplace_back( m_associatedTracks.at( index ) );
530 workVertex.Chi2PerTrk.emplace_back( AlgConsts::chi2PerTrackInitValue );
531 }
532
533 auto& vertexPos = workVertex.vertex;
534
535 m_fitSvc->setApproximateVertex( vertexPos.x(), vertexPos.y(), vertexPos.z(), istate );
536
537 ATH_MSG_VERBOSE( " >>> refitVertex: ListBaseTracks.size = " << ListBaseTracks.size()
538 << ", #selectedBaseTracks = " << workVertex.selectedTrackIndices.size()
539 << ", #assocTracks = " << workVertex.associatedTrackIndices.size() );
540 for( const auto *trk : ListBaseTracks ) {
541 ATH_MSG_VERBOSE( " >>> refitVertex: track index = " << trk->index() );
542 }
543
544 ATH_MSG_VERBOSE( " >>> refitVertex: m_fitSvc is reset." );
545
546 Amg::Vector3D initVertex;
547
548 StatusCode sc = m_fitSvc->VKalVrtFitFast( ListBaseTracks, initVertex, istate );/* Fast crude estimation */
549 if(sc.isFailure()) ATH_MSG_DEBUG(" >>> refitVertex: fast crude estimation failed.");
550 ATH_MSG_VERBOSE( " >>> refitVertex: Fast VKalVrtFit succeeded. vertex (r,z) = (" << initVertex.perp() << ", " << initVertex.z() << ", " << ")" );
551
552 if( vtxVtxDistance( initVertex, vertexPos ) > 10. ) {
553
554 m_fitSvc->setApproximateVertex( vertexPos.x(), vertexPos.y(), vertexPos.z(), istate );
555
556 } else {
557
558 m_fitSvc->setApproximateVertex( initVertex.x(), initVertex.y(), initVertex.z(), istate );
559
560 }
561
562 ATH_MSG_VERBOSE( " >>> refitVertex: approx vertex is set. Now going to perform fitting..." );
563
564 StatusCode SC=m_fitSvc->VKalVrtFit(ListBaseTracks,dummyNeutrals,
565 workVertex.vertex,
566 workVertex.vertexMom,
567 workVertex.Charge,
568 workVertex.vertexCov,
569 workVertex.Chi2PerTrk,
570 workVertex.TrkAtVrt,
571 workVertex.Chi2,
572 istate);
573
574 auto& cov = workVertex.vertexCov;
575
576 if(SC.isFailure()) ATH_MSG_DEBUG(" >>> refitVertex: SC in refitVertex returned failure ");
577 ATH_MSG_VERBOSE(" >>> refitVertex "<<SC<<", "<<ListBaseTracks.size()<<","<<workVertex.Chi2PerTrk.size());
578 ATH_MSG_VERBOSE( " >>> refitVertex: succeeded in fitting. New vertex pos (r,z) = (" << vertexPos.perp() << ", " << vertexPos.z() << ")" );
579 ATH_MSG_VERBOSE( " >>> refitVertex: New vertex cov = (" << cov.at(0) << ", " << cov.at(1) << ", " << cov.at(2) << ", " << cov.at(3) << ", " << cov.at(4) << ", " << cov.at(5) << ")" );
580
581 return SC;
582 }
583
584 //____________________________________________________________________________________________________
585 StatusCode VrtSecInclusive::refitVertexWithSuggestion( WrkVrt& workVertex, const Amg::Vector3D& suggestedPosition )
586 {
587 std::unique_ptr<Trk::IVKalState> state = m_fitSvc->makeState();
588 return refitVertexWithSuggestion (workVertex, suggestedPosition, *state);
589 }
590
591
592 //____________________________________________________________________________________________________
594 const Amg::Vector3D& suggestedPosition,
595 Trk::IVKalState& istate )
596 {
597
598 //
599 vector<const xAOD::NeutralParticle*> dummyNeutrals;
600
601 int nth = workVertex.selectedTrackIndices.size();
602
603 if(nth<2) {
604 workVertex.isGood = false;
605 return StatusCode::SUCCESS;
606 }
607
608 vector<const xAOD::TrackParticle*> ListBaseTracks;
609
610 workVertex.Chi2PerTrk.clear();
611
612 for( const auto& index : workVertex.selectedTrackIndices ) {
613 ListBaseTracks.emplace_back( m_selectedTracks.at( index ) );
614 workVertex.Chi2PerTrk.emplace_back( AlgConsts::chi2PerTrackInitValue );
615 }
616
617 for( const auto& index : workVertex.associatedTrackIndices ) {
618 ListBaseTracks.emplace_back( m_associatedTracks.at( index ) );
619 workVertex.Chi2PerTrk.emplace_back( AlgConsts::chi2PerTrackInitValue );
620 }
621
622 auto& vertexPos = workVertex.vertex;
623
624 m_fitSvc->setApproximateVertex( suggestedPosition.x(), suggestedPosition.y(), suggestedPosition.z(), istate );
625
626 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ <<": ListBaseTracks.size = " << ListBaseTracks.size()
627 << ", #selectedBaseTracks = " << workVertex.selectedTrackIndices.size()
628 << ", #assocTracks = " << workVertex.associatedTrackIndices.size() );
629 for( const auto *trk : ListBaseTracks ) {
630 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ << ": track index = " << trk->index() );
631 }
632
633 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ << ": m_fitSvc is reset." );
634
635 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ << ": approx vertex is set. Now going to perform fitting..." );
636
637 StatusCode SC=m_fitSvc->VKalVrtFit(ListBaseTracks,dummyNeutrals,
638 workVertex.vertex,
639 workVertex.vertexMom,
640 workVertex.Charge,
641 workVertex.vertexCov,
642 workVertex.Chi2PerTrk,
643 workVertex.TrkAtVrt,
644 workVertex.Chi2,
645 istate);
646
647 auto& cov = workVertex.vertexCov;
648
649 if( SC.isFailure() ) ATH_MSG_VERBOSE(" >>> " << __FUNCTION__ << ": SC in refitVertex returned failure ");
650 ATH_MSG_VERBOSE(" >>> " << __FUNCTION__ << ": "<<SC<<", "<<ListBaseTracks.size()<<","<<workVertex.Chi2PerTrk.size());
651
652 if( SC.isSuccess() ) {
653 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ << ": succeeded in fitting. New vertex pos = (" << vertexPos.x() << ", " << vertexPos.y() << ", " << vertexPos.z() << ")" );
654 ATH_MSG_VERBOSE( " >>> " << __FUNCTION__ << ": New vertex cov = (" << cov.at(0) << ", " << cov.at(1) << ", " << cov.at(2) << ", " << cov.at(3) << ", " << cov.at(4) << ", " << cov.at(5) << ")" );
655 }
656
657 return SC;
658 }
659
660
661 //____________________________________________________________________________________________________
662 size_t VrtSecInclusive::nTrkCommon( std::vector<WrkVrt> *workVerticesContainer, const std::pair<unsigned, unsigned>& pairIndex)
663 {
664 //
665 // Number of common tracks for 2 vertices
666 //
667
668 auto& trackIndices1 = workVerticesContainer->at( pairIndex.first ).selectedTrackIndices;
669 auto& trackIndices2 = workVerticesContainer->at( pairIndex.second ).selectedTrackIndices;
670
671 if( trackIndices1.size() < 2 ) return 0;
672 if( trackIndices2.size() < 2 ) return 0;
673
674 size_t nTrkCom = 0;
675
676 for( auto& index : trackIndices1 ) {
677 if( std::find(trackIndices2.begin(),trackIndices2.end(), index) != trackIndices2.end()) nTrkCom++;
678 }
679
680 return nTrkCom;
681 }
682
683
684
685 //____________________________________________________________________________________________________
687
688 declareProperty("GeoModel", m_jp.geoModel = VKalVrtAthena::GeoModel::Run2 );
689
690 declareProperty("TrackLocation", m_jp.TrackLocation = "InDetTrackParticles" );
691 declareProperty("MuonLocation", m_jp.MuonLocation = "Muons" );
692 declareProperty("ElectronLocation", m_jp.ElectronLocation = "Electrons" );
693 declareProperty("PrimVrtLocation", m_jp.PrimVrtLocation = "PrimaryVertices" );
694 declareProperty("McParticleContainer", m_jp.truthParticleContainerName = "TruthParticles" );
695 declareProperty("MCEventContainer", m_jp.mcEventContainerName = "TruthEvents" );
696 declareProperty("AugmentingVersionString", m_jp.augVerString = "_VSI" );
697 declareProperty("TruthParticleFilter", m_jp.truthParticleFilter = "Rhadron" ); // Either "", "Kshort", "Rhadron", "HNL", "HadInt", "Bhadron"
698
699 declareProperty("All2trkVerticesContainerName", m_jp.all2trksVerticesContainerName = "All2TrksVertices" );
700 declareProperty("SecondaryVerticesContainerName", m_jp.secondaryVerticesContainerName = "SecondaryVertices" );
701
702 declareProperty("FillHist", m_jp.FillHist = false );
703 declareProperty("FillNtuple", m_jp.FillNtuple = false );
704 declareProperty("FillIntermediateVertices", m_jp.FillIntermediateVertices = false );
705 declareProperty("DoIntersectionPos", m_jp.doIntersectionPos = false );
706 declareProperty("DoMapToLocal", m_jp.doMapToLocal = false );
707 declareProperty("DoTruth", m_jp.doTruth = false );
708 declareProperty("DoPVcompatibility", m_jp.doPVcompatibilityCut = true );
709 declareProperty("DoTightPVcompatibility", m_jp.doTightPVcompatibilityCut = false );
710 declareProperty("RemoveFake2TrkVrt", m_jp.removeFakeVrt = true );
711 declareProperty("DoDelayedFakeReject", m_jp.removeFakeVrtLate = false );
712 declareProperty("CheckHitPatternStrategy", m_checkPatternStrategy = "Classical" ); // Either Classical or Extrapolation
713 declareProperty("MCTrackResolution", m_jp.mcTrkResolution = 0.06 ); // see getTruth for explanation
714 declareProperty("TruthTrkLen", m_jp.TruthTrkLen = 1000 ); // in [mm]
715 declareProperty("ExtrapPV", m_jp.extrapPV = false ); // Leave false. only for testing
716 declareProperty("PassThroughTrackSelection", m_jp.passThroughTrackSelection = false );
717 declareProperty("DoFastMode", m_jp.doFastMode = false );
718
719
720 declareProperty("DoTwoTrSoftBtag", m_jp.doTwoTrSoftBtag = false );
721 declareProperty("TwoTrVrtAngleCut", m_jp.twoTrVrtAngleCut = -10 );
722 declareProperty("TwoTrVrtMinDistFromPVCut", m_jp.twoTrVrtMinDistFromPV = 0. );
723
724 declareProperty("TruncateListOfWorkingVertices", m_jp.truncateWrkVertices = true );
725 declareProperty("MaxNumberOfWorkingVertices", m_jp.maxWrkVertices = 1500 );
726
727 // default values are set upstream - check top of file
728 declareProperty("do_PVvetoCut", m_jp.do_PVvetoCut = true );
729 declareProperty("do_d0Cut", m_jp.do_d0Cut = true );
730 declareProperty("do_z0Cut", m_jp.do_z0Cut = true );
731 declareProperty("do_d0errCut", m_jp.do_d0errCut = false );
732 declareProperty("do_z0errCut", m_jp.do_z0errCut = false );
733 declareProperty("do_d0signifCut", m_jp.do_d0signifCut = false );
734 declareProperty("do_z0signifCut", m_jp.do_z0signifCut = false );
735
736 declareProperty("ImpactWrtBL", m_jp.ImpactWrtBL = true ); // false option is going to be deprecated
737 declareProperty("a0TrkPVDstMinCut", m_jp.d0TrkPVDstMinCut = 0. ); // in [mm]
738 declareProperty("a0TrkPVDstMaxCut", m_jp.d0TrkPVDstMaxCut = 1000. ); // in [mm]
739 declareProperty("a0TrkPVSignifCut", m_jp.d0TrkPVSignifCut = 0. ); // in [mm]
740 declareProperty("twoTrkVtxFormingD0Cut", m_jp.twoTrkVtxFormingD0Cut = 1. ); // in [mm]
741 declareProperty("zTrkPVDstMinCut", m_jp.z0TrkPVDstMinCut = 0. ); // in [mm]
742 declareProperty("zTrkPVDstMaxCut", m_jp.z0TrkPVDstMaxCut = 1000. ); // in [mm]
743 declareProperty("zTrkPVSignifCut", m_jp.z0TrkPVSignifCut = 0. ); // in unit of sigma
744 declareProperty("TrkA0ErrCut", m_jp.d0TrkErrorCut = 10000 ); // in [mm]
745 declareProperty("TrkZErrCut", m_jp.z0TrkErrorCut = 20000 ); // in [mm]
746
747 declareProperty("SelTrkMaxCutoff", m_jp.SelTrkMaxCutoff = 50 ); // max number of tracks
748 declareProperty("TrkPtCut", m_jp.TrkPtCut = 1000. ); // low pT threshold. in [MeV]
749 declareProperty("TrkChi2Cut", m_jp.TrkChi2Cut = 3. ); // in terms of chi2 / ndof
750 declareProperty("PVcompatibilityCut", m_jp.pvCompatibilityCut = -20. ); // in [mm]
751 declareProperty("SelVrtChi2Cut", m_jp.SelVrtChi2Cut = 4.5 ); // in terms of chi2 / ndof
752
753 declareProperty("CutSctHits", m_jp.CutSctHits = 0 );
754 declareProperty("CutPixelHits", m_jp.CutPixelHits = 0 );
755 declareProperty("CutSiHits", m_jp.CutSiHits = 0 );
756 declareProperty("DoSAloneTRT", m_jp.SAloneTRT = false ); // SAlone = "standalone"
757 declareProperty("CutBLayHits", m_jp.CutBLayHits = 0 );
758 declareProperty("CutSharedHits", m_jp.CutSharedHits = 0 );
759 declareProperty("doTRTPixCut", m_jp.doTRTPixCut = false ); // mode for R-hadron displaced vertex
760 declareProperty("CutTRTHits", m_jp.CutTRTHits = 0 );
761 declareProperty("CutTightSCTHits", m_jp.CutTightSCTHits = 7 );
762 declareProperty("CutTightTRTHits", m_jp.CutTightTRTHits = 20 );
763
764 declareProperty("TrkExtrapolator", m_jp.trkExtrapolator = 2 );
765
766 declareProperty("doReassembleVertices", m_jp.doReassembleVertices = false );
767 declareProperty("doMergeByShuffling", m_jp.doMergeByShuffling = false );
768 declareProperty("doSuggestedRefitOnMerging", m_jp.doSuggestedRefitOnMerging = true ); // sub-option of doMergeByShuffling-1
769 declareProperty("doMagnetMerging", m_jp.doMagnetMerging = true ); // sub-option of doMergeByShuffling-2
770 declareProperty("doWildMerging", m_jp.doWildMerging = true ); // sub-option of doMergeByShuffling-3
771 declareProperty("doMergeFinalVerticesDistance", m_jp.doMergeFinalVerticesDistance = false );
772 declareProperty("doAssociateNonSelectedTracks", m_jp.doAssociateNonSelectedTracks = false );
773 declareProperty("doFinalImproveChi2", m_jp.doFinalImproveChi2 = false );
774
775 declareProperty("VertexMergeCut", m_jp.VertexMergeCut = 3 );
776 declareProperty("TrackDetachCut", m_jp.TrackDetachCut = 6 );
777 declareProperty("associateMinDistanceToPV", m_jp.associateMinDistanceToPV = 0.5 );
778 declareProperty("associateMaxD0Signif", m_jp.associateMaxD0Signif = 5. ); // wrt. DV in unit of sigma
779 declareProperty("associateMaxZ0Signif", m_jp.associateMaxZ0Signif = 5. ); // wrt. DV in unit of sigma
780 declareProperty("associatePtCut", m_jp.associatePtCut = 0. ); // in [MeV]
781 declareProperty("associateChi2Cut", m_jp.associateChi2Cut = 20. );
782 declareProperty("reassembleMaxImpactParameterD0", m_jp.reassembleMaxImpactParameterD0 = 1. ); // wrt. DV in [mm]
783 declareProperty("reassembleMaxImpactParameterZ0", m_jp.reassembleMaxImpactParameterZ0 = 5. ); // wrt. DV in [mm]
784 declareProperty("mergeByShufflingMaxSignificance", m_jp.mergeByShufflingMaxSignificance = 100. ); // in unit of sigma
785 declareProperty("mergeByShufflingAllowance", m_jp.mergeByShufflingAllowance = 4. ); // in unit of sigma
786 declareProperty("VertexMergeFinalDistCut", m_jp.VertexMergeFinalDistCut = 1. ); // in [mm]
787 declareProperty("VertexMergeFinalDistScaling", m_jp.VertexMergeFinalDistScaling = 0. ); // in [1/mm]
788 declareProperty("improveChi2ProbThreshold", m_jp.improveChi2ProbThreshold = 1.e-4 );
789
790 // A test implementation for muon vertices
791 declareProperty("doSelectTracksFromMuons", m_jp.doSelectTracksFromMuons = false );
792 declareProperty("doRemoveCaloTaggedMuons", m_jp.doRemoveCaloTaggedMuons = false );
793 declareProperty("doSelectTracksFromElectrons", m_jp.doSelectTracksFromElectrons = false );
794 declareProperty("doSelectIDAndGSFTracks", m_jp.doSelectIDAndGSFTracks = false );
795 declareProperty("doRemoveNonLeptonVertices", m_jp.doRemoveNonLeptonVertices = false );
796
797 // Disappearing track vertices
798 declareProperty("doDisappearingTrackVertexing", m_jp.doDisappearingTrackVertexing = false );
799 declareProperty("twoTrVrtMaxPerigeeDist", m_jp.twoTrVrtMaxPerigeeDist = 50 ); // in [mm]
800 declareProperty("twoTrVrtMinRadius", m_jp.twoTrVrtMinRadius = 50 ); // in [mm]
801
802
803
804 // Select tracks with additonal LRT Cuts (inspiried by Run 3 LRT optimization studies)
805 declareProperty("doSelectTracksWithLRTCuts", m_jp.doSelectTracksWithLRTCuts = false );
806
807 // Additional dressing option
808 declareProperty("doAugmentDVimpactParametersToMuons", m_jp.doAugmentDVimpactParametersToMuons = false );
809 declareProperty("doAugmentDVimpactParametersToElectrons", m_jp.doAugmentDVimpactParametersToElectrons = false );
810
811 // Additional ToolHandles
812 declareProperty("VertexFitterTool", m_fitSvc, " Private TrkVKalVrtFitter" );
813 declareProperty("Extrapolator", m_extrapolator );
814 declareProperty("TrackToVertexTool", m_trackToVertexTool );
815 declareProperty("TrackToVertexIPEstimatorTool", m_trackToVertexIPEstimatorTool );
816 declareProperty("VertexMapper", m_vertexMapper );
817 declareProperty("TruthToTrack", m_truthToTrack );
818
819 }
820
821
822
823
824 //____________________________________________________________________________________________________
826
827 //--------------------------------------------------------
828 // Primary vertex extraction from TES
829 //
830
831 ATH_CHECK( evtStore()->retrieve( m_primaryVertices, "PrimaryVertices") );
832
833 if( m_jp.FillNtuple ) m_ntupleVars->get<unsigned int>( "NumPV" ) = 0;
834 m_thePV = nullptr;
835
836 ATH_MSG_DEBUG( "processPrimaryVertices(): pv_size = " << m_primaryVertices->size() );
837
838 // Loop over PV container and get number of tracks of each PV
839
840 for( const auto *vertex : *m_primaryVertices ) {
841
842 // Hide (2015-04-21): xAOD::Vertex may contain several types of vertices
843 // e.g. if VertexType==NoVtx, this is a dummy vertex.
844 // We only need to extract primary vertices here, and skip otherwise.
845
846 if( xAOD::VxType::PriVtx != vertex->vertexType() ) continue;
847
848 // Not considering pile-up; pick-up the first PV
849 m_thePV = vertex;
850
851 if( m_jp.FillNtuple ) {
852
853 if( 0 == m_ntupleVars->get<unsigned int>( "NumPV" ) ) {
854
855 m_ntupleVars->get<double>( "PVX" ) = vertex->x();
856 m_ntupleVars->get<double>( "PVY" ) = vertex->y();
857 m_ntupleVars->get<double>( "PVZ" ) = vertex->z();
858 m_ntupleVars->get<unsigned int>( "PVType" ) = vertex->vertexType();
859
860 // number of tracks associated to the PV
861 m_ntupleVars->get<unsigned int>( "NTrksPV" ) = vertex->nTrackParticles();
862 }
863
864 m_ntupleVars->get<unsigned int>( "NumPV" )++;
865
866 m_ntupleVars->get< vector<int> > ( "NdofTrksPV" ) .emplace_back( vertex->numberDoF() );
867 m_ntupleVars->get< vector<double> >( "PVZpile" ) .emplace_back( vertex->position().z() );
868 }
869
870 ATH_MSG_DEBUG("PrimVertex x/y/z/nDOF "
871 << vertex->x() << ","
872 << vertex->y() << ","
873 << vertex->z() << ","
874 << vertex->numberDoF() );
875
876 }
877
878 // Use the dummy PV if no PV is composed
879 if( !m_thePV ) {
880 ATH_MSG_DEBUG("No Reconstructed PV was found. Using the dummy PV instead.");
881 for( const auto *vertex : *m_primaryVertices ) {
882 if( xAOD::VxType::NoVtx != vertex->vertexType() ) continue;
883
884 if( m_jp.FillNtuple ) {
885 // Not considering pile-up; pick-up the first PV
886 if( 0 == m_ntupleVars->get<unsigned int>( "NumPV" ) ) {
887 m_thePV = vertex;
888
889 m_ntupleVars->get<double>( "PVX" ) = vertex->x();
890 m_ntupleVars->get<double>( "PVY" ) = vertex->y();
891 m_ntupleVars->get<double>( "PVZ" ) = vertex->z();
892 m_ntupleVars->get<unsigned int>( "PVType" ) = vertex->vertexType();
893
894 // number of tracks associated to the PV
895 m_ntupleVars->get<unsigned int>( "NTrksPV" ) = vertex->nTrackParticles();
896 }
897 }
898 }
899 }
900
901 // if thePV is null, the PV is not found.
902 if( !m_thePV ) {
903 ATH_MSG_DEBUG("No PV is found in the PV collection!");
904 return StatusCode::FAILURE;
905 }
906
907 ATH_MSG_DEBUG(" Primary vertex successful. thePV = " << m_thePV );
908
909 return StatusCode::SUCCESS;
910 }
911
912
913 //____________________________________________________________________________________________________
914 void VrtSecInclusive::trackClassification(std::vector<WrkVrt> *workVerticesContainer, std::map<long int, std::vector<long int> >& trackToVertexMap )
915 {
916 // Fill TrkInVrt with vertex IDs of each track
917
918 trackToVertexMap.clear();
919
920 for( size_t iv = 0; iv<workVerticesContainer->size(); iv++ ) {
921
922 WrkVrt& vertex = workVerticesContainer->at(iv);
923
924 auto& trackIndices = vertex.selectedTrackIndices;
925 if( !vertex.isGood ) continue;
926 if( trackIndices.size() < 2 ) continue;
927
928 for( auto& index : trackIndices ) {
929 trackToVertexMap[index].emplace_back( iv );
930 }
931 }
932
933 for( auto& pair: trackToVertexMap ) {
934 std::string msg = Form("track index [%ld]: vertices = (", pair.first);
935 for( auto& v : pair.second ) {
936 msg += Form("%ld, ", v);
937 }
938 msg += ")";
939 if( pair.second.size() >=2 ) ATH_MSG_VERBOSE(" >> " << __FUNCTION__ << ": " << msg );
940 }
941
942 }
943
944
945 //____________________________________________________________________________________________________
946 double VrtSecInclusive::findWorstChi2ofMaximallySharedTrack(std::vector<WrkVrt> *workVerticesContainer,
947 std::map<long int, std::vector<long int> >& trackToVertexMap,
948 long int & maxSharedTrack,
949 long int & worstMatchingVertex)
950 {
951
952 double worstChi2 = AlgConsts::invalidFloat;
953
954 // Find the track index that has the largest shared vertices
955 auto maxSharedTrackToVertices = std::max_element( trackToVertexMap.begin(), trackToVertexMap.end(), []( auto& p1, auto& p2 ) { return p1.second.size() < p2.second.size(); } );
956
957 if( maxSharedTrackToVertices == trackToVertexMap.end() ) return worstChi2;
958
959 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": max-shared track index = " << maxSharedTrackToVertices->first << ", number of shared vertices = " << maxSharedTrackToVertices->second.size() );
960
961 if( maxSharedTrackToVertices->second.size() < 2 ) return worstChi2;
962
963 // map of vertex index and the chi2 of the track for the maxSharedTrack
964 std::map<long int, double> vrtChi2Map;
965
966 // loop over vertices for the largest shared track
967 for( auto& iv : maxSharedTrackToVertices->second ) {
968 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": loop over vertices: vertex index " << iv );
969
970 auto& wrkvrt = workVerticesContainer->at( iv );
971 auto& trackIndices = wrkvrt.selectedTrackIndices;
972
973 // find the index of the track
974 auto index = std::find_if( trackIndices.begin(), trackIndices.end(), [&]( auto& index ) { return index == maxSharedTrackToVertices->first; } );
975 if( index == trackIndices.end() ) {
976 ATH_MSG_WARNING(" >> " << __FUNCTION__ << ": index not found (algorithm inconsistent)" );
977 return worstChi2;
978 }
979
980 auto& chi2 = wrkvrt.Chi2PerTrk.at( index - trackIndices.begin() );
981
982 vrtChi2Map.emplace( std::pair<long int, double>(iv, chi2) );
983 }
984
985 auto worstVrtChi2Pair = std::max_element( vrtChi2Map.begin(), vrtChi2Map.end(), []( auto& p1, auto& p2 ) { return p1.second < p2.second; } );
986
987 if( worstVrtChi2Pair == vrtChi2Map.end() ) {
988 ATH_MSG_WARNING(" >> " << __FUNCTION__ << ": max_element of vrtChi2Map not found" );
989 return worstChi2;
990 }
991
992 maxSharedTrack = maxSharedTrackToVertices->first;
993 worstMatchingVertex = worstVrtChi2Pair->first;
994 worstChi2 = worstVrtChi2Pair->second;
995
996 return worstChi2;
997 }
998
999
1000 //____________________________________________________________________________________________________
1001 void VrtSecInclusive::printWrkSet(const std::vector<WrkVrt> *workVerticesContainer, const std::string& name)
1002 {
1003 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": ===============================================================" );
1004 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": " << name << ": #vertices = " << workVerticesContainer->size() );
1005
1006 std::set<const xAOD::TrackParticle*> usedTracks;
1007
1008 auto concatenateIndicesToString = []( auto indices, auto& collection ) -> std::string {
1009 if( 0 == indices.size() ) return "";
1010 return std::accumulate( std::next(indices.begin()), indices.end(), std::to_string( indices.at(0) ),
1011 [&collection]( const std::string& str, auto& index ) { return str + ", " + std::to_string( collection.at(index)->index() ); } );
1012 };
1013
1014 std::map<const xAOD::TruthVertex*, bool> previous;
1015
1016 for( auto& pair : m_matchMap ) { previous.emplace( pair.first, pair.second ); }
1017
1018 m_matchMap.clear();
1019 for( const auto* truthVertex : m_tracingTruthVertices ) { m_matchMap.emplace( truthVertex, false ); }
1020
1021 for(size_t iv=0; iv<workVerticesContainer->size(); iv++) {
1022 const auto& wrkvrt = workVerticesContainer->at(iv);
1023
1024 if( wrkvrt.nTracksTotal() < 2 ) continue;
1025
1026 std::string sels = concatenateIndicesToString( wrkvrt.selectedTrackIndices, m_selectedTracks );
1027 std::string assocs = concatenateIndicesToString( wrkvrt.associatedTrackIndices, m_associatedTracks );
1028
1029 for( const auto& index : wrkvrt.selectedTrackIndices ) { usedTracks.insert( m_selectedTracks.at(index) ); }
1030 for( const auto& index : wrkvrt.associatedTrackIndices ) { usedTracks.insert( m_associatedTracks.at(index) ); }
1031
1032 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": " << name << " vertex [" << iv << "]: " << &wrkvrt
1033 << ", isGood = " << (wrkvrt.isGood? "true" : "false")
1034 << ", #ntrks(tot, sel, assoc) = (" << wrkvrt.nTracksTotal() << ", " << wrkvrt.selectedTrackIndices.size() << ", " << wrkvrt.associatedTrackIndices.size() << "), "
1035 << ", chi2/ndof = " << wrkvrt.fitQuality()
1036 << ", (r, z) = (" << wrkvrt.vertex.perp()
1037 << ", " << wrkvrt.vertex.z() << ")"
1038 << ", sels = { " << sels << " }"
1039 << ", assocs = { " << assocs << " }" );
1040
1041 // Truth match condition
1042 for( const auto* truthVertex : m_tracingTruthVertices ) {
1043
1044
1045 Amg::Vector3D vTruth( truthVertex->x(), truthVertex->y(), truthVertex->z() );
1046 Amg::Vector3D vReco ( wrkvrt.vertex.x(), wrkvrt.vertex.y(), wrkvrt.vertex.z() );
1047
1048 const auto distance = vReco - vTruth;
1049
1050 AmgSymMatrix(3) cov;
1051 cov.fillSymmetric( 0, 0, wrkvrt.vertexCov.at(0) );
1052 cov.fillSymmetric( 1, 0, wrkvrt.vertexCov.at(1) );
1053 cov.fillSymmetric( 1, 1, wrkvrt.vertexCov.at(2) );
1054 cov.fillSymmetric( 2, 0, wrkvrt.vertexCov.at(3) );
1055 cov.fillSymmetric( 2, 1, wrkvrt.vertexCov.at(4) );
1056 cov.fillSymmetric( 2, 2, wrkvrt.vertexCov.at(5) );
1057
1058 const double s2 = distance.transpose() * cov.inverse() * distance;
1059
1060 if( distance.norm() < 2.0 || s2 < 100. ) m_matchMap.at( truthVertex ) = true;
1061
1062 }
1063
1064 }
1065
1066 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": number of used tracks = " << usedTracks.size() );
1067
1068 if( !previous.empty() && previous.size() == m_matchMap.size() ) {
1069 for( auto& pair : m_matchMap ) {
1070 if( previous.find( pair.first ) == previous.end() ) continue;
1071 if( pair.second != previous.at( pair.first ) ) {
1072 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": Match flag has changed: (r, z) = (" << pair.first->perp() << ", " << pair.first->z() << ")" );
1073 }
1074 }
1075 }
1076
1077 if( m_jp.FillHist ) {
1078 for( auto& pair : m_matchMap ) {
1079 if( pair.second ) m_hists["nMatchedTruths"]->Fill( m_vertexingAlgorithmStep+2, pair.first->perp() );
1080 }
1081 }
1082
1083 std::string msg;
1084 for( const auto* trk : usedTracks ) { msg += Form("%ld, ", trk->index() ); }
1085
1086 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": used tracks = " << msg );
1087 ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": ===============================================================" );
1088
1089 }
1090
1091
1092 //____________________________________________________________________________________________________
1094 summary.numIBLHits = 0;
1095 summary.numBLayerHits = 0;
1096 summary.numPixelHits = 0;
1097 summary.numSctHits = 0;
1098 summary.numTrtHits = 0;
1099
1100 trk->summaryValue( summary.numIBLHits, xAOD::numberOfInnermostPixelLayerHits );
1101 trk->summaryValue( summary.numBLayerHits, xAOD::numberOfNextToInnermostPixelLayerHits );
1102 trk->summaryValue( summary.numPixelHits, xAOD::numberOfPixelHits );
1103 trk->summaryValue( summary.numSctHits, xAOD::numberOfSCTHits );
1104 trk->summaryValue( summary.numTrtHits, xAOD::numberOfTRTHits );
1105 }
1106
1107
1108
1109 //____________________________________________________________________________________________________
1111
1112 auto* pattern = new ExtrapolatedPattern;
1113 const EventContext& ctx = Gaudi::Hive::currentContext();
1114 std::vector<std::unique_ptr<Trk::TrackParameters>> paramsVector =
1115 m_extrapolator->extrapolateBlindly(ctx, trk->perigeeParameters(), direction);
1116
1118
1119 auto nDisabled = 0;
1120
1121 for( auto& params : paramsVector ) {
1122
1123 const TVector3 position( params->position().x(), params->position().y(), params->position().z() );
1124
1125 if( prevPos == position ) {
1126 continue;
1127 }
1128
1129 prevPos = position;
1130
1131 const auto* detElement = params->associatedSurface().associatedDetectorElement();
1132
1133 if( detElement ) {
1134
1135 enum { Pixel = 1, SCT = 2 };
1136
1137 const auto& id = detElement->identify();
1138 Flag good = false;
1139
1140 if( m_atlasId->is_pixel(id) ) {
1141
1142 auto idHash = m_pixelId->wafer_hash( id );
1143 good = m_pixelCondSummaryTool->isGood( idHash, ctx );
1144
1145 pattern->emplace_back( std::make_tuple( position, Pixel, m_pixelId->barrel_ec(id), m_pixelId->layer_disk(id), good ) );
1146
1147 } else if( m_atlasId->is_sct(id) ) {
1148
1149 auto idHash = m_sctId->wafer_hash( id );
1150 good = m_sctCondSummaryTool->isGood( idHash, ctx );
1151
1152 pattern->emplace_back( std::make_tuple( position, SCT, m_sctId->barrel_ec(id), m_sctId->layer_disk(id), good ) );
1153
1154 }
1155
1156 if( !pattern->empty() ) {
1157
1158 ATH_MSG_VERBOSE(" >> " << __FUNCTION__ << ", track " << trk << ": position = (" << position.Perp() << ", " << position.z() << ", " << position.Phi() << "), detElement ID = " << id << ", good = " << good
1159 << ": (det, bec, layer) = (" << std::get<1>( pattern->back() ) << ", " << std::get<2>( pattern->back() ) << ", " << std::get<3>( pattern->back() ) << ")" );
1160
1161 if( !good ) nDisabled++;
1162 }
1163
1164 }
1165
1166 }
1167
1168 if( m_jp.FillHist ) {
1169 m_hists["disabledCount"]->Fill( nDisabled );
1170 }
1171
1172 return pattern;
1173
1174 }
1175
1176
1177 //____________________________________________________________________________________________________
1179 {
1180
1181 if( m_extrapolatedPatternBank.find( trk ) == m_extrapolatedPatternBank.end() ) {
1182
1183 std::unique_ptr<ExtrapolatedPattern> exPattern_along( extrapolatedPattern( trk, Trk::alongMomentum ) );
1184 std::unique_ptr<ExtrapolatedPattern> exPattern_oppos( extrapolatedPattern( trk, Trk::oppositeMomentum ) );
1185
1186 m_extrapolatedPatternBank.emplace( trk, std::make_pair( std::move(exPattern_along), std::move(exPattern_oppos) ) );
1187
1188 }
1189
1190 auto& exPattern = m_extrapolatedPatternBank.at( trk );
1191
1192 using LayerCombination = std::vector<int>;
1193
1194 std::map<LayerCombination, unsigned> layerMap;
1195 if( layerMap.empty() ) {
1196 layerMap[ { 1, 0, 0 } ] = Trk::pixelBarrel0;
1197 layerMap[ { 1, 0, 1 } ] = Trk::pixelBarrel1;
1198 layerMap[ { 1, 0, 2 } ] = Trk::pixelBarrel2;
1199 layerMap[ { 1, 0, 3 } ] = Trk::pixelBarrel3;
1200
1201 layerMap[ { 1, 2, 0 } ] = Trk::pixelEndCap0;
1202 layerMap[ { 1, 2, 1 } ] = Trk::pixelEndCap1;
1203 layerMap[ { 1, 2, 2 } ] = Trk::pixelEndCap2;
1204 layerMap[ { 1,-2, 0 } ] = Trk::pixelEndCap0;
1205 layerMap[ { 1,-2, 1 } ] = Trk::pixelEndCap1;
1206 layerMap[ { 1,-2, 2 } ] = Trk::pixelEndCap2;
1207
1208 layerMap[ { 2, 0, 0 } ] = Trk::sctBarrel0;
1209 layerMap[ { 2, 0, 1 } ] = Trk::sctBarrel1;
1210 layerMap[ { 2, 0, 2 } ] = Trk::sctBarrel2;
1211 layerMap[ { 2, 0, 3 } ] = Trk::sctBarrel3;
1212
1213 layerMap[ { 2, 2, 0 } ] = Trk::sctEndCap0;
1214 layerMap[ { 2, 2, 1 } ] = Trk::sctEndCap1;
1215 layerMap[ { 2, 2, 2 } ] = Trk::sctEndCap2;
1216 layerMap[ { 2, 2, 3 } ] = Trk::sctEndCap3;
1217 layerMap[ { 2, 2, 4 } ] = Trk::sctEndCap4;
1218 layerMap[ { 2, 2, 5 } ] = Trk::sctEndCap5;
1219 layerMap[ { 2, 2, 6 } ] = Trk::sctEndCap6;
1220 layerMap[ { 2, 2, 7 } ] = Trk::sctEndCap7;
1221 layerMap[ { 2, 2, 8 } ] = Trk::sctEndCap8;
1222 layerMap[ { 2,-2, 0 } ] = Trk::sctEndCap0;
1223 layerMap[ { 2,-2, 1 } ] = Trk::sctEndCap1;
1224 layerMap[ { 2,-2, 2 } ] = Trk::sctEndCap2;
1225 layerMap[ { 2,-2, 3 } ] = Trk::sctEndCap3;
1226 layerMap[ { 2,-2, 4 } ] = Trk::sctEndCap4;
1227 layerMap[ { 2,-2, 5 } ] = Trk::sctEndCap5;
1228 layerMap[ { 2,-2, 6 } ] = Trk::sctEndCap6;
1229 layerMap[ { 2,-2, 7 } ] = Trk::sctEndCap7;
1230 layerMap[ { 2,-2, 8 } ] = Trk::sctEndCap8;
1231 }
1232
1233 enum { position=0, detector=1, bec=2, layer=3, isGood=4 };
1234
1235 // Lambda!
1236 auto getDetectorType = [&]( const ExtrapolatedPoint& point ) -> unsigned {
1237
1238 const LayerCombination comb { std::get<detector>( point ), std::get<bec>( point ), std::get<layer>( point ) };
1239
1240 for( auto& pair : layerMap ) {
1241 if( pair.first == comb ) {
1242 return pair.second;
1243 }
1244 }
1245
1247 };
1248
1249 enum { kShouldNotHaveHit, kShouldHaveHit, kMayHaveHit };
1250 std::vector<unsigned> expectedHitPattern(Trk::numberOfDetectorTypes, kShouldNotHaveHit);
1251
1252 auto minExpectedRadius = AlgConsts::maxValue;
1253
1254 // Loop over extrapolated points (along direction)
1255 auto& exPattern_along = *( exPattern.first );
1256
1257 for( auto itr = exPattern_along.begin(); itr != exPattern_along.end(); ++itr ) {
1258 if( std::next( itr ) == exPattern_along.end() ) continue;
1259
1260 const auto& point = *itr;
1261 const auto& nextPoint = *( std::next( itr ) );
1262
1263 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": isGood = " << std::get<isGood>( point ) );
1264
1265 const auto& thisPos = std::get<position>( point );
1266 const auto& nextPos = std::get<position>( nextPoint );
1267
1268 auto sectionVector = nextPos - thisPos;
1269 auto vertexVector = TVector3( vertex.x(), vertex.y(), vertex.z() ) - thisPos;
1270
1271
1272 const auto& detectorType = getDetectorType( point );
1273
1274 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": detType = " << detectorType );
1275
1276 if( detectorType == AlgConsts::invalidUnsigned ) continue;
1277 if( detectorType >= Trk::numberOfDetectorTypes ) continue;
1278
1279 // if the vertex is nearby (within 10 mm), the hit may be presnet ("X")
1280 if( vertexVector.Mag() < 10. ) {
1281 expectedHitPattern.at( detectorType ) = kMayHaveHit;
1282 continue;
1283 }
1284
1285 // if the front-end module is not active, then the hit is not expected,
1286 // which means the hit may be present
1287 if( !static_cast<bool>(std::get<isGood>( point )) ) {
1288 expectedHitPattern.at( detectorType ) = kMayHaveHit;
1289 continue;
1290 }
1291
1292 // if the inner product of the above two vectors is positive,
1293 // then point is inner than the vertex.
1294 // Else, the point is outer than the vertex and expect to have hits
1295 // when the track is originated from the vertex.
1296
1297 if( sectionVector.Mag() == 0. ) continue;
1298
1299 ATH_MSG_VERBOSE( " > " << __FUNCTION__
1300 << ": hitPos = (" << thisPos.Perp() << ", " << thisPos.z() << ", " << thisPos.Phi() << ")"
1301 << ", sectionVec = (" << sectionVector.Perp() << ", " << sectionVector.z() << ", " << sectionVector.Phi() << ")"
1302 << ", vertexVec = (" << vertexVector.Perp() << ", " << vertexVector.z() << ", " << vertexVector.Phi() << ")"
1303 << ", cos(s,v) = " << sectionVector * vertexVector / ( sectionVector.Mag() * vertexVector.Mag() + AlgConsts::infinitesimal ) );
1304
1305 if( sectionVector * vertexVector > 0. ) continue;
1306
1307 if( minExpectedRadius > thisPos.Perp() ) minExpectedRadius = thisPos.Perp();
1308
1309 // now, the hit is expected to present.
1310
1311 expectedHitPattern.at( detectorType ) = kShouldHaveHit;
1312 }
1313
1314 // Loop over extrapolated points (opposite direction)
1315 auto& exPattern_oppos = *( exPattern.second );
1316 bool oppositeFlag = false;
1317
1318 for( auto itr = exPattern_oppos.begin(); itr != exPattern_oppos.end(); ++itr ) {
1319 if( std::next( itr ) == exPattern_oppos.end() ) continue;
1320
1321 const auto& point = *itr;
1322 const auto& nextPoint = *( std::next( itr ) );
1323
1324 const auto& thisPos = std::get<position>( point );
1325 const auto& nextPos = std::get<position>( nextPoint );
1326
1327 auto sectionVector = nextPos - thisPos;
1328 auto vertexVector = TVector3( vertex.x(), vertex.y(), vertex.z() ) - thisPos;
1329
1330 const auto& detectorType = getDetectorType( point );
1331
1332 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": detType = " << detectorType );
1333
1334 ATH_MSG_DEBUG( " > " << __FUNCTION__
1335 << ": hitPos = (" << thisPos.Perp() << ", " << thisPos.z() << ", " << thisPos.Phi() << ")"
1336 << ", vertex = (" << vertex.perp() << ", " << vertex.z() << ", " << vertex.phi() << ")"
1337 << ", cos(s,v) = " << sectionVector * vertexVector / ( sectionVector.Mag() * vertexVector.Mag() + AlgConsts::infinitesimal ) );
1338
1339 if( detectorType == AlgConsts::invalidUnsigned ) continue;
1340 if( detectorType >= Trk::numberOfDetectorTypes ) continue;
1341
1342 if( sectionVector * vertexVector < 0. ) {
1343 oppositeFlag = true;
1344 }
1345 }
1346
1347 // If the first expected point's radius is smaller than the vertex radius,
1348 // it's the case that the vertex was reconstructed in the opposite phi-direction
1349 // to the track outgoing direction. Such a case should be rejected.
1350 // bool oppositeFlag = ( minExpectedRadius < vertex.perp() );
1351
1352 std::string msg = "Expected hit pattern: ";
1353 for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1354 msg += Form("%s", expectedHitPattern.at(i) < kMayHaveHit? Form("%u", expectedHitPattern.at(i)) : "X" );
1355 }
1356 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
1357
1358 msg = "Recorded hit pattern: ";
1359 for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1360 msg += Form("%u", ( trk->hitPattern() >> i ) & 1 );
1361 }
1362 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
1363
1364 std::vector<unsigned> matchedLayers;
1365
1366 for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1367 const unsigned recordedPattern = ( (trk->hitPattern()>>i) & 1 );
1368
1369 if( expectedHitPattern.at(i) == kMayHaveHit ) {
1370 matchedLayers.emplace_back( i );
1371 } else if( expectedHitPattern.at(i) == kShouldHaveHit ) {
1372 if( expectedHitPattern.at(i) != recordedPattern ) {
1373 break;
1374 } else {
1375 matchedLayers.emplace_back( i );
1376 }
1377 } else {
1378 if( expectedHitPattern.at(i) != recordedPattern ) {
1379 break;
1380 }
1381 }
1382
1383 }
1384
1385 uint8_t PixelHits = 0;
1386 uint8_t SctHits = 0;
1387 uint8_t TRTHits = 0;
1388 if( !(trk->summaryValue( PixelHits, xAOD::numberOfPixelHits ) ) ) PixelHits =0;
1389 if( !(trk->summaryValue( SctHits, xAOD::numberOfSCTHits ) ) ) SctHits =0;
1390 if( !(trk->summaryValue( TRTHits, xAOD::numberOfTRTHits ) ) ) TRTHits =0;
1391
1392 auto dphi = trk->phi() - vertex.phi();
1393 while( dphi > TMath::Pi() ) dphi -= TMath::TwoPi();
1394 while( dphi < -TMath::Pi() ) dphi += TMath::TwoPi();
1395
1396 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": vtx phi = " << vertex.phi() << ", track phi = " << trk->phi() << ", dphi = " << dphi
1397 << ", oppositeFlag = " << oppositeFlag
1398 << ", nPixelHits = " << static_cast<int>(PixelHits)
1399 << ", nSCTHits = " << static_cast<int>(SctHits)
1400 << ", nTRTHits = " << static_cast<int>(TRTHits)
1401 << ", nMatchedLayers = " << matchedLayers.size() );
1402
1403 if( PixelHits == 0 && vertex.perp() > 300. ) {
1404 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": vertex r > 300 mm, w/o no pixel hits" );
1405 }
1406
1407
1408 // Requires the first 2 layers with the hit matches.
1409 if( matchedLayers.size() < 2 ) return false;
1410
1411 // In case of the first matched layer is not within pixel barrel, requires the first 4 layers with the hit match
1412 if( matchedLayers.at(0) >= Trk::pixelEndCap0 ) {
1413 if( matchedLayers.size() < 4 ) return false;
1414 }
1415
1416 // Sometimes the vertex is reconstructed at the opposite phi direction.
1417 // In this case, the pattern match may pass.
1418 // This can be avoided by requiring that the
1419 if( oppositeFlag ) return false;
1420
1421 // The following condition should apply for vertices outer than IBL.
1422 if( false /*matchedLayers.at(0) > Trk::pixelBarrel0*/ ) {
1423
1424 // If the dphi (defined above) is opposite, reject.
1425 if( fabs( dphi ) > TMath::Pi()/2.0 ) return false;
1426
1427 // If the track is not within the forward hemisphere to the vertex, reject.
1428 TVector3 trkP; trkP.SetPtEtaPhi( trk->pt(), trk->eta(), trk->phi() );
1429 TVector3 vtx; vtx.SetXYZ( vertex.x(), vertex.y(), vertex.z() );
1430 if( trkP.Dot( vtx ) < 0. ) return false;
1431
1432 }
1433
1434 return true;
1435 }
1436
1437
1438 //____________________________________________________________________________________________________
1439 bool VrtSecInclusive::patternCheckRun2( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1440
1441 //
1442 // rough guesses for active layers:
1443 // BeamPipe: 23.5-24.3
1444 // IBL: 31.0-38.4
1445 // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1446 // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1447 //
1448
1449 const double rad = vertex.perp();
1450 const double absz = fabs( vertex.z() );
1451
1452 // vertex area classification
1453 enum vertexArea {
1454 insideBeamPipe,
1455
1456 insidePixelBarrel0,
1457 aroundPixelBarrel0,
1458
1459 outsidePixelBarrel0_and_insidePixelBarrel1,
1460 aroundPixelBarrel1,
1461
1462 outsidePixelBarrel1_and_insidePixelBarrel2,
1463 aroundPixelBarrel2,
1464
1465 outsidePixelBarrel2_and_insidePixelBarrel3,
1466 aroundPixelBarrel3,
1467
1468 outsidePixelBarrel3_and_insideSctBarrel0,
1469 aroundSctBarrel0,
1470
1471 outsideSctBarrel0_and_insideSctBarrel1,
1472 aroundSctBarrel1,
1473 };
1474
1475 // Mutually exclusive vertex position pattern
1476 int vertex_pattern = 0;
1477 if( rad < 23.50 ) {
1478 vertex_pattern = insideBeamPipe;
1479
1480 } else if( rad < 31.0 && absz < 331.5 ) {
1481 vertex_pattern = insidePixelBarrel0;
1482
1483 } else if( rad < 38.4 && absz < 331.5 ) {
1484 vertex_pattern = aroundPixelBarrel0;
1485
1486 } else if( rad < 47.7 && absz < 400.5 ) {
1487 vertex_pattern = outsidePixelBarrel0_and_insidePixelBarrel1;
1488
1489 } else if( rad < 54.4 && absz < 400.5 ) {
1490 vertex_pattern = aroundPixelBarrel1;
1491
1492 } else if( rad < 85.5 && absz < 400.5 ) {
1493 vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1494
1495 } else if( rad < 92.2 && absz < 400.5 ) {
1496 vertex_pattern = aroundPixelBarrel2;
1497
1498 } else if( rad < 119.3 && absz < 400.5 ) {
1499 vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1500
1501 } else if( rad < 126.1 && absz < 400.5 ) {
1502 vertex_pattern = aroundPixelBarrel3;
1503
1504 } else if( rad < 290 && absz < 749.0 ) {
1505 vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1506
1507 } else if( rad < 315 && absz < 749.0 ) {
1508 vertex_pattern = aroundSctBarrel0;
1509
1510 } else if( rad < 360 && absz < 749.0 ) {
1511 vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1512
1513 } else if( rad < 390 && absz < 749.0 ) {
1514 vertex_pattern = aroundSctBarrel1;
1515
1516 } else {
1517 }
1518
1519 unsigned nPixelLayers { 0 };
1520 {
1521 if ( pattern & (1 << Trk::pixelBarrel0) ) nPixelLayers++;
1522 if ( pattern & (1 << Trk::pixelBarrel1) ) nPixelLayers++;
1523 if ( pattern & (1 << Trk::pixelBarrel2) ) nPixelLayers++;
1524 if ( pattern & (1 << Trk::pixelBarrel3) ) nPixelLayers++;
1525 if ( pattern & (1 << Trk::pixelEndCap0) ) nPixelLayers++;
1526 if ( pattern & (1 << Trk::pixelEndCap1) ) nPixelLayers++;
1527 if ( pattern & (1 << Trk::pixelEndCap2) ) nPixelLayers++;
1528 }
1529
1531 if( vertex_pattern == insideBeamPipe ) {
1532
1533 if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1534 if( nPixelLayers < 3 ) return false;
1535
1536
1537 } else if( vertex_pattern == insidePixelBarrel0 ) {
1538
1539 if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1540 if( nPixelLayers < 3 ) return false;
1541 }
1542
1543
1544 else if( vertex_pattern == aroundPixelBarrel0 ) {
1545
1546 // require nothing for PixelBarrel0
1547 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1548 if( nPixelLayers < 2 ) return false;
1549 }
1550
1551
1552 else if( vertex_pattern == outsidePixelBarrel0_and_insidePixelBarrel1 ) {
1553
1554 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1555 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1556 if( nPixelLayers < 2 ) return false;
1557 }
1558
1559
1560 else if( vertex_pattern == aroundPixelBarrel1 ) {
1561
1562 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1563 // require nothing for PixelBarrel
1564 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1565 if( nPixelLayers < 2 ) return false;
1566 }
1567
1568
1569 else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1570
1571 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1572 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1573 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1574 if( nPixelLayers < 2 ) return false;
1575 }
1576
1577
1578 else if( vertex_pattern == aroundPixelBarrel2 ) {
1579
1580 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1581 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1582 // require nothing for PixelBarrel2
1583 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1584 }
1585
1586
1587 else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1588
1589 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1590 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1591 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1592 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1593 }
1594
1595 else if( vertex_pattern == aroundPixelBarrel3 ) {
1596
1597 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1598 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1599 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1600 // require nothing for PixelBarrel3
1601 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1602 }
1603
1604
1605 else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1606
1607 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1608 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1609 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1610 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1611 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1612 }
1613
1614
1615 else if( vertex_pattern == aroundSctBarrel0 ) {
1616
1617 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1618 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1619 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1620 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1621 // require nothing for SctBarrel0
1622 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1623 }
1624
1625
1626 else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1627
1628 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1629 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1630 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1631 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1632 if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1633 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1634 }
1635
1636
1637 else if( vertex_pattern == aroundSctBarrel1 ) {
1638 if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1639 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1640 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1641 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1642 if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1643 // require nothing for SctBarrel1
1644 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1645 }
1647
1648 return true;
1649
1650 }
1651
1652 //____________________________________________________________________________________________________
1653 bool VrtSecInclusive::patternCheckRun2OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1654
1655 //
1656 // rough guesses for active layers:
1657 // BeamPipe: 23.5-24.3
1658 // IBL: 31.0-38.4
1659 // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1660 // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1661 //
1662
1663 const double rad = vertex.perp();
1664 const double absz = fabs( vertex.z() );
1665
1666 // vertex area classification
1667 enum vertexArea {
1668 insideBeamPipe,
1669
1670 insidePixelBarrel0,
1671 aroundPixelBarrel0,
1672
1673 outsidePixelBarrel0_and_insidePixelBarrel1,
1674 aroundPixelBarrel1,
1675
1676 outsidePixelBarrel1_and_insidePixelBarrel2,
1677 aroundPixelBarrel2,
1678
1679 outsidePixelBarrel2_and_insidePixelBarrel3,
1680 aroundPixelBarrel3,
1681
1682 outsidePixelBarrel3_and_insideSctBarrel0,
1683 aroundSctBarrel0,
1684
1685 outsideSctBarrel0_and_insideSctBarrel1,
1686 aroundSctBarrel1,
1687 };
1688
1689 // Mutually exclusive vertex position pattern
1690 int vertex_pattern = 0;
1691 if( rad < 23.50 ) {
1692 vertex_pattern = insideBeamPipe;
1693
1694 } else if( rad < 31.0 && absz < 331.5 ) {
1695 vertex_pattern = insidePixelBarrel0;
1696
1697 } else if( rad < 38.4 && absz < 331.5 ) {
1698 vertex_pattern = aroundPixelBarrel0;
1699
1700 } else if( rad < 47.7 && absz < 400.5 ) {
1701 vertex_pattern = outsidePixelBarrel0_and_insidePixelBarrel1;
1702
1703 } else if( rad < 54.4 && absz < 400.5 ) {
1704 vertex_pattern = aroundPixelBarrel1;
1705
1706 } else if( rad < 85.5 && absz < 400.5 ) {
1707 vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1708
1709 } else if( rad < 92.2 && absz < 400.5 ) {
1710 vertex_pattern = aroundPixelBarrel2;
1711
1712 } else if( rad < 119.3 && absz < 400.5 ) {
1713 vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1714
1715 } else if( rad < 126.1 && absz < 400.5 ) {
1716 vertex_pattern = aroundPixelBarrel3;
1717
1718 } else if( rad < 290 && absz < 749.0 ) {
1719 vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1720
1721 } else if( rad < 315 && absz < 749.0 ) {
1722 vertex_pattern = aroundSctBarrel0;
1723
1724 } else if( rad < 360 && absz < 749.0 ) {
1725 vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1726
1727 } else if( rad < 390 && absz < 749.0 ) {
1728 vertex_pattern = aroundSctBarrel1;
1729
1730 } else {
1731 }
1732
1733
1734 unsigned nPixelLayers { 0 };
1735 {
1736 nPixelLayers += ( pattern & (1 << Trk::pixelBarrel0) );
1737 nPixelLayers += ( pattern & (1 << Trk::pixelBarrel1) );
1738 nPixelLayers += ( pattern & (1 << Trk::pixelBarrel2) );
1739 nPixelLayers += ( pattern & (1 << Trk::pixelBarrel3) );
1740 nPixelLayers += ( pattern & (1 << Trk::pixelEndCap0) );
1741 nPixelLayers += ( pattern & (1 << Trk::pixelEndCap1) );
1742 nPixelLayers += ( pattern & (1 << Trk::pixelEndCap2) );
1743 }
1744
1746 if( vertex_pattern == insideBeamPipe ) {
1747
1748 if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1749 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1750 if( nPixelLayers < 3 ) return false;
1751
1752 } else if( vertex_pattern == insidePixelBarrel0 ) {
1753
1754 if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1755 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1756 if( nPixelLayers < 3 ) return false;
1757
1758 }
1759
1760
1761 else if( vertex_pattern == aroundPixelBarrel0 ) {
1762
1763 // require nothing for PixelBarrel0
1764 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1765 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1766 if( nPixelLayers < 3 ) return false;
1767 }
1768
1769
1770 else if( vertex_pattern == outsidePixelBarrel0_and_insidePixelBarrel1 ) {
1771
1772 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1773 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1774 if( nPixelLayers < 3 ) return false;
1775 }
1776
1777
1778 else if( vertex_pattern == aroundPixelBarrel1 ) {
1779
1780 // require nothing for PixelBarrel1
1781 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1782 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1783 if( nPixelLayers < 2 ) return false;
1784 }
1785
1786
1787 else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1788
1789 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1790 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1791 if( nPixelLayers < 2 ) return false;
1792 }
1793
1794
1795 else if( vertex_pattern == aroundPixelBarrel2 ) {
1796
1797 // require nothing for PixelBarrel2
1798 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1799 }
1800
1801
1802 else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1803
1804 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1805 }
1806
1807 else if( vertex_pattern == aroundPixelBarrel3 ) {
1808
1809 // require nothing for PixelBarrel3
1810 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1811 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1812 }
1813
1814
1815 else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1816
1817 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1818 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1819 }
1820
1821
1822 else if( vertex_pattern == aroundSctBarrel0 ) {
1823
1824 // require nothing for SctBarrel0
1825 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1826 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1827 }
1828
1829
1830 else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1831
1832 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1833 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1834 }
1835
1836
1837 else if( vertex_pattern == aroundSctBarrel1 ) {
1838 // require nothing for SctBarrel1
1839 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1840 if( ! (pattern & (1<<Trk::sctBarrel3)) ) return false;
1841 }
1843
1844 return true;
1845
1846 }
1847
1848 //____________________________________________________________________________________________________
1849 bool VrtSecInclusive::patternCheckRun1( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1850 //
1851 // rough guesses for active layers:
1852 // BeamPipe: 25.0
1853 // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1854 // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1855 //
1856
1857 const double rad = vertex.perp();
1858 const double absz = fabs( vertex.z() );
1859
1860 // vertex area classification
1861 enum vertexArea {
1862 insideBeamPipe,
1863
1864 insidePixelBarrel1,
1865 aroundPixelBarrel1,
1866
1867 outsidePixelBarrel1_and_insidePixelBarrel2,
1868 aroundPixelBarrel2,
1869
1870 outsidePixelBarrel2_and_insidePixelBarrel3,
1871 aroundPixelBarrel3,
1872
1873 outsidePixelBarrel3_and_insideSctBarrel0,
1874 aroundSctBarrel0,
1875
1876 outsideSctBarrel0_and_insideSctBarrel1,
1877 aroundSctBarrel1,
1878 };
1879
1880 // Mutually exclusive vertex position pattern
1881 Int_t vertex_pattern = 0;
1882 if( rad < 25.00 ) {
1883 vertex_pattern = insideBeamPipe;
1884
1885 } else if( rad < 47.7 && absz < 400.5 ) {
1886 vertex_pattern = insidePixelBarrel1;
1887
1888 } else if( rad < 54.4 && absz < 400.5 ) {
1889 vertex_pattern = aroundPixelBarrel1;
1890
1891 } else if( rad < 85.5 && absz < 400.5 ) {
1892 vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1893
1894 } else if( rad < 92.2 && absz < 400.5 ) {
1895 vertex_pattern = aroundPixelBarrel2;
1896
1897 } else if( rad < 119.3 && absz < 400.5 ) {
1898 vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1899
1900 } else if( rad < 126.1 && absz < 400.5 ) {
1901 vertex_pattern = aroundPixelBarrel3;
1902
1903 } else if( rad < 290 && absz < 749.0 ) {
1904 vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1905
1906 } else if( rad < 315 && absz < 749.0 ) {
1907 vertex_pattern = aroundSctBarrel0;
1908
1909 } else if( rad < 360 && absz < 749.0 ) {
1910 vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1911
1912 } else if( rad < 390 && absz < 749.0 ) {
1913 vertex_pattern = aroundSctBarrel1;
1914
1915 } else {
1916 }
1917
1918
1920 if( vertex_pattern == insideBeamPipe ) {
1921
1922 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1923
1924 }
1925
1926
1927 else if( vertex_pattern == insidePixelBarrel1 ) {
1928
1929 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1930 }
1931
1932
1933 else if( vertex_pattern == aroundPixelBarrel1 ) {
1934
1935 // require nothing for PixelBarrel1
1936 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1937 }
1938
1939
1940 else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1941
1942 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1943 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1944 }
1945
1946
1947 else if( vertex_pattern == aroundPixelBarrel2 ) {
1948
1949 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1950 // require nothing for PixelBarrel2
1951 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1952 }
1953
1954
1955 else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1956
1957 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1958 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1959 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1960 }
1961
1962 else if( vertex_pattern == aroundPixelBarrel3 ) {
1963
1964 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1965 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1966 // require nothing for PixelBarrel3
1967 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1968 }
1969
1970
1971 else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1972
1973 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1974 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1975 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1976 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1977 }
1978
1979
1980 else if( vertex_pattern == aroundSctBarrel0 ) {
1981
1982 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1983 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1984 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1985 // require nothing for SctBarrel0
1986 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1987 }
1988
1989
1990 else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1991
1992 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1993 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1994 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1995 if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1996 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1997 }
1998
1999
2000 else if( vertex_pattern == aroundSctBarrel1 ) {
2001 if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
2002 if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
2003 if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
2004 if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
2005 // require nothing for SctBarrel1
2006 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
2007 }
2009
2010 return true;
2011 }
2012
2013 //____________________________________________________________________________________________________
2014 bool VrtSecInclusive::patternCheckRun1OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2015 //
2016 // rough guesses for active layers:
2017 // BeamPipe: 25.0
2018 // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
2019 // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
2020 //
2021
2022 const double rad = vertex.perp();
2023 const double absz = fabs( vertex.z() );
2024
2025 // vertex area classification
2026 enum vertexArea {
2027 insideBeamPipe,
2028
2029 insidePixelBarrel1,
2030 aroundPixelBarrel1,
2031
2032 outsidePixelBarrel1_and_insidePixelBarrel2,
2033 aroundPixelBarrel2,
2034
2035 outsidePixelBarrel2_and_insidePixelBarrel3,
2036 aroundPixelBarrel3,
2037
2038 outsidePixelBarrel3_and_insideSctBarrel0,
2039 aroundSctBarrel0,
2040
2041 outsideSctBarrel0_and_insideSctBarrel1,
2042 aroundSctBarrel1,
2043 };
2044
2045 // Mutually exclusive vertex position pattern
2046 Int_t vertex_pattern = 0;
2047 if( rad < 25.00 ) {
2048 vertex_pattern = insideBeamPipe;
2049
2050 } else if( rad < 47.7 && absz < 400.5 ) {
2051 vertex_pattern = insidePixelBarrel1;
2052
2053 } else if( rad < 54.4 && absz < 400.5 ) {
2054 vertex_pattern = aroundPixelBarrel1;
2055
2056 } else if( rad < 85.5 && absz < 400.5 ) {
2057 vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
2058
2059 } else if( rad < 92.2 && absz < 400.5 ) {
2060 vertex_pattern = aroundPixelBarrel2;
2061
2062 } else if( rad < 119.3 && absz < 400.5 ) {
2063 vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
2064
2065 } else if( rad < 126.1 && absz < 400.5 ) {
2066 vertex_pattern = aroundPixelBarrel3;
2067
2068 } else if( rad < 290 && absz < 749.0 ) {
2069 vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
2070
2071 } else if( rad < 315 && absz < 749.0 ) {
2072 vertex_pattern = aroundSctBarrel0;
2073
2074 } else if( rad < 360 && absz < 749.0 ) {
2075 vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
2076
2077 } else if( rad < 390 && absz < 749.0 ) {
2078 vertex_pattern = aroundSctBarrel1;
2079
2080 } else {
2081 }
2082
2083
2085 if( vertex_pattern == insideBeamPipe ) {
2086
2087 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
2088
2089 }
2090
2091
2092 else if( vertex_pattern == insidePixelBarrel1 ) {
2093
2094 if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
2095 }
2096
2097
2098 else if( vertex_pattern == aroundPixelBarrel1 ) {
2099
2100 // require nothing for PixelBarrel1
2101 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
2102 }
2103
2104
2105 else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
2106
2107 if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
2108 }
2109
2110
2111 else if( vertex_pattern == aroundPixelBarrel2 ) {
2112
2113 // require nothing for PixelBarrel2
2114 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
2115 }
2116
2117
2118 else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
2119
2120 if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
2121 }
2122
2123 else if( vertex_pattern == aroundPixelBarrel3 ) {
2124
2125 // require nothing for PixelBarrel3
2126 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
2127 }
2128
2129
2130 else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
2131
2132 if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
2133 }
2134
2135
2136 else if( vertex_pattern == aroundSctBarrel0 ) {
2137
2138 // require nothing for SctBarrel0
2139 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
2140 }
2141
2142
2143 else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
2144
2145 if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
2146 }
2147
2148
2149 else if( vertex_pattern == aroundSctBarrel1 ) {
2150 // require nothing for SctBarrel1
2151 if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
2152 }
2154
2155 return true;
2156 }
2157
2158 //____________________________________________________________________________________________________
2159 bool VrtSecInclusive::patternCheck( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2160 bool flag = false;
2161
2162 if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run2 ) {
2163 flag = patternCheckRun2( pattern, vertex );
2164 } else if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run1 ) {
2165 flag = patternCheckRun1( pattern, vertex );
2166 }
2167
2168 return flag;
2169 }
2170
2171 //____________________________________________________________________________________________________
2172 bool VrtSecInclusive::patternCheckOuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2173 bool flag = false;
2174
2175 if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run2 ) {
2176 flag = patternCheckRun2OuterOnly( pattern, vertex );
2177 } else if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run1 ) {
2178 flag = patternCheckRun1OuterOnly( pattern, vertex );
2179 }
2180
2181 return flag;
2182 }
2183
2184 //____________________________________________________________________________________________________
2186 {
2187
2188 const uint32_t pattern = trk->hitPattern();
2189
2190 return patternCheck( pattern, vertex );
2191
2192 }
2193
2194
2195 //____________________________________________________________________________________________________
2197 {
2198
2199 const uint32_t pattern = trk->hitPattern();
2200
2201 return patternCheckOuterOnly( pattern, vertex );
2202
2203 }
2204
2205
2206 //____________________________________________________________________________________________________
2208 {
2209
2210 if( m_extrapolatedPatternBank.find( trk ) == m_extrapolatedPatternBank.end() ) {
2211
2212 std::unique_ptr<ExtrapolatedPattern> exPattern_along( extrapolatedPattern( trk, Trk::alongMomentum ) );
2213
2214 m_extrapolatedPatternBank.emplace( trk, std::make_pair( std::move(exPattern_along), nullptr ) );
2215
2216 }
2217
2218 if( vertex.perp() < 31.0 ) {
2219 double dphi = trk->phi() - vertex.phi();
2220 while( dphi > TMath::Pi() ) { dphi -= TMath::TwoPi(); }
2221 while( dphi < -TMath::Pi() ) { dphi += TMath::TwoPi(); }
2222 if( cos(dphi) < -0.8 ) return false;
2223 }
2224
2225 auto& exPattern = m_extrapolatedPatternBank.at( trk );
2226
2227 using LayerCombination = std::vector<int>;
2228
2229 std::map<LayerCombination, unsigned> layerMap;
2230 if( layerMap.empty() ) {
2231 layerMap[ { 1, 0, 0 } ] = Trk::pixelBarrel0;
2232 layerMap[ { 1, 0, 1 } ] = Trk::pixelBarrel1;
2233 layerMap[ { 1, 0, 2 } ] = Trk::pixelBarrel2;
2234 layerMap[ { 1, 0, 3 } ] = Trk::pixelBarrel3;
2235
2236 layerMap[ { 1, 2, 0 } ] = Trk::pixelEndCap0;
2237 layerMap[ { 1, 2, 1 } ] = Trk::pixelEndCap1;
2238 layerMap[ { 1, 2, 2 } ] = Trk::pixelEndCap2;
2239 layerMap[ { 1,-2, 0 } ] = Trk::pixelEndCap0;
2240 layerMap[ { 1,-2, 1 } ] = Trk::pixelEndCap1;
2241 layerMap[ { 1,-2, 2 } ] = Trk::pixelEndCap2;
2242
2243 layerMap[ { 2, 0, 0 } ] = Trk::sctBarrel0;
2244 layerMap[ { 2, 0, 1 } ] = Trk::sctBarrel1;
2245 layerMap[ { 2, 0, 2 } ] = Trk::sctBarrel2;
2246 layerMap[ { 2, 0, 3 } ] = Trk::sctBarrel3;
2247
2248 layerMap[ { 2, 2, 0 } ] = Trk::sctEndCap0;
2249 layerMap[ { 2, 2, 1 } ] = Trk::sctEndCap1;
2250 layerMap[ { 2, 2, 2 } ] = Trk::sctEndCap2;
2251 layerMap[ { 2, 2, 3 } ] = Trk::sctEndCap3;
2252 layerMap[ { 2, 2, 4 } ] = Trk::sctEndCap4;
2253 layerMap[ { 2, 2, 5 } ] = Trk::sctEndCap5;
2254 layerMap[ { 2, 2, 6 } ] = Trk::sctEndCap6;
2255 layerMap[ { 2, 2, 7 } ] = Trk::sctEndCap7;
2256 layerMap[ { 2, 2, 8 } ] = Trk::sctEndCap8;
2257 layerMap[ { 2,-2, 0 } ] = Trk::sctEndCap0;
2258 layerMap[ { 2,-2, 1 } ] = Trk::sctEndCap1;
2259 layerMap[ { 2,-2, 2 } ] = Trk::sctEndCap2;
2260 layerMap[ { 2,-2, 3 } ] = Trk::sctEndCap3;
2261 layerMap[ { 2,-2, 4 } ] = Trk::sctEndCap4;
2262 layerMap[ { 2,-2, 5 } ] = Trk::sctEndCap5;
2263 layerMap[ { 2,-2, 6 } ] = Trk::sctEndCap6;
2264 layerMap[ { 2,-2, 7 } ] = Trk::sctEndCap7;
2265 layerMap[ { 2,-2, 8 } ] = Trk::sctEndCap8;
2266 }
2267
2268 enum { position=0, detector=1, bec=2, layer=3, isGood=4 };
2269
2270 // Lambda!
2271 auto getDetectorType = [&]( const ExtrapolatedPoint& point ) -> unsigned {
2272
2273 const LayerCombination comb { std::get<detector>( point ), std::get<bec>( point ), std::get<layer>( point ) };
2274
2275 for( auto& pair : layerMap ) {
2276 if( pair.first == comb ) {
2277 return pair.second;
2278 }
2279 }
2280
2282 };
2283
2284 uint32_t disabledPattern { 0 };
2285
2286 // Loop over extrapolated points (along direction)
2287 auto& exPattern_along = *( exPattern.first );
2288
2289 for( auto itr = exPattern_along.begin(); itr != exPattern_along.end(); ++itr ) {
2290 if( std::next( itr ) == exPattern_along.end() ) continue;
2291
2292 const auto& point = *itr;
2293
2294 ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": isGood = " << std::get<isGood>( point ) );
2295
2296 if( !std::get<isGood>( point ) ) {
2297 const auto& detectorType = getDetectorType( point );
2298 disabledPattern += (1 << detectorType);
2299 }
2300 }
2301
2302 uint32_t hitPattern = trk->hitPattern();
2303 uint32_t modifiedPattern = disabledPattern | hitPattern;
2304
2305 std::string msg = "Disabled hit pattern: ";
2306 for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
2307 msg += Form("%u", ( disabledPattern >> i ) & 1 );
2308 }
2309 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2310
2311 msg = "Recorded hit pattern: ";
2312 for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
2313 msg += Form("%u", ( hitPattern >> i ) & 1 );
2314 }
2315 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2316
2317 return patternCheck( modifiedPattern, vertex );
2318
2319 }
2320
2321
2322 //____________________________________________________________________________________________________
2324 const xAOD::TrackParticle *itrk,
2325 const xAOD::TrackParticle *jtrk )
2326 {
2327
2328 const bool& check_itrk = ( this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( itrk, FitVertex );
2329 const bool& check_jtrk = ( this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( jtrk, FitVertex );
2330
2331 return ( check_itrk && check_jtrk );
2332
2333 }
2334
2335
2336 //____________________________________________________________________________________________________
2338 {
2339
2340 const auto& vertex = wrkvrt.vertex;
2341
2342 std::map< std::deque<long int>*, std::vector<const xAOD::TrackParticle*>& > indexMap
2343 {
2345 };
2346
2347 for( auto& pair : indexMap ) {
2348
2349 auto* indices = pair.first;
2350 auto& tracks = pair.second;
2351
2352 auto newEnd = std::remove_if( indices->begin(), indices->end(),
2353 [&]( auto& index ) {
2354 bool isConsistent = (this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( tracks.at(index), vertex );
2355 return !isConsistent;
2356 } );
2357
2358 indices->erase( newEnd, indices->end() );
2359
2360 }
2361
2362 }
2363
2364
2365 //____________________________________________________________________________________________________
2367
2368 const xAOD::EventInfo* eventInfo { nullptr };
2369 const xAOD::TruthParticleContainer* truthParticles { nullptr };
2370 const xAOD::TruthVertexContainer* truthVertices { nullptr };
2371
2372 auto sc0 = evtStore()->retrieve( eventInfo, "EventInfo" );
2373 if( sc0.isFailure() ) { return; }
2374
2375 if( !eventInfo->eventType( xAOD::EventInfo::IS_SIMULATION ) ) {
2376 return;
2377 }
2378
2379 auto sc1 = evtStore()->retrieve( truthParticles, "TruthParticles" );
2380 if( sc1.isFailure() ) { return; }
2381
2382 auto sc2 = evtStore()->retrieve( truthVertices, "TruthVertices" );
2383 if( sc2.isFailure() ) { return; }
2384
2385 if( !truthParticles ) { return; }
2386 if( !truthVertices ) { return; }
2387
2388 m_tracingTruthVertices.clear();
2389
2390 std::vector<const xAOD::TruthParticle*> truthSvTracks;
2391
2392 // truth particle selection functions
2393
2394 auto selectNone = [](const xAOD::TruthVertex*) ->bool { return false; };
2395
2396 auto selectRhadron = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2397 if( truthVertex->nIncomingParticles() != 1 ) return false;
2398 if( !truthVertex->incomingParticle(0) ) return false;
2399 if( abs(truthVertex->incomingParticle(0)->pdgId()) < 1000000 ) return false;
2400 if( abs(truthVertex->incomingParticle(0)->pdgId()) > 1000000000 ) return false; // Nuclear codes, e.g. deuteron
2401 // neutralino in daughters
2402 bool hasNeutralino = false;
2403 for( unsigned ip = 0; ip < truthVertex->nOutgoingParticles(); ip++ ) {
2404 const auto* p = truthVertex->outgoingParticle(ip);
2405 if( abs( p->pdgId() ) == 1000022 ) {
2406 hasNeutralino = true;
2407 break;
2408 }
2409 }
2410 return hasNeutralino;
2411 };
2412
2413 auto selectHNL = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2414 if( truthVertex->nIncomingParticles() != 1 ) return false;
2415 if( !truthVertex->incomingParticle(0) ) return false;
2416 if( abs(truthVertex->incomingParticle(0)->pdgId()) != 50 ) return false;
2417 return true;
2418 };
2419
2420 auto selectHiggs = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2421 if( truthVertex->nIncomingParticles() != 1 ) return false;
2422 if( !truthVertex->incomingParticle(0) ) return false;
2423 if( abs(truthVertex->incomingParticle(0)->pdgId()) != 36 ) return false;
2424 return true;
2425 };
2426
2427 auto selectKshort = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2428 if( truthVertex->nIncomingParticles() != 1 ) return false;
2429 if( !truthVertex->incomingParticle(0) ) return false;
2430 if( abs(truthVertex->incomingParticle(0)->pdgId()) != 310 ) return false;
2431 return true;
2432 };
2433
2434 auto selectBhadron = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2435 if( truthVertex->nIncomingParticles() != 1 ) return false;
2436 if( !truthVertex->incomingParticle(0) ) return false;
2437 if( abs(truthVertex->incomingParticle(0)->pdgId()) <= 500 || abs(truthVertex->incomingParticle(0)->pdgId()) >= 600 ) return false;
2438 return true;
2439 };
2440
2441 auto selectHadInt = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2442 if( truthVertex->nIncomingParticles() != 1 ) return false;
2443 if( !truthVertex->incomingParticle(0) ) return false;
2444
2445 const auto* parent = truthVertex->incomingParticle(0);
2446 if( parent->isLepton() ) return false;
2447
2448 TLorentzVector p4sum_in;
2449 TLorentzVector p4sum_out;
2450 for( unsigned ip = 0; ip < truthVertex->nIncomingParticles(); ip++ ) {
2451 const auto* particle = truthVertex->incomingParticle(ip);
2452 TLorentzVector p4; p4.SetPtEtaPhiM( particle->pt(), particle->eta(), particle->phi(), particle->m() );
2453 p4sum_in += p4;
2454 }
2455 for( unsigned ip = 0; ip < truthVertex->nOutgoingParticles(); ip++ ) {
2456 const auto* particle = truthVertex->outgoingParticle(ip);
2457 TLorentzVector p4; p4.SetPtEtaPhiM( particle->pt(), particle->eta(), particle->phi(), particle->m() );
2458 p4sum_out += p4;
2459 }
2460 return p4sum_out.E() - p4sum_in.E() >= 100.;
2461 };
2462
2463
2464
2465 using ParticleSelectFunc = bool (*)(const xAOD::TruthVertex*);
2466 std::map<std::string, ParticleSelectFunc> selectFuncs { { "", selectNone },
2467 { "Kshort", selectKshort },
2468 { "Bhadron", selectBhadron },
2469 { "Rhadron", selectRhadron },
2470 { "HNL", selectHNL },
2471 { "Higgs", selectHiggs },
2472 { "HadInt", selectHadInt } };
2473
2474
2475 if( selectFuncs.find( m_jp.truthParticleFilter ) == selectFuncs.end() ) {
2476 ATH_MSG_WARNING( " > " << __FUNCTION__ << ": invalid function specification: " << m_jp.truthParticleFilter );
2477 return;
2478 }
2479
2480 auto selectFunc = selectFuncs.at( m_jp.truthParticleFilter );
2481
2482 // loop over truth vertices
2483 for( const auto *truthVertex : *truthVertices ) {
2484 if( selectFunc( truthVertex ) ) {
2485 m_tracingTruthVertices.emplace_back( truthVertex );
2486 std::string msg;
2487 msg += Form("pdgId = %d, (r, z) = (%.2f, %.2f), ", truthVertex->incomingParticle(0)->pdgId(), truthVertex->perp(), truthVertex->z());
2488 msg += Form("nOutgoing = %lu, ", truthVertex->nOutgoingParticles() );
2489 msg += Form("mass = %.3f GeV, pt = %.3f GeV", truthVertex->incomingParticle(0)->m()/1.e3, truthVertex->incomingParticle(0)->pt()/1.e3 );
2490 ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2491 }
2492 }
2493
2494 if( m_jp.FillHist ) {
2495 for( const auto* truthVertex : m_tracingTruthVertices ) {
2496 m_hists["nMatchedTruths"]->Fill( 0., truthVertex->perp() );
2497 }
2498 }
2499
2500 }
2501
2502} // end of namespace VKalVrtAthena
2503
2504
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define minValue(current, test)
#define AmgSymMatrix(dim)
static Double_t sc
#define max(a, b)
Definition cfImp.cxx:41
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
MsgStream & msg() const
size_t index() const
Return the index of this element within its container.
bool patternCheck(const uint32_t &pattern, const Amg::Vector3D &vertex)
std::map< std::string, PatternStrategyFunc > m_patternStrategyFuncs
StatusCode refitVertexWithSuggestion(WrkVrt &, const Amg::Vector3D &)
refit the vertex with suggestion
bool passedFakeReject(const Amg::Vector3D &FitVertex, const xAOD::TrackParticle *itrk, const xAOD::TrackParticle *jtrk)
Flag false if the consistituent tracks are not consistent with the vertex position.
bool checkTrackHitPatternToVertex(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
A classical method with hard-coded geometry.
static bool patternCheckRun1(const uint32_t &pattern, const Amg::Vector3D &vertex)
bool checkTrackHitPatternToVertexByExtrapolationAssist(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
New method with track extrapolation.
void printWrkSet(const std::vector< WrkVrt > *WrkVrtSet, const std::string &name)
print the contents of reconstructed vertices
static bool patternCheckRun2(const uint32_t &pattern, const Amg::Vector3D &vertex)
bool checkTrackHitPatternToVertexByExtrapolation(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
New method with track extrapolation.
ToolHandle< IInDetConditionsTool > m_pixelCondSummaryTool
Condition service.
void trackClassification(std::vector< WrkVrt > *, std::map< long int, std::vector< long int > > &)
ExtrapolatedPattern * extrapolatedPattern(const xAOD::TrackParticle *, enum Trk::PropDirection)
double findWorstChi2ofMaximallySharedTrack(std::vector< WrkVrt > *, std::map< long int, std::vector< long int > > &, long int &, long int &)
std::unique_ptr< NtupleVars > m_ntupleVars
static bool patternCheckRun1OuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
void removeInconsistentTracks(WrkVrt &)
Remove inconsistent tracks from vertices.
ToolHandle< IInDetConditionsTool > m_sctCondSummaryTool
ToolHandle< Trk::IVertexMapper > m_vertexMapper
bool checkTrackHitPatternToVertexOuterOnly(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
A classical method with hard-coded geometry.
ToolHandle< Trk::IExtrapolator > m_extrapolator
ToolHandle< Trk::ITrackToVertexIPEstimator > m_trackToVertexIPEstimatorTool
std::map< std::string, TH1 * > m_hists
static void fillTrackSummary(track_summary &summary, const xAOD::TrackParticle *trk)
retrieve the track hit information
double improveVertexChi2(WrkVrt &)
attempt to improve the vertex chi2 by removing the most-outlier track one by one until the vertex chi...
double distanceBetweenVertices(const WrkVrt &, const WrkVrt &) const
calculate the physical distance
static size_t nTrkCommon(std::vector< WrkVrt > *WrkVrtSet, const std::pair< unsigned, unsigned > &pairIndex)
returns the number of tracks commonly present in both vertices
static double findMinVerticesNextPair(std::vector< WrkVrt > *, std::pair< unsigned, unsigned > &)
returns the next pair of vertices that give next-to-minimum distance significance
std::tuple< const TVector3, Detector, Bec, Layer, Flag > ExtrapolatedPoint
double findMinVerticesPair(std::vector< WrkVrt > *, std::pair< unsigned, unsigned > &, const AlgForVerticesPair &)
returns the pair of vertices that give minimum in terms of some observable (e.g.
const AtlasDetectorID * m_atlasId
static void removeTrackFromVertex(std::vector< WrkVrt > *, std::vector< std::deque< long int > > *, const long int &, const long int &)
ToolHandle< Reco::ITrackToVertex > m_trackToVertexTool
get a handle on the Track to Vertex tool
std::vector< ExtrapolatedPoint > ExtrapolatedPattern
const xAOD::VertexContainer * m_primaryVertices
StatusCode mergeVertices(WrkVrt &destination, WrkVrt &source)
the 2nd vertex is merged into the 1st vertex.
std::vector< const xAOD::TrackParticle * > m_associatedTracks
ToolHandle< Trk::ITrkVKalVrtFitter > m_fitSvc
std::vector< const xAOD::TrackParticle * > m_selectedTracks
StatusCode disassembleVertex(std::vector< WrkVrt > *, const unsigned &vertexIndex)
double(VrtSecInclusive::*)(const WrkVrt &, const WrkVrt &) const AlgForVerticesPair
bool patternCheckOuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
struct VKalVrtAthena::VrtSecInclusive::track_summary_properties track_summary
ToolHandle< Trk::ITruthToTrack > m_truthToTrack
static bool patternCheckRun2OuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
std::vector< const xAOD::TruthVertex * > m_tracingTruthVertices
double significanceBetweenVertices(const WrkVrt &, const WrkVrt &) const
calculate the significance (Mahalanobis distance) between two reconstructed vertices
std::map< const xAOD::TruthVertex *, bool > m_matchMap
bool eventType(EventType type) const
Check for one particular bitmask value.
@ IS_SIMULATION
true: simulation, false: data
float z0() const
Returns the parameter.
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
float d0() const
Returns the parameter.
uint32_t hitPattern() const
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
double chi2(TH1 *h0, TH1 *h1)
std::string algorithm
Definition hcg.cxx:85
Eigen::Matrix< double, 3, 1 > Vector3D
PropDirection
PropDirection, enum for direction of the propagation.
@ oppositeMomentum
@ alongMomentum
@ pixelBarrel0
there are three or four pixel barrel layers (R1/R2)
@ pixelEndCap0
three pixel discs (on each side)
const xAOD::TrackParticle * getLeptonTrackParticle(const xAOD::Muon *muon, const unsigned &trackType)
bool isAssociatedToVertices(const xAOD::TrackParticle *trk, const xAOD::VertexContainer *vertices)
double vtxVtxDistance(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
void genSequence(const xAOD::Muon *, std::vector< unsigned > &trackTypes)
Definition index.py:1
STL namespace.
DataModel_detail::iterator< DVL > unique(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of unique for DataVector/List.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.
const xAOD::TrackParticle * getOriginalTrackParticleFromGSF(const xAOD::TrackParticle *trkPar)
Helper function for getting the "Original" Track Particle (i.e before GSF) via the GSF Track Particle...
@ PriVtx
Primary vertex.
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
EventInfo_v1 EventInfo
Definition of the latest event info version.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthVertexContainer_v1 TruthVertexContainer
Declare the latest version of the truth vertex container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Muon_v1 Muon
Reference the current persistent version:
@ numberOfTRTHits
number of TRT hits [unit8_t].
@ numberOfNextToInnermostPixelLayerHits
these are the hits in the 1st pixel barrel layer
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.
Electron_v1 Electron
Definition of the current "egamma version".
double Chi2
VKalVrt fit covariance.
std::vector< double > vertexCov
VKalVrt fit vertex 4-momentum.
std::vector< double > Chi2PerTrk
VKalVrt fit chi2 result.
std::deque< long int > selectedTrackIndices
flagged true for good vertex candidates
std::deque< long int > associatedTrackIndices
list if indices in TrackParticleContainer for selectedBaseTracks
long int Charge
list of VKalVrt fit chi2 for each track
Amg::Vector3D vertex
list if indices in TrackParticleContainer for associatedTracks
TLorentzVector vertexMom
VKalVrt fit vertex position.
double closestWrkVrtValue
stores the index of the closest WrkVrt in std::vector<WrkVrt>
std::vector< std::vector< double > > TrkAtVrt
total charge of the vertex
unsigned long closestWrkVrtIndex
list of track parameters wrt the reconstructed vertex