ATLAS Offline Software
Loading...
Searching...
No Matches
Sel2TrkVertices.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
7// Header include
10#include "MVAUtils/BDT.h"
11//-------------------------------------------------
12// Other stuff
15
16#include "TMath.h"
17#include "TH1.h"
18#include "TH2D.h"
19//
20
21
22namespace Rec{
23
24
25//
26//
27//--------------------------------------------------------
28// Template routine for 2track secondary vertices selection
29//
30
31 void NewVrtSecInclusiveTool::select2TrVrt(const EventContext& ctx,
32 std::vector<const xAOD::TrackParticle*> & selectedTracks,
33 const xAOD::Vertex & primVrt,
34 std::map<long int,std::vector<double>> & goodVrt,
35 compatibilityGraph_t& compatibilityGraph )
36 const
37 {
38 std::vector<const xAOD::NeutralParticle*> neutralPartDummy(0);
39 std::vector<const xAOD::TrackParticle*> tracksForFit(2,nullptr);
40 std::vector<double> impact,impactError;
41 std::vector<double> inpMass(2,m_massPi);
42 long int Charge;
43 int i,j;
44 StatusCode sc;
45 Vrt2Tr tmpVrt;
46 std::vector<Vrt2Tr> all2TrVrt(0);
47 TLorentzVector PSum2T;
48 Amg::Vector3D iniVrt(0.,0.,0.);
49//
50 int NTracks = (int) (selectedTracks.size());
51//
52// impact parameters with sign calculations
53//
54 uint8_t nPixelHits;
55 double signifR=0.,signifZ=0.;
56 std::vector<int> nPixHits(NTracks,0);
57 std::vector<double> trackSignif(NTracks),dRdZratio(NTracks);
58 for (i=0; i<NTracks; i++) {
59 m_fitSvc->VKalGetImpact(ctx, selectedTracks[i], primVrt.position(), 1, impact, impactError);
60 signifR = impact[0]/ sqrt(impactError[0]);
61 signifZ = impact[1]/ sqrt(impactError[2]);
62 trackSignif[i] = sqrt( signifR*signifR + signifZ*signifZ);
63 dRdZratio[i] = std::abs(signifR/signifZ);
64 if( !(selectedTracks[i]->summaryValue(nPixelHits,xAOD::numberOfPixelHits)) ) nPixelHits=0;
65 nPixHits[i]=nPixelHits;
66 if(m_fillHist){
67 Hists& h = getHists();
68 h.m_hb_impactR->Fill( signifR, m_w_1);
69 h.m_hb_impactZ->Fill( signifZ, m_w_1);
70 h.m_hb_impactRZ->Fill(signifR, signifZ, m_w_1);
71 h.m_hb_impact->Fill( trackSignif[i], m_w_1);
72 if( i<DevTuple::maxNTrk){
73 Hists& h = getHists();
74 h.m_curTup->pttrk[i]=selectedTracks[i]->pt();
75 h.m_curTup->d0trk[i]=selectedTracks[i]->d0();
76 h.m_curTup->Sig3D[i]=trackSignif[i];
77 h.m_curTup->idHF[i] =getIdHF(selectedTracks[i]);
78 h.m_curTup->dRdZrat[i] =dRdZratio[i];
79 h.m_curTup->displaced[i]=isDisplaced(selectedTracks[i]);
80
81 uint8_t TRTHits;
82 if( !(selectedTracks[i]->summaryValue( TRTHits,xAOD::numberOfTRTHits))) TRTHits=0;
83 h.m_curTup->trkTRT[i] =TRTHits;
84 h.m_curTup->etatrk[i] =selectedTracks[i]->eta();
85 }
86 }
87 }
88
89 if( m_fillHist ){
90 Hists& h = getHists();
91 h.m_curTup->nTrk=NTracks < DevTuple::maxNTrk ? NTracks : DevTuple::maxNTrk ;
92 h.m_curTup->n2Vrt=0;
93 }
94
95 std::vector<std::vector<std::tuple<int,float>>> trkCount(NTracks);
96 std::unique_ptr<Trk::IVKalState> state = m_fitSvc->makeState(ctx);
97 m_fitSvc->setMassInputParticles( inpMass, *state ); // Use pion masses for fit
98 for (i=0; i<NTracks-1; i++) {
99 if(trackSignif[i]<m_trkSigCut || dRdZratio[i]<m_dRdZRatioCut )continue;
100 for (j=i+1; j<NTracks; j++) {
101 if(trackSignif[j]<m_trkSigCut || dRdZratio[j]<m_dRdZRatioCut )continue;
102 PSum2T=selectedTracks[i]->p4()+selectedTracks[j]->p4();
103 if( std::abs(selectedTracks[i]->eta()-selectedTracks[j]->eta())==0 &&
104 std::abs(selectedTracks[i]->phi()-selectedTracks[j]->phi())==0 &&
105 std::abs(selectedTracks[i]->pt() -selectedTracks[j]->pt())==0 ) continue; //remove duplicated tracks
106 tracksForFit[0]=selectedTracks[i];
107 tracksForFit[1]=selectedTracks[j];
108 double minDZ=0.;
109 sc=m_fitSvc->VKalVrtFitFast(tracksForFit,tmpVrt.fitVertex,minDZ,*state); /* Fast crude estimation*/
110 if( sc.isFailure() ) { /* No initial estimation */
111 iniVrt=primVrt.position();
112 } else {
113 double cosMomVrtDir = projSV_PV(tmpVrt.fitVertex,primVrt,PSum2T);
114 if( cosMomVrtDir>0. ) iniVrt=tmpVrt.fitVertex; /* Good initial estimation */
115 else iniVrt=primVrt.position();
116 }
117 if(nPixHits[i]>0 && nPixHits[j]>0){ if(minDZ> m_fastZSVCut) continue; } // Drop SV candidates with big Z track-track distance.
118 else{ if(minDZ>2.*m_fastZSVCut) continue; } // Drop SV candidates with big Z track-track distance.
119 m_fitSvc->setApproximateVertex(iniVrt.x(), iniVrt.y(), iniVrt.z(), *state);
120 sc=m_fitSvc->VKalVrtFit(tracksForFit, neutralPartDummy, tmpVrt.fitVertex, tmpVrt.momentum, Charge,
121 tmpVrt.errorMatrix, tmpVrt.chi2PerTrk, tmpVrt.trkAtVrt, tmpVrt.chi2, *state, false );
122 if(sc.isFailure()) continue; /* No fit */
123 if(tmpVrt.chi2>20.) continue; /*Too bad Chi2 for nDoF=1 fit*/
126 float vQuality;
127 xAOD::Vertex testV;
128 testV.makePrivateStore();
129 testV.setPosition(tmpVrt.fitVertex);
130 std::vector<float> testVcov(tmpVrt.errorMatrix.begin(),tmpVrt.errorMatrix.end());
131 testV.setCovariance(testVcov);
132 testV.setFitQuality(tmpVrt.chi2,1.);
133 bool acceptV=m_ini_v2trselector->isgood(ctx,
134 std::make_pair(selectedTracks[i],selectedTracks[j]), testV,
135 std::make_pair(momAtVrt(tmpVrt.trkAtVrt[0]),momAtVrt(tmpVrt.trkAtVrt[1])), primVrt, vQuality);
136 if(!acceptV) continue; // Main 2-track vertex selection
137
138
139//Check close material layer. Not validated in detail, use with care.
140 double vrtR=tmpVrt.fitVertex.perp();
141 double dstMatSignif=1.e4;
142 if(m_removeTrkMatSignif>0. && vrtR>20.){
143 double vrtRErr=vrtRadiusError(tmpVrt.fitVertex,tmpVrt.errorMatrix );
144 if(vrtR<30.){ dstMatSignif=std::abs(vrtR-m_beampipeR)/vrtRErr;}
145 else { dstMatSignif=distToMatLayerSignificance(tmpVrt);} //Material in Pixel volume
146 if(dstMatSignif<m_removeTrkMatSignif)continue;
147 }
148//
149// Debugging and BDT
150 double minPtT = std::min(tracksForFit[0]->pt(),tracksForFit[1]->pt());
151 if( m_fillHist ){
152 float ihitR = selectedTracks[i]->radiusOfFirstHit();
153 float jhitR = selectedTracks[j]->radiusOfFirstHit();
154 int ihitIBL = getIBLHit(selectedTracks[i]);
155 int jhitIBL = getIBLHit(selectedTracks[j]);
156 int ihitBL = getBLHit (selectedTracks[i]);
157 int jhitBL = getBLHit (selectedTracks[j]);
158 double Prob2v=TMath::Prob(tmpVrt.chi2,1);
159 double cosSVPV=projSV_PV(tmpVrt.fitVertex, primVrt, tmpVrt.momentum);
160 ROOT::Math::PxPyPzMVector SVPV(tmpVrt.fitVertex.x()-primVrt.x(),
161 tmpVrt.fitVertex.y()-primVrt.y(),
162 tmpVrt.fitVertex.z()-primVrt.z(), 10.);
163
164 Hists& h = getHists();
165
166 if(Charge==0){h.m_hb_massPiPi->Fill(tmpVrt.momentum.M(),1.);}
167 h.m_hb_cosSVMom->Fill(cosSVPV,1.);
168 h.m_hb_etaSV->Fill(SVPV.Eta(),1.);
169
170 double Sig3D=0.,Sig2D=0., Dist2D=0.;
171 int idisk1=0,idisk2=0,idisk3=0,jdisk1=0,jdisk2=0,jdisk3=0;
172 int sumIBLHits = std::max(ihitIBL,0)+std::max(jhitIBL,0);
173 int sumBLHits = std::max(ihitBL, 0)+std::max(jhitBL, 0);
174 getPixelDiscs(selectedTracks[i],idisk1,idisk2,idisk3);
175 getPixelDiscs(selectedTracks[j],jdisk1,jdisk2,jdisk3);
176 vrtVrtDist(primVrt, tmpVrt.fitVertex, tmpVrt.errorMatrix, Sig3D);
177 Dist2D=vrtVrtDist2D(primVrt, tmpVrt.fitVertex, tmpVrt.errorMatrix, Sig2D);
178 int barVrt1=getProdVrtBarcode(tracksForFit[0],0.1); // FIXME barcode-based
179 int barVrt2=getProdVrtBarcode(tracksForFit[1],0.1); // FIXME barcode-based
180 h.m_hb_signif3D->Fill(Sig3D,1.);
181 h.m_curTup->VrtTrkHF [h.m_curTup->n2Vrt] = getIdHF(tracksForFit[0])+ getIdHF(tracksForFit[1]);
182 h.m_curTup->VrtTrkI [h.m_curTup->n2Vrt] = getG4Inter(tracksForFit[0])+ getG4Inter(tracksForFit[1]);
183 h.m_curTup->VrtTrueBar[h.m_curTup->n2Vrt] = (barVrt1 && barVrt2 && barVrt1==barVrt2) ? 1 : 0;
184 h.m_curTup->VrtTrueNear[h.m_curTup->n2Vrt] = checkTrue2TrVrt(tracksForFit[0],tracksForFit[1],0.1);
185 h.m_curTup->VrtCh [h.m_curTup->n2Vrt] = Charge;
186 h.m_curTup->VrtProb [h.m_curTup->n2Vrt] = Prob2v;
187 h.m_curTup->VrtSig3D [h.m_curTup->n2Vrt] = Sig3D;
188 h.m_curTup->VrtSig2D [h.m_curTup->n2Vrt] = Sig2D;
189 h.m_curTup->VrtDist2D[h.m_curTup->n2Vrt] = vrtR<20. ? Dist2D : vrtR;
190 h.m_curTup->VrtM [h.m_curTup->n2Vrt] = tmpVrt.momentum.M();
191 h.m_curTup->VrtZ [h.m_curTup->n2Vrt] = tmpVrt.fitVertex.z();
192 h.m_curTup->VrtPt [h.m_curTup->n2Vrt] = tmpVrt.momentum.Pt();
193 h.m_curTup->VrtEta [h.m_curTup->n2Vrt] = SVPV.Eta();
194 h.m_curTup->VrtIBL [h.m_curTup->n2Vrt] = sumIBLHits;
195 h.m_curTup->VrtBL [h.m_curTup->n2Vrt] = sumBLHits;
196 h.m_curTup->VrtCosSPM[h.m_curTup->n2Vrt] = cosSVPV;
197 h.m_curTup->VMinPtT [h.m_curTup->n2Vrt] = minPtT;
198 h.m_curTup->VMinS3DT [h.m_curTup->n2Vrt] = std::min(trackSignif[i],trackSignif[j]);
199 h.m_curTup->VMaxS3DT [h.m_curTup->n2Vrt] = std::max(trackSignif[i],trackSignif[j]);
200 h.m_curTup->VrtBDT [h.m_curTup->n2Vrt] = vQuality;
201 h.m_curTup->VrtHR1 [h.m_curTup->n2Vrt] = ihitR;
202 h.m_curTup->VrtHR2 [h.m_curTup->n2Vrt] = jhitR;
203 h.m_curTup->VrtDZ [h.m_curTup->n2Vrt] = minDZ;
204 h.m_curTup->VrtDisk [h.m_curTup->n2Vrt] = idisk1+10*idisk2+20*idisk3+30*jdisk1+40*jdisk2+50*jdisk3;
205 h.m_curTup->VSigMat [h.m_curTup->n2Vrt] = dstMatSignif;
206 h.m_curTup->VrtIT [h.m_curTup->n2Vrt] = i;
207 h.m_curTup->VrtJT [h.m_curTup->n2Vrt] = j;
208 if(h.m_curTup->n2Vrt<DevTuple::maxNVrt-1)h.m_curTup->n2Vrt++;
209 }
210//
211//--- Save good candidate for multi-vertex fit
212//
213 add_edge(i,j,compatibilityGraph);
214 goodVrt[NTracks*i+j]=std::vector<double>{tmpVrt.fitVertex.x(),tmpVrt.fitVertex.y(),tmpVrt.fitVertex.z()};
215 trkCount[i].emplace_back(j,vQuality); trkCount[j].emplace_back(i,vQuality);
216 }
217 }
218 //=== Resolve -!----!- case to speed up cluster finding
219 for(int t=0; t<NTracks; t++){
220 if(trkCount[t].size()==2){
221 i=std::get<0>(trkCount[t][0]);
222 j=std::get<0>(trkCount[t][1]);
223 if(trkCount[i].size()==1 && trkCount[j].size()==1 ){
224 if( std::get<1>(trkCount[t][0]) < std::get<1>(trkCount[t][1]) ) {
225 remove_edge(t,i,compatibilityGraph);
226 if(t<i)goodVrt.erase(NTracks*t+i); else goodVrt.erase(NTracks*i+t);
227 trkCount[i].clear();
228 trkCount[t].erase(trkCount[t].begin()+0);
229 } else {
230 remove_edge(t,j,compatibilityGraph);
231 if(t<j)goodVrt.erase(NTracks*t+j); else goodVrt.erase(NTracks*j+t);
232 trkCount[j].clear();
233 trkCount[t].erase(trkCount[t].begin()+1);
234 }
235 }
236 }
237 }
238
239 }
240
241
242} //end of namespace
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
static Double_t sc
size_t size() const
Number of registered mappings.
Header file for AthHistogramAlgorithm.
ToolHandle< Rec::ITwoTrackVertexSelector > m_ini_v2trselector
static double vrtVrtDist2D(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr, double &signif)
static int getProdVrtBarcode(const xAOD::TrackParticle *TP, float resolLimit=0.1)
static double projSV_PV(const Amg::Vector3D &SV, const xAOD::Vertex &PV, const TLorentzVector &Direction)
boost::adjacency_list< boost::listS, boost::vecS, boost::undirectedS > compatibilityGraph_t
void select2TrVrt(const EventContext &ctx, std::vector< const xAOD::TrackParticle * > &SelectedTracks, const xAOD::Vertex &primVrt, std::map< long int, std::vector< double > > &vrt, compatibilityGraph_t &compatibilityGraph) const
ToolHandle< Trk::TrkVKalVrtFitter > m_fitSvc
static void getPixelDiscs(const xAOD::TrackParticle *Part, int &d0Hit, int &d1Hit, int &d2Hit)
Gaudi::Property< float > m_trkSigCut
ROOT::Math::PxPyPzEVector momAtVrt(const std::vector< double > &inpTrk) const
Gaudi::Property< float > m_removeTrkMatSignif
Gaudi::Property< float > m_fastZSVCut
Gaudi::Property< float > m_dRdZRatioCut
static double vrtVrtDist(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr, double &signif)
Gaudi::Property< bool > m_fillHist
static double vrtRadiusError(const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr)
static bool checkTrue2TrVrt(const xAOD::TrackParticle *TP1, const xAOD::TrackParticle *TP2, float nearCut=0.1)
Gaudi::Property< float > m_beampipeR
float z() const
Returns the z position.
float y() const
Returns the y position.
void setCovariance(const std::vector< float > &value)
Sets the covariance matrix as a simple vector of values.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
void setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
const Amg::Vector3D & position() const
Returns the 3-pos.
float x() const
Returns the x position.
Eigen::Matrix< double, 3, 1 > Vector3D
Gaudi Tools.
Vertex_v1 Vertex
Define the latest version of the vertex class.
@ numberOfTRTHits
number of TRT hits [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
std::vector< std::vector< double > > trkAtVrt