ATLAS Offline Software
Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // 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 
29 using namespace std;
30 
31 namespace 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" );
78  return AlgConsts::maxValue;
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  //____________________________________________________________________________________________________
265  double VrtSecInclusive::improveVertexChi2( WrkVrt& vertex )
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 
385  double minValue = AlgConsts::maxValue;
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 
423  double minValue = AlgConsts::maxValue;
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  //____________________________________________________________________________________________________
452  StatusCode VrtSecInclusive::mergeVertices( WrkVrt& v1, WrkVrt& v2 )
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  //____________________________________________________________________________________________________
593  StatusCode VrtSecInclusive::refitVertexWithSuggestion( WrkVrt& workVertex,
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 = "" );
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  // Select tracks with additonal LRT Cuts (inspiried by Run 3 LRT optimization studies)
798  declareProperty("doSelectTracksWithLRTCuts", m_jp.doSelectTracksWithLRTCuts = false );
799 
800  // Additional dressing option
801  declareProperty("doAugmentDVimpactParametersToMuons", m_jp.doAugmentDVimpactParametersToMuons = false );
802  declareProperty("doAugmentDVimpactParametersToElectrons", m_jp.doAugmentDVimpactParametersToElectrons = false );
803 
804  // Additional ToolHandles
805  declareProperty("VertexFitterTool", m_fitSvc, " Private TrkVKalVrtFitter" );
806  declareProperty("Extrapolator", m_extrapolator );
807  declareProperty("TrackToVertexTool", m_trackToVertexTool );
808  declareProperty("TrackToVertexIPEstimatorTool", m_trackToVertexIPEstimatorTool );
809  declareProperty("VertexMapper", m_vertexMapper );
810  declareProperty("TruthToTrack", m_truthToTrack );
811 
812  }
813 
814 
815 
816 
817  //____________________________________________________________________________________________________
818  StatusCode VrtSecInclusive::processPrimaryVertices() {
819 
820  //--------------------------------------------------------
821  // Primary vertex extraction from TES
822  //
823 
824  ATH_CHECK( evtStore()->retrieve( m_primaryVertices, "PrimaryVertices") );
825 
826  if( m_jp.FillNtuple ) m_ntupleVars->get<unsigned int>( "NumPV" ) = 0;
827  m_thePV = nullptr;
828 
829  ATH_MSG_DEBUG( "processPrimaryVertices(): pv_size = " << m_primaryVertices->size() );
830 
831  // Loop over PV container and get number of tracks of each PV
832 
833  for( const auto *vertex : *m_primaryVertices ) {
834 
835  // Hide (2015-04-21): xAOD::Vertex may contain several types of vertices
836  // e.g. if VertexType==NoVtx, this is a dummy vertex.
837  // We only need to extract primary vertices here, and skip otherwise.
838 
839  if( xAOD::VxType::PriVtx != vertex->vertexType() ) continue;
840 
841  // Not considering pile-up; pick-up the first PV
842  m_thePV = vertex;
843 
844  if( m_jp.FillNtuple ) {
845 
846  if( 0 == m_ntupleVars->get<unsigned int>( "NumPV" ) ) {
847 
848  m_ntupleVars->get<double>( "PVX" ) = vertex->x();
849  m_ntupleVars->get<double>( "PVY" ) = vertex->y();
850  m_ntupleVars->get<double>( "PVZ" ) = vertex->z();
851  m_ntupleVars->get<unsigned int>( "PVType" ) = vertex->vertexType();
852 
853  // number of tracks associated to the PV
854  m_ntupleVars->get<unsigned int>( "NTrksPV" ) = vertex->nTrackParticles();
855  }
856 
857  m_ntupleVars->get<unsigned int>( "NumPV" )++;
858 
859  m_ntupleVars->get< vector<int> > ( "NdofTrksPV" ) .emplace_back( vertex->numberDoF() );
860  m_ntupleVars->get< vector<double> >( "PVZpile" ) .emplace_back( vertex->position().z() );
861  }
862 
863  ATH_MSG_DEBUG("PrimVertex x/y/z/nDOF "
864  << vertex->x() << ","
865  << vertex->y() << ","
866  << vertex->z() << ","
867  << vertex->numberDoF() );
868 
869  }
870 
871  // Use the dummy PV if no PV is composed
872  if( !m_thePV ) {
873  ATH_MSG_DEBUG("No Reconstructed PV was found. Using the dummy PV instead.");
874  for( const auto *vertex : *m_primaryVertices ) {
875  if( xAOD::VxType::NoVtx != vertex->vertexType() ) continue;
876 
877  if( m_jp.FillNtuple ) {
878  // Not considering pile-up; pick-up the first PV
879  if( 0 == m_ntupleVars->get<unsigned int>( "NumPV" ) ) {
880  m_thePV = vertex;
881 
882  m_ntupleVars->get<double>( "PVX" ) = vertex->x();
883  m_ntupleVars->get<double>( "PVY" ) = vertex->y();
884  m_ntupleVars->get<double>( "PVZ" ) = vertex->z();
885  m_ntupleVars->get<unsigned int>( "PVType" ) = vertex->vertexType();
886 
887  // number of tracks associated to the PV
888  m_ntupleVars->get<unsigned int>( "NTrksPV" ) = vertex->nTrackParticles();
889  }
890  }
891  }
892  }
893 
894  // if thePV is null, the PV is not found.
895  if( !m_thePV ) {
896  ATH_MSG_DEBUG("No PV is found in the PV collection!");
897  return StatusCode::FAILURE;
898  }
899 
900  ATH_MSG_DEBUG(" Primary vertex successful. thePV = " << m_thePV );
901 
902  return StatusCode::SUCCESS;
903  }
904 
905 
906  //____________________________________________________________________________________________________
907  void VrtSecInclusive::trackClassification(std::vector<WrkVrt> *workVerticesContainer, std::map<long int, std::vector<long int> >& trackToVertexMap )
908  {
909  // Fill TrkInVrt with vertex IDs of each track
910 
911  trackToVertexMap.clear();
912 
913  for( size_t iv = 0; iv<workVerticesContainer->size(); iv++ ) {
914 
915  WrkVrt& vertex = workVerticesContainer->at(iv);
916 
917  auto& trackIndices = vertex.selectedTrackIndices;
918  if( !vertex.isGood ) continue;
919  if( trackIndices.size() < 2 ) continue;
920 
921  for( auto& index : trackIndices ) {
922  trackToVertexMap[index].emplace_back( iv );
923  }
924  }
925 
926  for( auto& pair: trackToVertexMap ) {
927  std::string msg = Form("track index [%ld]: vertices = (", pair.first);
928  for( auto& v : pair.second ) {
929  msg += Form("%ld, ", v);
930  }
931  msg += ")";
932  if( pair.second.size() >=2 ) ATH_MSG_VERBOSE(" >> " << __FUNCTION__ << ": " << msg );
933  }
934 
935  }
936 
937 
938  //____________________________________________________________________________________________________
939  double VrtSecInclusive::findWorstChi2ofMaximallySharedTrack(std::vector<WrkVrt> *workVerticesContainer,
940  std::map<long int, std::vector<long int> >& trackToVertexMap,
941  long int & maxSharedTrack,
942  long int & worstMatchingVertex)
943  {
944 
945  double worstChi2 = AlgConsts::invalidFloat;
946 
947  // Find the track index that has the largest shared vertices
948  auto maxSharedTrackToVertices = std::max_element( trackToVertexMap.begin(), trackToVertexMap.end(), []( auto& p1, auto& p2 ) { return p1.second.size() < p2.second.size(); } );
949 
950  if( maxSharedTrackToVertices == trackToVertexMap.end() ) return worstChi2;
951 
952  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": max-shared track index = " << maxSharedTrackToVertices->first << ", number of shared vertices = " << maxSharedTrackToVertices->second.size() );
953 
954  if( maxSharedTrackToVertices->second.size() < 2 ) return worstChi2;
955 
956  // map of vertex index and the chi2 of the track for the maxSharedTrack
957  std::map<long int, double> vrtChi2Map;
958 
959  // loop over vertices for the largest shared track
960  for( auto& iv : maxSharedTrackToVertices->second ) {
961  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": loop over vertices: vertex index " << iv );
962 
963  auto& wrkvrt = workVerticesContainer->at( iv );
964  auto& trackIndices = wrkvrt.selectedTrackIndices;
965 
966  // find the index of the track
967  auto index = std::find_if( trackIndices.begin(), trackIndices.end(), [&]( auto& index ) { return index == maxSharedTrackToVertices->first; } );
968  if( index == trackIndices.end() ) {
969  ATH_MSG_WARNING(" >> " << __FUNCTION__ << ": index not found (algorithm inconsistent)" );
970  return worstChi2;
971  }
972 
973  auto& chi2 = wrkvrt.Chi2PerTrk.at( index - trackIndices.begin() );
974 
975  vrtChi2Map.emplace( std::pair<long int, double>(iv, chi2) );
976  }
977 
978  auto worstVrtChi2Pair = std::max_element( vrtChi2Map.begin(), vrtChi2Map.end(), []( auto& p1, auto& p2 ) { return p1.second < p2.second; } );
979 
980  if( worstVrtChi2Pair == vrtChi2Map.end() ) {
981  ATH_MSG_WARNING(" >> " << __FUNCTION__ << ": max_element of vrtChi2Map not found" );
982  return worstChi2;
983  }
984 
985  maxSharedTrack = maxSharedTrackToVertices->first;
986  worstMatchingVertex = worstVrtChi2Pair->first;
987  worstChi2 = worstVrtChi2Pair->second;
988 
989  return worstChi2;
990  }
991 
992 
993  //____________________________________________________________________________________________________
994  void VrtSecInclusive::printWrkSet(const std::vector<WrkVrt> *workVerticesContainer, const std::string& name)
995  {
996  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": ===============================================================" );
997  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": " << name << ": #vertices = " << workVerticesContainer->size() );
998 
999  std::set<const xAOD::TrackParticle*> usedTracks;
1000 
1001  auto concatenateIndicesToString = []( auto indices, auto& collection ) -> std::string {
1002  if( 0 == indices.size() ) return "";
1003  return std::accumulate( std::next(indices.begin()), indices.end(), std::to_string( indices.at(0) ),
1004  [&collection]( const std::string& str, auto& index ) { return str + ", " + std::to_string( collection.at(index)->index() ); } );
1005  };
1006 
1007  std::map<const xAOD::TruthVertex*, bool> previous;
1008 
1009  for( auto& pair : m_matchMap ) { previous.emplace( pair.first, pair.second ); }
1010 
1011  m_matchMap.clear();
1012  for( const auto* truthVertex : m_tracingTruthVertices ) { m_matchMap.emplace( truthVertex, false ); }
1013 
1014  for(size_t iv=0; iv<workVerticesContainer->size(); iv++) {
1015  const auto& wrkvrt = workVerticesContainer->at(iv);
1016 
1017  if( wrkvrt.nTracksTotal() < 2 ) continue;
1018 
1019  std::string sels = concatenateIndicesToString( wrkvrt.selectedTrackIndices, *m_selectedTracks );
1020  std::string assocs = concatenateIndicesToString( wrkvrt.associatedTrackIndices, *m_associatedTracks );
1021 
1022  for( const auto& index : wrkvrt.selectedTrackIndices ) { usedTracks.insert( m_selectedTracks->at(index) ); }
1023  for( const auto& index : wrkvrt.associatedTrackIndices ) { usedTracks.insert( m_associatedTracks->at(index) ); }
1024 
1025  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": " << name << " vertex [" << iv << "]: " << &wrkvrt
1026  << ", isGood = " << (wrkvrt.isGood? "true" : "false")
1027  << ", #ntrks(tot, sel, assoc) = (" << wrkvrt.nTracksTotal() << ", " << wrkvrt.selectedTrackIndices.size() << ", " << wrkvrt.associatedTrackIndices.size() << "), "
1028  << ", chi2/ndof = " << wrkvrt.fitQuality()
1029  << ", (r, z) = (" << wrkvrt.vertex.perp()
1030  << ", " << wrkvrt.vertex.z() << ")"
1031  << ", sels = { " << sels << " }"
1032  << ", assocs = { " << assocs << " }" );
1033 
1034  // Truth match condition
1035  for( const auto* truthVertex : m_tracingTruthVertices ) {
1036 
1037 
1038  Amg::Vector3D vTruth( truthVertex->x(), truthVertex->y(), truthVertex->z() );
1039  Amg::Vector3D vReco ( wrkvrt.vertex.x(), wrkvrt.vertex.y(), wrkvrt.vertex.z() );
1040 
1041  const auto distance = vReco - vTruth;
1042 
1043  AmgSymMatrix(3) cov;
1044  cov.fillSymmetric( 0, 0, wrkvrt.vertexCov.at(0) );
1045  cov.fillSymmetric( 1, 0, wrkvrt.vertexCov.at(1) );
1046  cov.fillSymmetric( 1, 1, wrkvrt.vertexCov.at(2) );
1047  cov.fillSymmetric( 2, 0, wrkvrt.vertexCov.at(3) );
1048  cov.fillSymmetric( 2, 1, wrkvrt.vertexCov.at(4) );
1049  cov.fillSymmetric( 2, 2, wrkvrt.vertexCov.at(5) );
1050 
1051  const double s2 = distance.transpose() * cov.inverse() * distance;
1052 
1053  if( distance.norm() < 2.0 || s2 < 100. ) m_matchMap.at( truthVertex ) = true;
1054 
1055  }
1056 
1057  }
1058 
1059  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": number of used tracks = " << usedTracks.size() );
1060 
1061  if( !previous.empty() && previous.size() == m_matchMap.size() ) {
1062  for( auto& pair : m_matchMap ) {
1063  if( previous.find( pair.first ) == previous.end() ) continue;
1064  if( pair.second != previous.at( pair.first ) ) {
1065  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": Match flag has changed: (r, z) = (" << pair.first->perp() << ", " << pair.first->z() << ")" );
1066  }
1067  }
1068  }
1069 
1070  if( m_jp.FillHist ) {
1071  for( auto& pair : m_matchMap ) {
1072  if( pair.second ) m_hists["nMatchedTruths"]->Fill( m_vertexingAlgorithmStep+2, pair.first->perp() );
1073  }
1074  }
1075 
1076  std::string msg;
1077  for( const auto* trk : usedTracks ) { msg += Form("%ld, ", trk->index() ); }
1078 
1079  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": used tracks = " << msg );
1080  ATH_MSG_DEBUG( " >> " << __FUNCTION__ << ": ===============================================================" );
1081 
1082  }
1083 
1084 
1085  //____________________________________________________________________________________________________
1086  void VrtSecInclusive::fillTrackSummary( track_summary& summary, const xAOD::TrackParticle *trk ) {
1087  summary.numIBLHits = 0;
1088  summary.numBLayerHits = 0;
1089  summary.numPixelHits = 0;
1090  summary.numSctHits = 0;
1091  summary.numTrtHits = 0;
1092 
1095  trk->summaryValue( summary.numPixelHits, xAOD::numberOfPixelHits );
1096  trk->summaryValue( summary.numSctHits, xAOD::numberOfSCTHits );
1097  trk->summaryValue( summary.numTrtHits, xAOD::numberOfTRTHits );
1098  }
1099 
1100 
1101 
1102  //____________________________________________________________________________________________________
1103  VrtSecInclusive::ExtrapolatedPattern* VrtSecInclusive::extrapolatedPattern( const xAOD::TrackParticle* trk, enum Trk::PropDirection direction ) {
1104 
1105  auto* pattern = new ExtrapolatedPattern;
1106  const EventContext& ctx = Gaudi::Hive::currentContext();
1107  std::vector<std::unique_ptr<Trk::TrackParameters>> paramsVector =
1108  m_extrapolator->extrapolateBlindly(ctx, trk->perigeeParameters(), direction);
1109 
1111 
1112  auto nDisabled = 0;
1113 
1114  for( auto& params : paramsVector ) {
1115 
1116  const TVector3 position( params->position().x(), params->position().y(), params->position().z() );
1117 
1118  if( prevPos == position ) {
1119  continue;
1120  }
1121 
1122  prevPos = position;
1123 
1124  const auto* detElement = params->associatedSurface().associatedDetectorElement();
1125 
1126  if( detElement ) {
1127 
1128  enum { Pixel = 1, SCT = 2 };
1129 
1130  const auto& id = detElement->identify();
1131  Flag good = false;
1132 
1133  if( m_atlasId->is_pixel(id) ) {
1134 
1135  auto idHash = m_pixelId->wafer_hash( id );
1136  good = m_pixelCondSummaryTool->isGood( idHash, ctx );
1137 
1138  pattern->emplace_back( std::make_tuple( position, Pixel, m_pixelId->barrel_ec(id), m_pixelId->layer_disk(id), good ) );
1139 
1140  } else if( m_atlasId->is_sct(id) ) {
1141 
1142  auto idHash = m_sctId->wafer_hash( id );
1143  good = m_sctCondSummaryTool->isGood( idHash, ctx );
1144 
1145  pattern->emplace_back( std::make_tuple( position, SCT, m_sctId->barrel_ec(id), m_sctId->layer_disk(id), good ) );
1146 
1147  }
1148 
1149  if( !pattern->empty() ) {
1150 
1151  ATH_MSG_VERBOSE(" >> " << __FUNCTION__ << ", track " << trk << ": position = (" << position.Perp() << ", " << position.z() << ", " << position.Phi() << "), detElement ID = " << id << ", good = " << good
1152  << ": (det, bec, layer) = (" << std::get<1>( pattern->back() ) << ", " << std::get<2>( pattern->back() ) << ", " << std::get<3>( pattern->back() ) << ")" );
1153 
1154  if( !good ) nDisabled++;
1155  }
1156 
1157  }
1158 
1159  }
1160 
1161  if( m_jp.FillHist ) {
1162  m_hists["disabledCount"]->Fill( nDisabled );
1163  }
1164 
1165  return pattern;
1166 
1167  }
1168 
1169 
1170  //____________________________________________________________________________________________________
1171  bool VrtSecInclusive::checkTrackHitPatternToVertexByExtrapolation( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex )
1172  {
1173 
1174  if( m_extrapolatedPatternBank.find( trk ) == m_extrapolatedPatternBank.end() ) {
1175 
1176  std::unique_ptr<ExtrapolatedPattern> exPattern_along( extrapolatedPattern( trk, Trk::alongMomentum ) );
1177  std::unique_ptr<ExtrapolatedPattern> exPattern_oppos( extrapolatedPattern( trk, Trk::oppositeMomentum ) );
1178 
1179  m_extrapolatedPatternBank.emplace( trk, std::make_pair( std::move(exPattern_along), std::move(exPattern_oppos) ) );
1180 
1181  }
1182 
1183  auto& exPattern = m_extrapolatedPatternBank.at( trk );
1184 
1185  using LayerCombination = std::vector<int>;
1186 
1187  std::map<LayerCombination, unsigned> layerMap;
1188  if( layerMap.empty() ) {
1189  layerMap[ { 1, 0, 0 } ] = Trk::pixelBarrel0;
1190  layerMap[ { 1, 0, 1 } ] = Trk::pixelBarrel1;
1191  layerMap[ { 1, 0, 2 } ] = Trk::pixelBarrel2;
1192  layerMap[ { 1, 0, 3 } ] = Trk::pixelBarrel3;
1193 
1194  layerMap[ { 1, 2, 0 } ] = Trk::pixelEndCap0;
1195  layerMap[ { 1, 2, 1 } ] = Trk::pixelEndCap1;
1196  layerMap[ { 1, 2, 2 } ] = Trk::pixelEndCap2;
1197  layerMap[ { 1,-2, 0 } ] = Trk::pixelEndCap0;
1198  layerMap[ { 1,-2, 1 } ] = Trk::pixelEndCap1;
1199  layerMap[ { 1,-2, 2 } ] = Trk::pixelEndCap2;
1200 
1201  layerMap[ { 2, 0, 0 } ] = Trk::sctBarrel0;
1202  layerMap[ { 2, 0, 1 } ] = Trk::sctBarrel1;
1203  layerMap[ { 2, 0, 2 } ] = Trk::sctBarrel2;
1204  layerMap[ { 2, 0, 3 } ] = Trk::sctBarrel3;
1205 
1206  layerMap[ { 2, 2, 0 } ] = Trk::sctEndCap0;
1207  layerMap[ { 2, 2, 1 } ] = Trk::sctEndCap1;
1208  layerMap[ { 2, 2, 2 } ] = Trk::sctEndCap2;
1209  layerMap[ { 2, 2, 3 } ] = Trk::sctEndCap3;
1210  layerMap[ { 2, 2, 4 } ] = Trk::sctEndCap4;
1211  layerMap[ { 2, 2, 5 } ] = Trk::sctEndCap5;
1212  layerMap[ { 2, 2, 6 } ] = Trk::sctEndCap6;
1213  layerMap[ { 2, 2, 7 } ] = Trk::sctEndCap7;
1214  layerMap[ { 2, 2, 8 } ] = Trk::sctEndCap8;
1215  layerMap[ { 2,-2, 0 } ] = Trk::sctEndCap0;
1216  layerMap[ { 2,-2, 1 } ] = Trk::sctEndCap1;
1217  layerMap[ { 2,-2, 2 } ] = Trk::sctEndCap2;
1218  layerMap[ { 2,-2, 3 } ] = Trk::sctEndCap3;
1219  layerMap[ { 2,-2, 4 } ] = Trk::sctEndCap4;
1220  layerMap[ { 2,-2, 5 } ] = Trk::sctEndCap5;
1221  layerMap[ { 2,-2, 6 } ] = Trk::sctEndCap6;
1222  layerMap[ { 2,-2, 7 } ] = Trk::sctEndCap7;
1223  layerMap[ { 2,-2, 8 } ] = Trk::sctEndCap8;
1224  }
1225 
1226  enum { position=0, detector=1, bec=2, layer=3, isGood=4 };
1227 
1228  // Lambda!
1229  auto getDetectorType = [&]( const ExtrapolatedPoint& point ) -> unsigned {
1230 
1231  const LayerCombination comb { std::get<detector>( point ), std::get<bec>( point ), std::get<layer>( point ) };
1232 
1233  for( auto& pair : layerMap ) {
1234  if( pair.first == comb ) {
1235  return pair.second;
1236  }
1237  }
1238 
1240  };
1241 
1242  enum { kShouldNotHaveHit, kShouldHaveHit, kMayHaveHit };
1243  std::vector<unsigned> expectedHitPattern(Trk::numberOfDetectorTypes, kShouldNotHaveHit);
1244 
1245  auto minExpectedRadius = AlgConsts::maxValue;
1246 
1247  // Loop over extrapolated points (along direction)
1248  auto& exPattern_along = *( exPattern.first );
1249 
1250  for( auto itr = exPattern_along.begin(); itr != exPattern_along.end(); ++itr ) {
1251  if( std::next( itr ) == exPattern_along.end() ) continue;
1252 
1253  const auto& point = *itr;
1254  const auto& nextPoint = *( std::next( itr ) );
1255 
1256  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": isGood = " << std::get<isGood>( point ) );
1257 
1258  const auto& thisPos = std::get<position>( point );
1259  const auto& nextPos = std::get<position>( nextPoint );
1260 
1261  auto sectionVector = nextPos - thisPos;
1262  auto vertexVector = TVector3( vertex.x(), vertex.y(), vertex.z() ) - thisPos;
1263 
1264 
1265  const auto& detectorType = getDetectorType( point );
1266 
1267  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": detType = " << detectorType );
1268 
1269  if( detectorType == AlgConsts::invalidUnsigned ) continue;
1270  if( detectorType >= Trk::numberOfDetectorTypes ) continue;
1271 
1272  // if the vertex is nearby (within 10 mm), the hit may be presnet ("X")
1273  if( vertexVector.Mag() < 10. ) {
1274  expectedHitPattern.at( detectorType ) = kMayHaveHit;
1275  continue;
1276  }
1277 
1278  // if the front-end module is not active, then the hit is not expected,
1279  // which means the hit may be present
1280  if( !static_cast<bool>(std::get<isGood>( point )) ) {
1281  expectedHitPattern.at( detectorType ) = kMayHaveHit;
1282  continue;
1283  }
1284 
1285  // if the inner product of the above two vectors is positive,
1286  // then point is inner than the vertex.
1287  // Else, the point is outer than the vertex and expect to have hits
1288  // when the track is originated from the vertex.
1289 
1290  if( sectionVector.Mag() == 0. ) continue;
1291 
1292  ATH_MSG_VERBOSE( " > " << __FUNCTION__
1293  << ": hitPos = (" << thisPos.Perp() << ", " << thisPos.z() << ", " << thisPos.Phi() << ")"
1294  << ", sectionVec = (" << sectionVector.Perp() << ", " << sectionVector.z() << ", " << sectionVector.Phi() << ")"
1295  << ", vertexVec = (" << vertexVector.Perp() << ", " << vertexVector.z() << ", " << vertexVector.Phi() << ")"
1296  << ", cos(s,v) = " << sectionVector * vertexVector / ( sectionVector.Mag() * vertexVector.Mag() + AlgConsts::infinitesimal ) );
1297 
1298  if( sectionVector * vertexVector > 0. ) continue;
1299 
1300  if( minExpectedRadius > thisPos.Perp() ) minExpectedRadius = thisPos.Perp();
1301 
1302  // now, the hit is expected to present.
1303 
1304  expectedHitPattern.at( detectorType ) = kShouldHaveHit;
1305  }
1306 
1307  // Loop over extrapolated points (opposite direction)
1308  auto& exPattern_oppos = *( exPattern.second );
1309  bool oppositeFlag = false;
1310 
1311  for( auto itr = exPattern_oppos.begin(); itr != exPattern_oppos.end(); ++itr ) {
1312  if( std::next( itr ) == exPattern_oppos.end() ) continue;
1313 
1314  const auto& point = *itr;
1315  const auto& nextPoint = *( std::next( itr ) );
1316 
1317  const auto& thisPos = std::get<position>( point );
1318  const auto& nextPos = std::get<position>( nextPoint );
1319 
1320  auto sectionVector = nextPos - thisPos;
1321  auto vertexVector = TVector3( vertex.x(), vertex.y(), vertex.z() ) - thisPos;
1322 
1323  const auto& detectorType = getDetectorType( point );
1324 
1325  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": detType = " << detectorType );
1326 
1327  ATH_MSG_DEBUG( " > " << __FUNCTION__
1328  << ": hitPos = (" << thisPos.Perp() << ", " << thisPos.z() << ", " << thisPos.Phi() << ")"
1329  << ", vertex = (" << vertex.perp() << ", " << vertex.z() << ", " << vertex.phi() << ")"
1330  << ", cos(s,v) = " << sectionVector * vertexVector / ( sectionVector.Mag() * vertexVector.Mag() + AlgConsts::infinitesimal ) );
1331 
1332  if( detectorType == AlgConsts::invalidUnsigned ) continue;
1333  if( detectorType >= Trk::numberOfDetectorTypes ) continue;
1334 
1335  if( sectionVector * vertexVector < 0. ) {
1336  oppositeFlag = true;
1337  }
1338  }
1339 
1340  // If the first expected point's radius is smaller than the vertex radius,
1341  // it's the case that the vertex was reconstructed in the opposite phi-direction
1342  // to the track outgoing direction. Such a case should be rejected.
1343  // bool oppositeFlag = ( minExpectedRadius < vertex.perp() );
1344 
1345  std::string msg = "Expected hit pattern: ";
1346  for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1347  msg += Form("%s", expectedHitPattern.at(i) < kMayHaveHit? Form("%u", expectedHitPattern.at(i)) : "X" );
1348  }
1349  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
1350 
1351  msg = "Recorded hit pattern: ";
1352  for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1353  msg += Form("%u", ( trk->hitPattern() >> i ) & 1 );
1354  }
1355  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
1356 
1357  std::vector<unsigned> matchedLayers;
1358 
1359  for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
1360  const unsigned recordedPattern = ( (trk->hitPattern()>>i) & 1 );
1361 
1362  if( expectedHitPattern.at(i) == kMayHaveHit ) {
1363  matchedLayers.emplace_back( i );
1364  } else if( expectedHitPattern.at(i) == kShouldHaveHit ) {
1365  if( expectedHitPattern.at(i) != recordedPattern ) {
1366  break;
1367  } else {
1368  matchedLayers.emplace_back( i );
1369  }
1370  } else {
1371  if( expectedHitPattern.at(i) != recordedPattern ) {
1372  break;
1373  }
1374  }
1375 
1376  }
1377 
1378  uint8_t PixelHits = 0;
1379  uint8_t SctHits = 0;
1380  uint8_t TRTHits = 0;
1381  if( !(trk->summaryValue( PixelHits, xAOD::numberOfPixelHits ) ) ) PixelHits =0;
1382  if( !(trk->summaryValue( SctHits, xAOD::numberOfSCTHits ) ) ) SctHits =0;
1383  if( !(trk->summaryValue( TRTHits, xAOD::numberOfTRTHits ) ) ) TRTHits =0;
1384 
1385  auto dphi = trk->phi() - vertex.phi();
1386  while( dphi > TMath::Pi() ) dphi -= TMath::TwoPi();
1387  while( dphi < -TMath::Pi() ) dphi += TMath::TwoPi();
1388 
1389  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": vtx phi = " << vertex.phi() << ", track phi = " << trk->phi() << ", dphi = " << dphi
1390  << ", oppositeFlag = " << oppositeFlag
1391  << ", nPixelHits = " << static_cast<int>(PixelHits)
1392  << ", nSCTHits = " << static_cast<int>(SctHits)
1393  << ", nTRTHits = " << static_cast<int>(TRTHits)
1394  << ", nMatchedLayers = " << matchedLayers.size() );
1395 
1396  if( PixelHits == 0 && vertex.perp() > 300. ) {
1397  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": vertex r > 300 mm, w/o no pixel hits" );
1398  }
1399 
1400 
1401  // Requires the first 2 layers with the hit matches.
1402  if( matchedLayers.size() < 2 ) return false;
1403 
1404  // In case of the first matched layer is not within pixel barrel, requires the first 4 layers with the hit match
1405  if( matchedLayers.at(0) >= Trk::pixelEndCap0 ) {
1406  if( matchedLayers.size() < 4 ) return false;
1407  }
1408 
1409  // Sometimes the vertex is reconstructed at the opposite phi direction.
1410  // In this case, the pattern match may pass.
1411  // This can be avoided by requiring that the
1412  if( oppositeFlag ) return false;
1413 
1414  // The following condition should apply for vertices outer than IBL.
1415  if( false /*matchedLayers.at(0) > Trk::pixelBarrel0*/ ) {
1416 
1417  // If the dphi (defined above) is opposite, reject.
1418  if( fabs( dphi ) > TMath::Pi()/2.0 ) return false;
1419 
1420  // If the track is not within the forward hemisphere to the vertex, reject.
1421  TVector3 trkP; trkP.SetPtEtaPhi( trk->pt(), trk->eta(), trk->phi() );
1422  TVector3 vtx; vtx.SetXYZ( vertex.x(), vertex.y(), vertex.z() );
1423  if( trkP.Dot( vtx ) < 0. ) return false;
1424 
1425  }
1426 
1427  return true;
1428  }
1429 
1430 
1431  //____________________________________________________________________________________________________
1432  bool VrtSecInclusive::patternCheckRun2( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1433 
1434  //
1435  // rough guesses for active layers:
1436  // BeamPipe: 23.5-24.3
1437  // IBL: 31.0-38.4
1438  // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1439  // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1440  //
1441 
1442  const double rad = vertex.perp();
1443  const double absz = fabs( vertex.z() );
1444 
1445  // vertex area classification
1446  enum vertexArea {
1447  insideBeamPipe,
1448 
1449  insidePixelBarrel0,
1450  aroundPixelBarrel0,
1451 
1452  outsidePixelBarrel0_and_insidePixelBarrel1,
1453  aroundPixelBarrel1,
1454 
1455  outsidePixelBarrel1_and_insidePixelBarrel2,
1456  aroundPixelBarrel2,
1457 
1458  outsidePixelBarrel2_and_insidePixelBarrel3,
1459  aroundPixelBarrel3,
1460 
1461  outsidePixelBarrel3_and_insideSctBarrel0,
1462  aroundSctBarrel0,
1463 
1464  outsideSctBarrel0_and_insideSctBarrel1,
1465  aroundSctBarrel1,
1466  };
1467 
1468  // Mutually exclusive vertex position pattern
1469  int vertex_pattern = 0;
1470  if( rad < 23.50 ) {
1471  vertex_pattern = insideBeamPipe;
1472 
1473  } else if( rad < 31.0 && absz < 331.5 ) {
1474  vertex_pattern = insidePixelBarrel0;
1475 
1476  } else if( rad < 38.4 && absz < 331.5 ) {
1477  vertex_pattern = aroundPixelBarrel0;
1478 
1479  } else if( rad < 47.7 && absz < 400.5 ) {
1480  vertex_pattern = outsidePixelBarrel0_and_insidePixelBarrel1;
1481 
1482  } else if( rad < 54.4 && absz < 400.5 ) {
1483  vertex_pattern = aroundPixelBarrel1;
1484 
1485  } else if( rad < 85.5 && absz < 400.5 ) {
1486  vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1487 
1488  } else if( rad < 92.2 && absz < 400.5 ) {
1489  vertex_pattern = aroundPixelBarrel2;
1490 
1491  } else if( rad < 119.3 && absz < 400.5 ) {
1492  vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1493 
1494  } else if( rad < 126.1 && absz < 400.5 ) {
1495  vertex_pattern = aroundPixelBarrel3;
1496 
1497  } else if( rad < 290 && absz < 749.0 ) {
1498  vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1499 
1500  } else if( rad < 315 && absz < 749.0 ) {
1501  vertex_pattern = aroundSctBarrel0;
1502 
1503  } else if( rad < 360 && absz < 749.0 ) {
1504  vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1505 
1506  } else if( rad < 390 && absz < 749.0 ) {
1507  vertex_pattern = aroundSctBarrel1;
1508 
1509  } else {
1510  }
1511 
1512  unsigned nPixelLayers { 0 };
1513  {
1514  if ( pattern & (1 << Trk::pixelBarrel0) ) nPixelLayers++;
1515  if ( pattern & (1 << Trk::pixelBarrel1) ) nPixelLayers++;
1516  if ( pattern & (1 << Trk::pixelBarrel2) ) nPixelLayers++;
1517  if ( pattern & (1 << Trk::pixelBarrel3) ) nPixelLayers++;
1518  if ( pattern & (1 << Trk::pixelEndCap0) ) nPixelLayers++;
1519  if ( pattern & (1 << Trk::pixelEndCap1) ) nPixelLayers++;
1520  if ( pattern & (1 << Trk::pixelEndCap2) ) nPixelLayers++;
1521  }
1522 
1524  if( vertex_pattern == insideBeamPipe ) {
1525 
1526  if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1527  if( nPixelLayers < 3 ) return false;
1528 
1529 
1530  } else if( vertex_pattern == insidePixelBarrel0 ) {
1531 
1532  if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1533  if( nPixelLayers < 3 ) return false;
1534  }
1535 
1536 
1537  else if( vertex_pattern == aroundPixelBarrel0 ) {
1538 
1539  // require nothing for PixelBarrel0
1540  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1541  if( nPixelLayers < 2 ) return false;
1542  }
1543 
1544 
1545  else if( vertex_pattern == outsidePixelBarrel0_and_insidePixelBarrel1 ) {
1546 
1547  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1548  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1549  if( nPixelLayers < 2 ) return false;
1550  }
1551 
1552 
1553  else if( vertex_pattern == aroundPixelBarrel1 ) {
1554 
1555  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1556  // require nothing for PixelBarrel
1557  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1558  if( nPixelLayers < 2 ) return false;
1559  }
1560 
1561 
1562  else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1563 
1564  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1565  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1566  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1567  if( nPixelLayers < 2 ) return false;
1568  }
1569 
1570 
1571  else if( vertex_pattern == aroundPixelBarrel2 ) {
1572 
1573  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1574  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1575  // require nothing for PixelBarrel2
1576  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1577  }
1578 
1579 
1580  else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1581 
1582  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1583  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1584  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1585  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1586  }
1587 
1588  else if( vertex_pattern == aroundPixelBarrel3 ) {
1589 
1590  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1591  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1592  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1593  // require nothing for PixelBarrel3
1594  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1595  }
1596 
1597 
1598  else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1599 
1600  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1601  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1602  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1603  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1604  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1605  }
1606 
1607 
1608  else if( vertex_pattern == aroundSctBarrel0 ) {
1609 
1610  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1611  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1612  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1613  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1614  // require nothing for SctBarrel0
1615  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1616  }
1617 
1618 
1619  else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1620 
1621  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1622  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1623  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1624  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1625  if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1626  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1627  }
1628 
1629 
1630  else if( vertex_pattern == aroundSctBarrel1 ) {
1631  if( (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1632  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1633  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1634  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1635  if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1636  // require nothing for SctBarrel1
1637  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1638  }
1640 
1641  return true;
1642 
1643  }
1644 
1645  //____________________________________________________________________________________________________
1646  bool VrtSecInclusive::patternCheckRun2OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1647 
1648  //
1649  // rough guesses for active layers:
1650  // BeamPipe: 23.5-24.3
1651  // IBL: 31.0-38.4
1652  // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1653  // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1654  //
1655 
1656  const double rad = vertex.perp();
1657  const double absz = fabs( vertex.z() );
1658 
1659  // vertex area classification
1660  enum vertexArea {
1661  insideBeamPipe,
1662 
1663  insidePixelBarrel0,
1664  aroundPixelBarrel0,
1665 
1666  outsidePixelBarrel0_and_insidePixelBarrel1,
1667  aroundPixelBarrel1,
1668 
1669  outsidePixelBarrel1_and_insidePixelBarrel2,
1670  aroundPixelBarrel2,
1671 
1672  outsidePixelBarrel2_and_insidePixelBarrel3,
1673  aroundPixelBarrel3,
1674 
1675  outsidePixelBarrel3_and_insideSctBarrel0,
1676  aroundSctBarrel0,
1677 
1678  outsideSctBarrel0_and_insideSctBarrel1,
1679  aroundSctBarrel1,
1680  };
1681 
1682  // Mutually exclusive vertex position pattern
1683  int vertex_pattern = 0;
1684  if( rad < 23.50 ) {
1685  vertex_pattern = insideBeamPipe;
1686 
1687  } else if( rad < 31.0 && absz < 331.5 ) {
1688  vertex_pattern = insidePixelBarrel0;
1689 
1690  } else if( rad < 38.4 && absz < 331.5 ) {
1691  vertex_pattern = aroundPixelBarrel0;
1692 
1693  } else if( rad < 47.7 && absz < 400.5 ) {
1694  vertex_pattern = outsidePixelBarrel0_and_insidePixelBarrel1;
1695 
1696  } else if( rad < 54.4 && absz < 400.5 ) {
1697  vertex_pattern = aroundPixelBarrel1;
1698 
1699  } else if( rad < 85.5 && absz < 400.5 ) {
1700  vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1701 
1702  } else if( rad < 92.2 && absz < 400.5 ) {
1703  vertex_pattern = aroundPixelBarrel2;
1704 
1705  } else if( rad < 119.3 && absz < 400.5 ) {
1706  vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1707 
1708  } else if( rad < 126.1 && absz < 400.5 ) {
1709  vertex_pattern = aroundPixelBarrel3;
1710 
1711  } else if( rad < 290 && absz < 749.0 ) {
1712  vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1713 
1714  } else if( rad < 315 && absz < 749.0 ) {
1715  vertex_pattern = aroundSctBarrel0;
1716 
1717  } else if( rad < 360 && absz < 749.0 ) {
1718  vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1719 
1720  } else if( rad < 390 && absz < 749.0 ) {
1721  vertex_pattern = aroundSctBarrel1;
1722 
1723  } else {
1724  }
1725 
1726 
1727  unsigned nPixelLayers { 0 };
1728  {
1729  nPixelLayers += ( pattern & (1 << Trk::pixelBarrel0) );
1730  nPixelLayers += ( pattern & (1 << Trk::pixelBarrel1) );
1731  nPixelLayers += ( pattern & (1 << Trk::pixelBarrel2) );
1732  nPixelLayers += ( pattern & (1 << Trk::pixelBarrel3) );
1733  nPixelLayers += ( pattern & (1 << Trk::pixelEndCap0) );
1734  nPixelLayers += ( pattern & (1 << Trk::pixelEndCap1) );
1735  nPixelLayers += ( pattern & (1 << Trk::pixelEndCap2) );
1736  }
1737 
1739  if( vertex_pattern == insideBeamPipe ) {
1740 
1741  if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1742  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1743  if( nPixelLayers < 3 ) return false;
1744 
1745  } else if( vertex_pattern == insidePixelBarrel0 ) {
1746 
1747  if( ! (pattern & (1<<Trk::pixelBarrel0)) ) return false;
1748  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1749  if( nPixelLayers < 3 ) return false;
1750 
1751  }
1752 
1753 
1754  else if( vertex_pattern == aroundPixelBarrel0 ) {
1755 
1756  // require nothing for PixelBarrel0
1757  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1758  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1759  if( nPixelLayers < 3 ) return false;
1760  }
1761 
1762 
1763  else if( vertex_pattern == outsidePixelBarrel0_and_insidePixelBarrel1 ) {
1764 
1765  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1766  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1767  if( nPixelLayers < 3 ) return false;
1768  }
1769 
1770 
1771  else if( vertex_pattern == aroundPixelBarrel1 ) {
1772 
1773  // require nothing for PixelBarrel1
1774  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1775  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1776  if( nPixelLayers < 2 ) return false;
1777  }
1778 
1779 
1780  else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1781 
1782  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1783  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1784  if( nPixelLayers < 2 ) return false;
1785  }
1786 
1787 
1788  else if( vertex_pattern == aroundPixelBarrel2 ) {
1789 
1790  // require nothing for PixelBarrel2
1791  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1792  }
1793 
1794 
1795  else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1796 
1797  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1798  }
1799 
1800  else if( vertex_pattern == aroundPixelBarrel3 ) {
1801 
1802  // require nothing for PixelBarrel3
1803  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1804  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1805  }
1806 
1807 
1808  else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1809 
1810  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1811  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1812  }
1813 
1814 
1815  else if( vertex_pattern == aroundSctBarrel0 ) {
1816 
1817  // require nothing for SctBarrel0
1818  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1819  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1820  }
1821 
1822 
1823  else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1824 
1825  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1826  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1827  }
1828 
1829 
1830  else if( vertex_pattern == aroundSctBarrel1 ) {
1831  // require nothing for SctBarrel1
1832  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
1833  if( ! (pattern & (1<<Trk::sctBarrel3)) ) return false;
1834  }
1836 
1837  return true;
1838 
1839  }
1840 
1841  //____________________________________________________________________________________________________
1842  bool VrtSecInclusive::patternCheckRun1( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
1843  //
1844  // rough guesses for active layers:
1845  // BeamPipe: 25.0
1846  // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1847  // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1848  //
1849 
1850  const double rad = vertex.perp();
1851  const double absz = fabs( vertex.z() );
1852 
1853  // vertex area classification
1854  enum vertexArea {
1855  insideBeamPipe,
1856 
1857  insidePixelBarrel1,
1858  aroundPixelBarrel1,
1859 
1860  outsidePixelBarrel1_and_insidePixelBarrel2,
1861  aroundPixelBarrel2,
1862 
1863  outsidePixelBarrel2_and_insidePixelBarrel3,
1864  aroundPixelBarrel3,
1865 
1866  outsidePixelBarrel3_and_insideSctBarrel0,
1867  aroundSctBarrel0,
1868 
1869  outsideSctBarrel0_and_insideSctBarrel1,
1870  aroundSctBarrel1,
1871  };
1872 
1873  // Mutually exclusive vertex position pattern
1874  Int_t vertex_pattern = 0;
1875  if( rad < 25.00 ) {
1876  vertex_pattern = insideBeamPipe;
1877 
1878  } else if( rad < 47.7 && absz < 400.5 ) {
1879  vertex_pattern = insidePixelBarrel1;
1880 
1881  } else if( rad < 54.4 && absz < 400.5 ) {
1882  vertex_pattern = aroundPixelBarrel1;
1883 
1884  } else if( rad < 85.5 && absz < 400.5 ) {
1885  vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1886 
1887  } else if( rad < 92.2 && absz < 400.5 ) {
1888  vertex_pattern = aroundPixelBarrel2;
1889 
1890  } else if( rad < 119.3 && absz < 400.5 ) {
1891  vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1892 
1893  } else if( rad < 126.1 && absz < 400.5 ) {
1894  vertex_pattern = aroundPixelBarrel3;
1895 
1896  } else if( rad < 290 && absz < 749.0 ) {
1897  vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1898 
1899  } else if( rad < 315 && absz < 749.0 ) {
1900  vertex_pattern = aroundSctBarrel0;
1901 
1902  } else if( rad < 360 && absz < 749.0 ) {
1903  vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1904 
1905  } else if( rad < 390 && absz < 749.0 ) {
1906  vertex_pattern = aroundSctBarrel1;
1907 
1908  } else {
1909  }
1910 
1911 
1913  if( vertex_pattern == insideBeamPipe ) {
1914 
1915  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1916 
1917  }
1918 
1919 
1920  else if( vertex_pattern == insidePixelBarrel1 ) {
1921 
1922  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1923  }
1924 
1925 
1926  else if( vertex_pattern == aroundPixelBarrel1 ) {
1927 
1928  // require nothing for PixelBarrel1
1929  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1930  }
1931 
1932 
1933  else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1934 
1935  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1936  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1937  }
1938 
1939 
1940  else if( vertex_pattern == aroundPixelBarrel2 ) {
1941 
1942  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1943  // require nothing for PixelBarrel2
1944  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1945  }
1946 
1947 
1948  else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1949 
1950  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1951  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1952  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1953  }
1954 
1955  else if( vertex_pattern == aroundPixelBarrel3 ) {
1956 
1957  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1958  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1959  // require nothing for PixelBarrel3
1960  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1961  }
1962 
1963 
1964  else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1965 
1966  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1967  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1968  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1969  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
1970  }
1971 
1972 
1973  else if( vertex_pattern == aroundSctBarrel0 ) {
1974 
1975  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1976  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1977  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1978  // require nothing for SctBarrel0
1979  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1980  }
1981 
1982 
1983  else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1984 
1985  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1986  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1987  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1988  if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1989  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
1990  }
1991 
1992 
1993  else if( vertex_pattern == aroundSctBarrel1 ) {
1994  if( (pattern & (1<<Trk::pixelBarrel1)) ) return false;
1995  if( (pattern & (1<<Trk::pixelBarrel2)) ) return false;
1996  if( (pattern & (1<<Trk::pixelBarrel3)) ) return false;
1997  if( (pattern & (1<<Trk::sctBarrel0)) ) return false;
1998  // require nothing for SctBarrel1
1999  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
2000  }
2002 
2003  return true;
2004  }
2005 
2006  //____________________________________________________________________________________________________
2007  bool VrtSecInclusive::patternCheckRun1OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2008  //
2009  // rough guesses for active layers:
2010  // BeamPipe: 25.0
2011  // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
2012  // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
2013  //
2014 
2015  const double rad = vertex.perp();
2016  const double absz = fabs( vertex.z() );
2017 
2018  // vertex area classification
2019  enum vertexArea {
2020  insideBeamPipe,
2021 
2022  insidePixelBarrel1,
2023  aroundPixelBarrel1,
2024 
2025  outsidePixelBarrel1_and_insidePixelBarrel2,
2026  aroundPixelBarrel2,
2027 
2028  outsidePixelBarrel2_and_insidePixelBarrel3,
2029  aroundPixelBarrel3,
2030 
2031  outsidePixelBarrel3_and_insideSctBarrel0,
2032  aroundSctBarrel0,
2033 
2034  outsideSctBarrel0_and_insideSctBarrel1,
2035  aroundSctBarrel1,
2036  };
2037 
2038  // Mutually exclusive vertex position pattern
2039  Int_t vertex_pattern = 0;
2040  if( rad < 25.00 ) {
2041  vertex_pattern = insideBeamPipe;
2042 
2043  } else if( rad < 47.7 && absz < 400.5 ) {
2044  vertex_pattern = insidePixelBarrel1;
2045 
2046  } else if( rad < 54.4 && absz < 400.5 ) {
2047  vertex_pattern = aroundPixelBarrel1;
2048 
2049  } else if( rad < 85.5 && absz < 400.5 ) {
2050  vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
2051 
2052  } else if( rad < 92.2 && absz < 400.5 ) {
2053  vertex_pattern = aroundPixelBarrel2;
2054 
2055  } else if( rad < 119.3 && absz < 400.5 ) {
2056  vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
2057 
2058  } else if( rad < 126.1 && absz < 400.5 ) {
2059  vertex_pattern = aroundPixelBarrel3;
2060 
2061  } else if( rad < 290 && absz < 749.0 ) {
2062  vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
2063 
2064  } else if( rad < 315 && absz < 749.0 ) {
2065  vertex_pattern = aroundSctBarrel0;
2066 
2067  } else if( rad < 360 && absz < 749.0 ) {
2068  vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
2069 
2070  } else if( rad < 390 && absz < 749.0 ) {
2071  vertex_pattern = aroundSctBarrel1;
2072 
2073  } else {
2074  }
2075 
2076 
2078  if( vertex_pattern == insideBeamPipe ) {
2079 
2080  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
2081 
2082  }
2083 
2084 
2085  else if( vertex_pattern == insidePixelBarrel1 ) {
2086 
2087  if( ! (pattern & (1<<Trk::pixelBarrel1)) ) return false;
2088  }
2089 
2090 
2091  else if( vertex_pattern == aroundPixelBarrel1 ) {
2092 
2093  // require nothing for PixelBarrel1
2094  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
2095  }
2096 
2097 
2098  else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
2099 
2100  if( ! (pattern & (1<<Trk::pixelBarrel2)) ) return false;
2101  }
2102 
2103 
2104  else if( vertex_pattern == aroundPixelBarrel2 ) {
2105 
2106  // require nothing for PixelBarrel2
2107  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
2108  }
2109 
2110 
2111  else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
2112 
2113  if( ! (pattern & (1<<Trk::pixelBarrel3)) ) return false;
2114  }
2115 
2116  else if( vertex_pattern == aroundPixelBarrel3 ) {
2117 
2118  // require nothing for PixelBarrel3
2119  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
2120  }
2121 
2122 
2123  else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
2124 
2125  if( ! (pattern & (1<<Trk::sctBarrel0)) ) return false;
2126  }
2127 
2128 
2129  else if( vertex_pattern == aroundSctBarrel0 ) {
2130 
2131  // require nothing for SctBarrel0
2132  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
2133  }
2134 
2135 
2136  else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
2137 
2138  if( ! (pattern & (1<<Trk::sctBarrel1)) ) return false;
2139  }
2140 
2141 
2142  else if( vertex_pattern == aroundSctBarrel1 ) {
2143  // require nothing for SctBarrel1
2144  if( ! (pattern & (1<<Trk::sctBarrel2)) ) return false;
2145  }
2147 
2148  return true;
2149  }
2150 
2151  //____________________________________________________________________________________________________
2152  bool VrtSecInclusive::patternCheck( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2153  bool flag = false;
2154 
2155  if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run2 ) {
2156  flag = patternCheckRun2( pattern, vertex );
2157  } else if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run1 ) {
2158  flag = patternCheckRun1( pattern, vertex );
2159  }
2160 
2161  return flag;
2162  }
2163 
2164  //____________________________________________________________________________________________________
2165  bool VrtSecInclusive::patternCheckOuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex ) {
2166  bool flag = false;
2167 
2168  if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run2 ) {
2169  flag = patternCheckRun2OuterOnly( pattern, vertex );
2170  } else if( m_jp.geoModel == VKalVrtAthena::GeoModel::Run1 ) {
2171  flag = patternCheckRun1OuterOnly( pattern, vertex );
2172  }
2173 
2174  return flag;
2175  }
2176 
2177  //____________________________________________________________________________________________________
2178  bool VrtSecInclusive::checkTrackHitPatternToVertex( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex )
2179  {
2180 
2181  const uint32_t pattern = trk->hitPattern();
2182 
2183  return patternCheck( pattern, vertex );
2184 
2185  }
2186 
2187 
2188  //____________________________________________________________________________________________________
2189  bool VrtSecInclusive::checkTrackHitPatternToVertexOuterOnly( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex )
2190  {
2191 
2192  const uint32_t pattern = trk->hitPattern();
2193 
2194  return patternCheckOuterOnly( pattern, vertex );
2195 
2196  }
2197 
2198 
2199  //____________________________________________________________________________________________________
2200  bool VrtSecInclusive::checkTrackHitPatternToVertexByExtrapolationAssist( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex )
2201  {
2202 
2203  if( m_extrapolatedPatternBank.find( trk ) == m_extrapolatedPatternBank.end() ) {
2204 
2205  std::unique_ptr<ExtrapolatedPattern> exPattern_along( extrapolatedPattern( trk, Trk::alongMomentum ) );
2206 
2207  m_extrapolatedPatternBank.emplace( trk, std::make_pair( std::move(exPattern_along), nullptr ) );
2208 
2209  }
2210 
2211  if( vertex.perp() < 31.0 ) {
2212  double dphi = trk->phi() - vertex.phi();
2213  while( dphi > TMath::Pi() ) { dphi -= TMath::TwoPi(); }
2214  while( dphi < -TMath::Pi() ) { dphi += TMath::TwoPi(); }
2215  if( cos(dphi) < -0.8 ) return false;
2216  }
2217 
2218  auto& exPattern = m_extrapolatedPatternBank.at( trk );
2219 
2220  using LayerCombination = std::vector<int>;
2221 
2222  std::map<LayerCombination, unsigned> layerMap;
2223  if( layerMap.empty() ) {
2224  layerMap[ { 1, 0, 0 } ] = Trk::pixelBarrel0;
2225  layerMap[ { 1, 0, 1 } ] = Trk::pixelBarrel1;
2226  layerMap[ { 1, 0, 2 } ] = Trk::pixelBarrel2;
2227  layerMap[ { 1, 0, 3 } ] = Trk::pixelBarrel3;
2228 
2229  layerMap[ { 1, 2, 0 } ] = Trk::pixelEndCap0;
2230  layerMap[ { 1, 2, 1 } ] = Trk::pixelEndCap1;
2231  layerMap[ { 1, 2, 2 } ] = Trk::pixelEndCap2;
2232  layerMap[ { 1,-2, 0 } ] = Trk::pixelEndCap0;
2233  layerMap[ { 1,-2, 1 } ] = Trk::pixelEndCap1;
2234  layerMap[ { 1,-2, 2 } ] = Trk::pixelEndCap2;
2235 
2236  layerMap[ { 2, 0, 0 } ] = Trk::sctBarrel0;
2237  layerMap[ { 2, 0, 1 } ] = Trk::sctBarrel1;
2238  layerMap[ { 2, 0, 2 } ] = Trk::sctBarrel2;
2239  layerMap[ { 2, 0, 3 } ] = Trk::sctBarrel3;
2240 
2241  layerMap[ { 2, 2, 0 } ] = Trk::sctEndCap0;
2242  layerMap[ { 2, 2, 1 } ] = Trk::sctEndCap1;
2243  layerMap[ { 2, 2, 2 } ] = Trk::sctEndCap2;
2244  layerMap[ { 2, 2, 3 } ] = Trk::sctEndCap3;
2245  layerMap[ { 2, 2, 4 } ] = Trk::sctEndCap4;
2246  layerMap[ { 2, 2, 5 } ] = Trk::sctEndCap5;
2247  layerMap[ { 2, 2, 6 } ] = Trk::sctEndCap6;
2248  layerMap[ { 2, 2, 7 } ] = Trk::sctEndCap7;
2249  layerMap[ { 2, 2, 8 } ] = Trk::sctEndCap8;
2250  layerMap[ { 2,-2, 0 } ] = Trk::sctEndCap0;
2251  layerMap[ { 2,-2, 1 } ] = Trk::sctEndCap1;
2252  layerMap[ { 2,-2, 2 } ] = Trk::sctEndCap2;
2253  layerMap[ { 2,-2, 3 } ] = Trk::sctEndCap3;
2254  layerMap[ { 2,-2, 4 } ] = Trk::sctEndCap4;
2255  layerMap[ { 2,-2, 5 } ] = Trk::sctEndCap5;
2256  layerMap[ { 2,-2, 6 } ] = Trk::sctEndCap6;
2257  layerMap[ { 2,-2, 7 } ] = Trk::sctEndCap7;
2258  layerMap[ { 2,-2, 8 } ] = Trk::sctEndCap8;
2259  }
2260 
2261  enum { position=0, detector=1, bec=2, layer=3, isGood=4 };
2262 
2263  // Lambda!
2264  auto getDetectorType = [&]( const ExtrapolatedPoint& point ) -> unsigned {
2265 
2266  const LayerCombination comb { std::get<detector>( point ), std::get<bec>( point ), std::get<layer>( point ) };
2267 
2268  for( auto& pair : layerMap ) {
2269  if( pair.first == comb ) {
2270  return pair.second;
2271  }
2272  }
2273 
2275  };
2276 
2277  uint32_t disabledPattern { 0 };
2278 
2279  // Loop over extrapolated points (along direction)
2280  auto& exPattern_along = *( exPattern.first );
2281 
2282  for( auto itr = exPattern_along.begin(); itr != exPattern_along.end(); ++itr ) {
2283  if( std::next( itr ) == exPattern_along.end() ) continue;
2284 
2285  const auto& point = *itr;
2286 
2287  ATH_MSG_VERBOSE( " > " << __FUNCTION__ << ": isGood = " << std::get<isGood>( point ) );
2288 
2289  if( !std::get<isGood>( point ) ) {
2290  const auto& detectorType = getDetectorType( point );
2291  disabledPattern += (1 << detectorType);
2292  }
2293  }
2294 
2295  uint32_t hitPattern = trk->hitPattern();
2296  uint32_t modifiedPattern = disabledPattern | hitPattern;
2297 
2298  std::string msg = "Disabled hit pattern: ";
2299  for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
2300  msg += Form("%u", ( disabledPattern >> i ) & 1 );
2301  }
2302  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2303 
2304  msg = "Recorded hit pattern: ";
2305  for( unsigned i=0; i<Trk::numberOfDetectorTypes; i++) {
2306  msg += Form("%u", ( hitPattern >> i ) & 1 );
2307  }
2308  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2309 
2310  return patternCheck( modifiedPattern, vertex );
2311 
2312  }
2313 
2314 
2315  //____________________________________________________________________________________________________
2316  bool VrtSecInclusive::passedFakeReject( const Amg::Vector3D& FitVertex,
2317  const xAOD::TrackParticle *itrk,
2318  const xAOD::TrackParticle *jtrk )
2319  {
2320 
2321  const bool& check_itrk = ( this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( itrk, FitVertex );
2322  const bool& check_jtrk = ( this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( jtrk, FitVertex );
2323 
2324  return ( check_itrk && check_jtrk );
2325 
2326  }
2327 
2328 
2329  //____________________________________________________________________________________________________
2330  void VrtSecInclusive::removeInconsistentTracks( WrkVrt& wrkvrt )
2331  {
2332 
2333  const auto& vertex = wrkvrt.vertex;
2334 
2335  std::map< std::deque<long int>*, std::vector<const xAOD::TrackParticle*>& > indexMap
2336  {
2337  { &(wrkvrt.selectedTrackIndices), *m_selectedTracks }, { &(wrkvrt.associatedTrackIndices), *m_associatedTracks }
2338  };
2339 
2340  for( auto& pair : indexMap ) {
2341 
2342  auto* indices = pair.first;
2343  auto& tracks = pair.second;
2344 
2345  auto newEnd = std::remove_if( indices->begin(), indices->end(),
2346  [&]( auto& index ) {
2347  bool isConsistent = (this->*m_patternStrategyFuncs[m_checkPatternStrategy] )( tracks.at(index), vertex );
2348  return !isConsistent;
2349  } );
2350 
2351  indices->erase( newEnd, indices->end() );
2352 
2353  }
2354 
2355  }
2356 
2357 
2358  //____________________________________________________________________________________________________
2359  void VrtSecInclusive::dumpTruthInformation() {
2360 
2361  const xAOD::EventInfo* eventInfo { nullptr };
2362  const xAOD::TruthParticleContainer* truthParticles { nullptr };
2363  const xAOD::TruthVertexContainer* truthVertices { nullptr };
2364 
2365  auto sc0 = evtStore()->retrieve( eventInfo, "EventInfo" );
2366  if( sc0.isFailure() ) { return; }
2367 
2368  if( !eventInfo->eventType( xAOD::EventInfo::IS_SIMULATION ) ) {
2369  return;
2370  }
2371 
2372  auto sc1 = evtStore()->retrieve( truthParticles, "TruthParticles" );
2373  if( sc1.isFailure() ) { return; }
2374 
2375  auto sc2 = evtStore()->retrieve( truthVertices, "TruthVertices" );
2376  if( sc2.isFailure() ) { return; }
2377 
2378  if( !truthParticles ) { return; }
2379  if( !truthVertices ) { return; }
2380 
2381  m_tracingTruthVertices.clear();
2382 
2383  std::vector<const xAOD::TruthParticle*> truthSvTracks;
2384 
2385  // truth particle selection functions
2386 
2387  auto selectNone = [](const xAOD::TruthVertex*) ->bool { return false; };
2388 
2389  auto selectRhadron = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2390  if( truthVertex->nIncomingParticles() != 1 ) return false;
2391  if( !truthVertex->incomingParticle(0) ) return false;
2392  if( abs(truthVertex->incomingParticle(0)->pdgId()) < 1000000 ) return false;
2393  if( abs(truthVertex->incomingParticle(0)->pdgId()) > 1000000000 ) return false; // Nuclear codes, e.g. deuteron
2394  // neutralino in daughters
2395  bool hasNeutralino = false;
2396  for( unsigned ip = 0; ip < truthVertex->nOutgoingParticles(); ip++ ) {
2397  const auto* p = truthVertex->outgoingParticle(ip);
2398  if( abs( p->pdgId() ) == 1000022 ) {
2399  hasNeutralino = true;
2400  break;
2401  }
2402  }
2403  return hasNeutralino;
2404  };
2405 
2406  auto selectHNL = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2407  if( truthVertex->nIncomingParticles() != 1 ) return false;
2408  if( !truthVertex->incomingParticle(0) ) return false;
2409  if( abs(truthVertex->incomingParticle(0)->pdgId()) != 50 ) return false;
2410  return true;
2411  };
2412 
2413  auto selectHiggs = [](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()) != 36 ) return false;
2417  return true;
2418  };
2419 
2420  auto selectKshort = [](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()) != 310 ) return false;
2424  return true;
2425  };
2426 
2427  auto selectBhadron = [](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()) <= 500 || abs(truthVertex->incomingParticle(0)->pdgId()) >= 600 ) return false;
2431  return true;
2432  };
2433 
2434  auto selectHadInt = [](const xAOD::TruthVertex* truthVertex ) -> bool {
2435  if( truthVertex->nIncomingParticles() != 1 ) return false;
2436  if( !truthVertex->incomingParticle(0) ) return false;
2437 
2438  const auto* parent = truthVertex->incomingParticle(0);
2439  if( parent->isLepton() ) return false;
2440 
2441  TLorentzVector p4sum_in;
2442  TLorentzVector p4sum_out;
2443  for( unsigned ip = 0; ip < truthVertex->nIncomingParticles(); ip++ ) {
2444  const auto* particle = truthVertex->incomingParticle(ip);
2445  TLorentzVector p4; p4.SetPtEtaPhiM( particle->pt(), particle->eta(), particle->phi(), particle->m() );
2446  p4sum_in += p4;
2447  }
2448  for( unsigned ip = 0; ip < truthVertex->nOutgoingParticles(); ip++ ) {
2449  const auto* particle = truthVertex->outgoingParticle(ip);
2450  TLorentzVector p4; p4.SetPtEtaPhiM( particle->pt(), particle->eta(), particle->phi(), particle->m() );
2451  p4sum_out += p4;
2452  }
2453  return p4sum_out.E() - p4sum_in.E() >= 100.;
2454  };
2455 
2456 
2457 
2458  using ParticleSelectFunc = bool (*)(const xAOD::TruthVertex*);
2459  std::map<std::string, ParticleSelectFunc> selectFuncs { { "", selectNone },
2460  { "Kshort", selectKshort },
2461  { "Bhadron", selectBhadron },
2462  { "Rhadron", selectRhadron },
2463  { "HNL", selectHNL },
2464  { "Higgs", selectHiggs },
2465  { "HadInt", selectHadInt } };
2466 
2467 
2468  if( selectFuncs.find( m_jp.truthParticleFilter ) == selectFuncs.end() ) {
2469  ATH_MSG_WARNING( " > " << __FUNCTION__ << ": invalid function specification: " << m_jp.truthParticleFilter );
2470  return;
2471  }
2472 
2473  auto selectFunc = selectFuncs.at( m_jp.truthParticleFilter );
2474 
2475  // loop over truth vertices
2476  for( const auto *truthVertex : *truthVertices ) {
2477  if( selectFunc( truthVertex ) ) {
2478  m_tracingTruthVertices.emplace_back( truthVertex );
2479  std::string msg;
2480  msg += Form("pdgId = %d, (r, z) = (%.2f, %.2f), ", truthVertex->incomingParticle(0)->pdgId(), truthVertex->perp(), truthVertex->z());
2481  msg += Form("nOutgoing = %lu, ", truthVertex->nOutgoingParticles() );
2482  msg += Form("mass = %.3f GeV, pt = %.3f GeV", truthVertex->incomingParticle(0)->m()/1.e3, truthVertex->incomingParticle(0)->pt()/1.e3 );
2483  ATH_MSG_DEBUG( " > " << __FUNCTION__ << ": " << msg );
2484  }
2485  }
2486 
2487  if( m_jp.FillHist ) {
2488  for( const auto* truthVertex : m_tracingTruthVertices ) {
2489  m_hists["nMatchedTruths"]->Fill( 0., truthVertex->perp() );
2490  }
2491  }
2492 
2493  }
2494 
2495 } // end of namespace VKalVrtAthena
2496 
2497 
python.PyKernel.retrieve
def retrieve(aClass, aKey=None)
Definition: PyKernel.py:110
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
Trk::pixelEndCap2
@ pixelEndCap2
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:239
mergePhysValFiles.pattern
pattern
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:26
xAOD::TrackParticle_v1::pt
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition: TrackParticle_v1.cxx:73
VKalVrtAthena::VrtSecInclusive::WrkVrt::vertex
Amg::Vector3D vertex
list if indices in TrackParticleContainer for associatedTracks
Definition: VrtSecInclusive.h:343
VKalVrtAthena::AlgConsts::infinitesimal
constexpr double infinitesimal
Definition: Reconstruction/VKalVrt/VrtSecInclusive/VrtSecInclusive/Constants.h:19
xAOD::muon
@ muon
Definition: TrackingPrimitives.h:195
VKalVrtAthena::VrtSecInclusive::WrkVrt::Charge
long int Charge
list of VKalVrt fit chi2 for each track
Definition: VrtSecInclusive.h:349
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
algorithm
std::string algorithm
Definition: hcg.cxx:82
VKalVrtAthena::VrtSecInclusive::WrkVrt::associatedTrackIndices
std::deque< long int > associatedTrackIndices
list if indices in TrackParticleContainer for selectedBaseTracks
Definition: VrtSecInclusive.h:342
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
maxValue
#define maxValue(current, test)
Definition: CompoundLayerMaterialCreator.h:22
max
#define max(a, b)
Definition: cfImp.cxx:41
VKalVrtAthena::VrtSecInclusive::WrkVrt::isGood
bool isGood
Definition: VrtSecInclusive.h:340
xAOD::uint8_t
uint8_t
Definition: Muon_v1.cxx:575
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
ParticleGun_SamplingFraction.bec
int bec
Definition: ParticleGun_SamplingFraction.py:89
Pixel
Definition: PixelFEUtils.h:16
Trk::sctEndCap2
@ sctEndCap2
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:248
xAOD::uint32_t
setEventNumber uint32_t
Definition: EventInfo_v1.cxx:127
Trk::oppositeMomentum
@ oppositeMomentum
Definition: PropDirection.h:21
index
Definition: index.py:1
xAOD::TrackParticle_v1::eta
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition: TrackParticle_v1.cxx:77
Trk::sctEndCap0
@ sctEndCap0
and 9 sct discs (on each side)
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:246
TruthVertexContainer.h
TruthParticleContainer.h
Trk::indices
std::pair< long int, long int > indices
Definition: AlSymMatBase.h:24
VKalVrtAthena::AlgConsts::chi2PerTrackInitValue
constexpr double chi2PerTrackInitValue
Definition: Reconstruction/VKalVrt/VrtSecInclusive/VrtSecInclusive/Constants.h:21
xAOD::TrackParticle_v1::summaryValue
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
Definition: TrackParticle_v1.cxx:736
accumulate
bool accumulate(AccumulateMap &map, std::vector< module_t > const &modules, FPGATrackSimMatrixAccumulator const &acc)
Accumulates an accumulator (e.g.
Definition: FPGATrackSimMatrixAccumulator.cxx:22
VKalVrtAthena::GeoModel::Run2
@ Run2
Definition: VrtSecInclusive.h:72
plotBeamSpotVxVal.cov
cov
Definition: plotBeamSpotVxVal.py:201
xAOD::TrackParticle_v1::z0
float z0() const
Returns the parameter.
Trk::sctBarrel0
@ sctBarrel0
four sct barrel layers
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:241
ElectronxAODHelpers.h
VKalVrtAthena::VrtSecInclusive::WrkVrt::vertexMom
TLorentzVector vertexMom
VKalVrt fit vertex position.
Definition: VrtSecInclusive.h:344
ParticleJetTools::declareProperties
void declareProperties(T &tool, LabelNames *n)
Definition: ParticleJetLabelCommon.h:105
VKalVrtAthena::VrtSecInclusive::track_summary_properties
Definition: VrtSecInclusive.h:492
athena.value
value
Definition: athena.py:122
xAOD::numberOfPixelHits
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
Definition: TrackingPrimitives.h:259
xAOD::numberOfTRTHits
@ numberOfTRTHits
number of TRT hits [unit8_t].
Definition: TrackingPrimitives.h:275
Trk::alongMomentum
@ alongMomentum
Definition: PropDirection.h:20
drawFromPickle.cos
cos
Definition: drawFromPickle.py:36
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
NtupleVars.h
Trk::sctEndCap5
@ sctEndCap5
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:251
xAOD::EventInfo_v1::IS_SIMULATION
@ IS_SIMULATION
true: simulation, false: data
Definition: EventInfo_v1.h:151
xAOD::VxType::NoVtx
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
Definition: TrackingPrimitives.h:570
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:52
std::sort
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
Definition: DVL_algorithms.h:554
xAOD::Muon_v1
Class describing a Muon.
Definition: Muon_v1.h:38
xAOD::TrackParticle_v1::d0
float d0() const
Returns the parameter.
AthenaPoolTestRead.sc
sc
Definition: AthenaPoolTestRead.py:27
Trk::sctBarrel3
@ sctBarrel3
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:244
mergePhysValFiles.end
end
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:93
TRT::Hit::detector
@ detector
Definition: HitInfo.h:78
VrtSecInclusive.h
InDet::ExclusiveOrigin::Primary
@ Primary
Definition: InDetTrackTruthOriginDefs.h:163
VKalVrtAthena
Definition: AANT_Tools.cxx:24
xAOD::TrackParticle_v1::hitPattern
uint32_t hitPattern() const
VKalVrtAthena::VrtSecInclusive::ExtrapolatedPoint
std::tuple< const TVector3, Detector, Bec, Layer, Flag > ExtrapolatedPoint
Definition: VrtSecInclusive.h:365
Trk::pixelEndCap0
@ pixelEndCap0
three pixel discs (on each side)
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:237
VKalVrtAthena::isAssociatedToVertices
bool isAssociatedToVertices(const xAOD::TrackParticle *trk, const xAOD::VertexContainer *vertices)
Definition: Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx:34
xAOD::TrackParticle_v1::perigeeParameters
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
Definition: TrackParticle_v1.cxx:485
Trk::PropDirection
PropDirection
Definition: PropDirection.h:19
VKalVrtAthena::genSequence
void genSequence(const xAOD::Electron *electron, std::vector< unsigned > &trackTypes)
Definition: Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx:94
fillPileUpNoiseLumi.next
next
Definition: fillPileUpNoiseLumi.py:52
VKalVrtAthena::VrtSecInclusive::WrkVrt::Chi2PerTrk
std::vector< double > Chi2PerTrk
VKalVrt fit chi2 result.
Definition: VrtSecInclusive.h:348
Trk::pixelBarrel2
@ pixelBarrel2
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:234
lumiFormat.i
int i
Definition: lumiFormat.py:92
VKalVrtAthena::GeoModel::Run1
@ Run1
Definition: VrtSecInclusive.h:72
Trk::pixelEndCap1
@ pixelEndCap1
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:238
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
TRT::Hit::layer
@ layer
Definition: HitInfo.h:79
master.flag
bool flag
Definition: master.py:29
xAOD::VxType::PriVtx
@ PriVtx
Primary vertex.
Definition: TrackingPrimitives.h:571
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:522
VKalVrtAthena::VrtSecInclusive::WrkVrt::vertexCov
std::vector< double > vertexCov
VKalVrt fit vertex 4-momentum.
Definition: VrtSecInclusive.h:345
test_pyathena.parent
parent
Definition: test_pyathena.py:15
Trk::sctEndCap3
@ sctEndCap3
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:249
VKalVrtAthena::VrtSecInclusive::WrkVrt::Chi2
double Chi2
VKalVrt fit covariance.
Definition: VrtSecInclusive.h:346
Trk::sctEndCap8
@ sctEndCap8
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:254
find_tgc_unfilled_channelids.ip
ip
Definition: find_tgc_unfilled_channelids.py:3
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
Trk::sctEndCap7
@ sctEndCap7
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:253
TRT_PAI_gasdata::SC
const float SC[NC]
Cross sections for Carbon.
Definition: TRT_PAI_gasdata.h:255
Trk::pixelBarrel3
@ pixelBarrel3
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:235
std::remove_if
std::reverse_iterator< DataModel_detail::iterator< DVL > > remove_if(typename std::reverse_iterator< DataModel_detail::iterator< DVL > > beg, typename std::reverse_iterator< DataModel_detail::iterator< DVL > > end, Predicate pred)
Specialization of remove_if for DataVector/List.
Definition: DVL_algorithms.h:114
histSizes.list
def list(name, path='/')
Definition: histSizes.py:38
TrackSummary.h
VKalVrtAthena::VrtSecInclusive::AlgForVerticesPair
double(VrtSecInclusive::*)(const WrkVrt &, const WrkVrt &) const AlgForVerticesPair
Definition: VrtSecInclusive.h:479
VKalVrtAthena::vtxVtxDistance
double vtxVtxDistance(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
Definition: Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx:55
VKalVrtAthena::VrtSecInclusive::WrkVrt::closestWrkVrtIndex
unsigned long closestWrkVrtIndex
list of track parameters wrt the reconstructed vertex
Definition: VrtSecInclusive.h:351
SG::AuxElement::index
size_t index() const
Return the index of this element within its container.
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
DataVector::clear
void clear()
Erase all the elements in the collection.
SCT
Definition: SCT_ChipUtils.h:14
Trk::pixelBarrel1
@ pixelBarrel1
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:233
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:41
xAOD::numberOfNextToInnermostPixelLayerHits
@ numberOfNextToInnermostPixelLayerHits
these are the hits in the 1st pixel barrel layer
Definition: TrackingPrimitives.h:248
VKalVrtAthena::VrtSecInclusive::Flag
int Flag
Definition: VrtSecInclusive.h:364
Trk::sctEndCap4
@ sctEndCap4
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:250
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
ActsTrk::to_string
std::string to_string(const DetectorType &type)
Definition: GeometryDefs.h:34
Trk::sctEndCap1
@ sctEndCap1
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:247
VKalVrtAthena::VrtSecInclusive::ExtrapolatedPattern
std::vector< ExtrapolatedPoint > ExtrapolatedPattern
Definition: VrtSecInclusive.h:366
xAOD::EgammaHelpers::getOriginalTrackParticleFromGSF
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...
Definition: ElectronxAODHelpers.cxx:22
VKalVrtAthena::VrtSecInclusive::WrkVrt::selectedTrackIndices
std::deque< long int > selectedTrackIndices
flagged true for good vertex candidates
Definition: VrtSecInclusive.h:341
VKalVrtAthena::VrtSecInclusive::WrkVrt
Definition: VrtSecInclusive.h:339
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
xAOD::Electron_v1
Definition: Electron_v1.h:34
Trk::sctBarrel2
@ sctBarrel2
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:243
VKalVrtAthena::getLeptonTrackParticle
const xAOD::TrackParticle * getLeptonTrackParticle(const xAOD::Electron *electron, const unsigned &trackType)
Definition: Reconstruction/VKalVrt/VrtSecInclusive/src/Utilities.cxx:106
EventInfo.h
xAOD::EventInfo_v1
Class describing the basic event information.
Definition: EventInfo_v1.h:43
ReadCellNoiseFromCoolCompare.v2
v2
Definition: ReadCellNoiseFromCoolCompare.py:364
VKalVrtAthena::VrtSecInclusive::WrkVrt::closestWrkVrtValue
double closestWrkVrtValue
stores the index of the closest WrkVrt in std::vector<WrkVrt>
Definition: VrtSecInclusive.h:352
python.PyAthena.v
v
Definition: PyAthena.py:157
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
DeMoScan.index
string index
Definition: DeMoScan.py:362
std::unique
DataModel_detail::iterator< DVL > unique(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of unique for DataVector/List.
Definition: DVL_algorithms.h:135
ReadBchFromCool.good
good
Definition: ReadBchFromCool.py:433
VKalVrtAthena::AlgConsts::invalidFloat
constexpr double invalidFloat
Definition: Reconstruction/VKalVrt/VrtSecInclusive/VrtSecInclusive/Constants.h:23
Trk::IVKalState
Definition: IVKalState.h:21
VKalVrtAthena::VrtSecInclusive::WrkVrt::fitQuality
double fitQuality() const
Definition: VrtSecInclusive.h:357
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
VKalVrtAthena::VrtSecInclusive::WrkVrt::TrkAtVrt
std::vector< std::vector< double > > TrkAtVrt
total charge of the vertex
Definition: VrtSecInclusive.h:350
Trk::numberOfDetectorTypes
@ numberOfDetectorTypes
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:263
ReadCellNoiseFromCoolCompare.s2
s2
Definition: ReadCellNoiseFromCoolCompare.py:379
xAOD::EgammaParameters::electron
@ electron
Definition: EgammaEnums.h:18
declareProperty
#define declareProperty(n, p, h)
Definition: BaseFakeBkgTool.cxx:15
xAOD::numberOfSCTHits
@ numberOfSCTHits
number of hits in SCT [unit8_t].
Definition: TrackingPrimitives.h:268
str
Definition: BTagTrackIpAccessor.cxx:11
Trk::pixelBarrel0
@ pixelBarrel0
there are three or four pixel barrel layers (R1/R2)
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:232
Trk::sctBarrel1
@ sctBarrel1
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:242
xAOD::TrackParticle_v1
Class describing a TrackParticle.
Definition: TrackParticle_v1.h:43
PowhegControl_ttFCNC_NLO.params
params
Definition: PowhegControl_ttFCNC_NLO.py:226
IntersectionPos.h
xAOD::bool
setBGCode setTAP setLVL2ErrorBits bool
Definition: TrigDecision_v1.cxx:60
minValue
#define minValue(current, test)
Definition: CompoundLayerMaterialCreator.h:21
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
python.AutoConfigFlags.msg
msg
Definition: AutoConfigFlags.py:7
VKalVrtAthena::AlgConsts::invalidUnsigned
constexpr unsigned invalidUnsigned
Definition: Reconstruction/VKalVrt/VrtSecInclusive/VrtSecInclusive/Constants.h:25
python.SystemOfUnits.rad
int rad
Definition: SystemOfUnits.py:111
xAOD::TrackParticle_v1::phi
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
xAOD::numberOfInnermostPixelLayerHits
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
Definition: TrackingPrimitives.h:237
Trk::sctEndCap6
@ sctEndCap6
Definition: Tracking/TrkEvent/TrkTrackSummary/TrkTrackSummary/TrackSummary.h:252
Trk::previous
@ previous
Definition: BinningData.h:32
SCT_Monitoring::summary
@ summary
Definition: SCT_MonitoringNumbers.h:65