ATLAS Offline Software
Reconstruction/VKalVrt/NewVrtSecInclusiveTool/src/CutTrk.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
7 
8 // Header include
11 
12 #include "TH1.h"
13 
14 //-------------------------------------------------
15 namespace Rec{
16 
17 
18 //==============================================================================================================
19 // xAOD based stuff
20 //
22  const xAOD::Vertex & primVrt)
23  const
24  {
25  std::vector<const xAOD::TrackParticle*>& inpTrk = xAODwrk->inpTrk;
26  std::vector<const xAOD::TrackParticle*>& selectedTracks = xAODwrk->listSelTracks;
27  std::vector<const xAOD::TrackParticle*>::const_iterator i_ntrk;
28  std::vector<double> impact,impactError;
29  std::multimap<double,const xAOD::TrackParticle*> orderedTrk;
30  if(m_fillHist){
31  Hists& h = getHists();
32  h.m_hb_ntrkInput->Fill( inpTrk.size()+0.1, m_w_1);
33  }
34  for (i_ntrk = inpTrk.begin(); i_ntrk < inpTrk.end(); ++i_ntrk) {
35 //
36 //-- Perigee in TrackParticle
37 //
38  if(m_fillHist){
39  Hists& h = getHists();
40  h.m_hb_trkSelect->Fill( 0., m_w_1);
41  }
42  if((*i_ntrk)->numberDoF() == 0) continue; //Protection
43  double trkChi2 = (*i_ntrk)->chiSquared() / (*i_ntrk)->numberDoF();
44  if(trkChi2 > m_cutChi2) continue;
45  if( (*i_ntrk)->pt() < m_cutPt) continue;
46  if(m_fillHist){
47  Hists& h = getHists();
48  h.m_hb_trkSelect->Fill( 1., m_w_1);
49  }
50 
51  const Trk::Perigee mPer=(*i_ntrk)->perigeeParameters() ;
52  const AmgSymMatrix(5) * locCov = mPer.covariance();
53  const double CovTrkMtx00 = (*locCov)(0,0);
54 
55  uint8_t PixelHits,SctHits,BLayHits,TRTHits;
56  if( !((*i_ntrk)->summaryValue(PixelHits,xAOD::numberOfPixelHits)) ) continue; // Track is
57  if( !((*i_ntrk)->summaryValue( SctHits,xAOD::numberOfSCTHits)) ) continue; // definitely
58  if( SctHits<3 ) continue; // bad
59  if( !((*i_ntrk)->summaryValue( TRTHits,xAOD::numberOfTRTHits)) ) continue;
60  if( !((*i_ntrk)->summaryValue(BLayHits,xAOD::numberOfInnermostPixelLayerHits))) BLayHits=0;
61  if(m_fillHist){
62  Hists& h = getHists();
63  h.m_hb_trkSelect->Fill( 2., m_w_1);
64  }
65 
66  Amg::Vector3D perigeePos=mPer.position();
67  double impactD0=sqrt( (perigeePos.x()-primVrt.x())*(perigeePos.x()-primVrt.x())
68  +(perigeePos.y()-primVrt.y())*(perigeePos.y()-primVrt.y()) );
69  double impactZ=perigeePos.z()-primVrt.z();
70  if(m_fillHist){
71  Hists& h = getHists();
72  h.m_hb_trkD0->Fill( impactD0, m_w_1);
73  h.m_hb_trkZ ->Fill( impactZ, m_w_1);
74  }
75  if(std::abs(impactZ)*std::sin((*i_ntrk)->theta())>m_cutZVrt) continue;
76  if(impactD0>m_cutD0Max) continue;
77  if(impactD0<m_cutD0Min) continue;
78  if(m_fillHist){
79  Hists& h = getHists();
80  h.m_hb_trkSelect->Fill( 3., m_w_1);
81  }
82 
83  double bX=xAODwrk->beamX + (perigeePos.z()-xAODwrk->beamZ)*xAODwrk->tanBeamTiltX;
84  double bY=xAODwrk->beamY + (perigeePos.z()-xAODwrk->beamZ)*xAODwrk->tanBeamTiltY;
85  double impactBeam=sqrt( (perigeePos.x()-bX)*(perigeePos.x()-bX) + (perigeePos.y()-bY)*(perigeePos.y()-bY));
86 //----Anti-pileup
87  double signifBeam = impactBeam / sqrt(CovTrkMtx00);
88  if(signifBeam < m_antiPileupSigRCut) continue; // cut against tracks from pileup vertices
89  if(m_fillHist){
90  Hists& h = getHists();
91  h.m_hb_trkSelect->Fill( 4., m_w_1);
92  }
93 //----
94  if(PixelHits < m_cutPixelHits) continue;
95  if(SctHits < m_cutSctHits) continue;
96  if((PixelHits+SctHits) < m_cutSiHits) continue;
97  if(BLayHits < m_cutBLayHits) continue;
98  if(std::abs((*i_ntrk)->eta())<1.9 && TRTHits < m_cutTRTHits) continue; //TRT hits must be present inside TRT
99  if(m_fillHist){
100  Hists& h = getHists();
101  h.m_hb_trkSelect->Fill( 5., m_w_1);
102  }
103 //----
104  orderedTrk.emplace(signifBeam,*i_ntrk);
105  selectedTracks.push_back(*i_ntrk);
106  }
107 //---- Order tracks according to ranks
108  std::map<double,const xAOD::TrackParticle*>::reverse_iterator rt=orderedTrk.rbegin();
109  selectedTracks.resize(orderedTrk.size());
110  for ( int cntt=0; rt!=orderedTrk.rend(); ++rt,++cntt) {selectedTracks[cntt]=(*rt).second;}
111  }
112 
113 }//end namespace
xAOD::Vertex_v1::x
float x() const
Returns the x position.
Rec::workVectorArrxAOD::beamY
double beamY
Definition: GNNVertexFitterTool.h:43
Rec::NewVrtSecInclusiveTool::selGoodTrkParticle
void selGoodTrkParticle(workVectorArrxAOD *xAODwrk, const xAOD::Vertex &primVrt) const
Definition: Reconstruction/VKalVrt/NewVrtSecInclusiveTool/src/CutTrk.cxx:21
Rec::NewVrtSecInclusiveTool::Hists
Definition: NewVrtSecInclusiveTool.h:101
xAOD::uint8_t
uint8_t
Definition: Muon_v1.cxx:557
Rec::NewVrtSecInclusiveTool::m_cutPixelHits
long int m_cutPixelHits
Definition: NewVrtSecInclusiveTool.h:141
Trk::ParametersT
Dummy class used to allow special convertors to be called for surfaces owned by a detector element.
Definition: EMErrorDetail.h:25
Rec::workVectorArrxAOD::listSelTracks
std::vector< const xAOD::TrackParticle * > listSelTracks
Definition: GNNVertexFitterTool.h:41
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
Rec::workVectorArrxAOD::tanBeamTiltX
double tanBeamTiltX
Definition: GNNVertexFitterTool.h:45
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:50
Rec::NewVrtSecInclusiveTool::m_antiPileupSigRCut
float m_antiPileupSigRCut
Definition: NewVrtSecInclusiveTool.h:159
TrkVKalVrtFitter.h
Rec::NewVrtSecInclusiveTool::m_cutPt
double m_cutPt
Definition: NewVrtSecInclusiveTool.h:146
Rec::NewVrtSecInclusiveTool::m_cutD0Min
double m_cutD0Min
Definition: NewVrtSecInclusiveTool.h:149
Rec::NewVrtSecInclusiveTool::m_cutD0Max
double m_cutD0Max
Definition: NewVrtSecInclusiveTool.h:148
Rec::NewVrtSecInclusiveTool::m_cutSiHits
long int m_cutSiHits
Definition: NewVrtSecInclusiveTool.h:143
Rec::NewVrtSecInclusiveTool::m_cutSctHits
long int m_cutSctHits
Definition: NewVrtSecInclusiveTool.h:140
Rec
Name: MuonSpContainer.h Package : offline/Reconstruction/MuonIdentification/muonEvent.
Definition: FakeTrackBuilder.h:10
Rec::workVectorArrxAOD::tanBeamTiltY
double tanBeamTiltY
Definition: GNNVertexFitterTool.h:46
xAOD::Vertex_v1::z
float z() const
Returns the z position.
Rec::NewVrtSecInclusiveTool::m_cutBLayHits
long int m_cutBLayHits
Definition: NewVrtSecInclusiveTool.h:144
Rec::workVectorArrxAOD
Definition: GNNVertexFitterTool.h:40
Rec::workVectorArrxAOD::inpTrk
std::vector< const xAOD::TrackParticle * > inpTrk
Definition: NewVrtSecInclusiveTool.h:67
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
NewVrtSecInclusiveTool.h
Rec::NewVrtSecInclusiveTool::m_w_1
double m_w_1
Definition: NewVrtSecInclusiveTool.h:99
h
xAOD::Vertex_v1
Class describing a Vertex.
Definition: Vertex_v1.h:42
Rec::NewVrtSecInclusiveTool::getHists
Hists & getHists() const
Definition: NewVrtSecInclusiveTool.cxx:361
xAOD::Vertex_v1::y
float y() const
Returns the y position.
xAOD::numberOfSCTHits
@ numberOfSCTHits
number of hits in SCT [unit8_t].
Definition: TrackingPrimitives.h:268
Rec::workVectorArrxAOD::beamX
double beamX
Definition: GNNVertexFitterTool.h:42
drawFromPickle.sin
sin
Definition: drawFromPickle.py:36
Rec::NewVrtSecInclusiveTool::m_fillHist
bool m_fillHist
Definition: NewVrtSecInclusiveTool.h:170
Rec::NewVrtSecInclusiveTool::m_cutChi2
double m_cutChi2
Definition: NewVrtSecInclusiveTool.h:150
Rec::NewVrtSecInclusiveTool::m_cutZVrt
double m_cutZVrt
Definition: NewVrtSecInclusiveTool.h:147
Rec::workVectorArrxAOD::beamZ
double beamZ
Definition: GNNVertexFitterTool.h:44
xAOD::numberOfInnermostPixelLayerHits
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
Definition: TrackingPrimitives.h:237
Rec::NewVrtSecInclusiveTool::m_cutTRTHits
long int m_cutTRTHits
Definition: NewVrtSecInclusiveTool.h:142