ATLAS Offline Software
Loading...
Searching...
No Matches
VrtSecMulti.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
7//
8// Header include
10//-------------------------------------------------
11// Other stuff
15#include "MVAUtils/BDT.h"
16
17#include "boost/graph/bron_kerbosch_all_cliques.hpp"
18#include "TMath.h"
19#include "TH1.h"
20
21#include <algorithm>
22#include <ranges>
23//
24
25
26namespace Rec{
27
28
29 std::vector<xAOD::Vertex*> NewVrtSecInclusiveTool::getVrtSecMulti( workVectorArrxAOD * xAODwrk,
30 const xAOD::Vertex & primVrt,
31 compatibilityGraph_t& compatibilityGraph )
32 const
33 {
34
35 const double probVrtMergeLimit=0.01;
36
37 int i,j;
38 int inpNPart=xAODwrk->inpTrk.size();
39 std::vector<const xAOD::NeutralParticle*> neutralPartDummy(0);
40 ATH_MSG_DEBUG( "getVrtSecMulti() called with NPart=" << inpNPart);
41
42 std::vector<xAOD::Vertex*> finalVertices(0);
43
44 if( inpNPart < 2 ) { return finalVertices;} // 0,1 track => nothing to do!
45
46 selGoodTrkParticle( xAODwrk, primVrt);
47 int nTracks = xAODwrk->listSelTracks.size();
48
49 if( nTracks < 2 ) { return finalVertices;} // 0,1 selected track => nothing to do!
50
51 ATH_MSG_DEBUG( "Number of selected tracks = " <<nTracks);
52
53 if(m_fillHist){
54 Hists& h = getHists();
55 h.m_hb_ntrksel->Fill( (double) nTracks, m_w_1);
56 }
57
58//
59// inpTrk[] - input track list
60// listSelTracks[] - list of good tracks in jet for vertex search
61//------------------------------------------------------------
62// Initial track list ready
63// Find 2track vertices
64//
65
66 std::map<long int,std::vector<double>> foundVrt2t;
67 select2TrVrt(xAODwrk->listSelTracks, primVrt, foundVrt2t, compatibilityGraph);
68
69//---
70 ATH_MSG_DEBUG(" Defined edges in the graph="<< num_edges(compatibilityGraph));
71 if(num_edges(compatibilityGraph)==0){ return finalVertices;} //No vertices!
72
73//
74// m_Incomp[] - main vector of pointers for multivertex search
75//-----------------------------------------------------------------------------------------------------
76// Secondary track list is ready
77// Creation of initial vertex set
78//
79
80
81 std::unique_ptr<std::vector<WrkVrt>> wrkVrtSet = std::make_unique<std::vector<WrkVrt>>();
82 WrkVrt newvrt; newvrt.Good=true;
83 std::unique_ptr<Trk::IVKalState> state = m_fitSvc->makeState();
84 StatusCode sc;
85 long int NPTR=0, nth=2; // VK nth=2 to speed up PGRAPH when it's used
86
87
88 std::vector< std::vector<int> > allCliques;
89 bron_kerbosch_all_cliques(compatibilityGraph, clique_visitor(allCliques));
90 for(int cq=0; cq<(int)allCliques.size();cq++){
91 newvrt.selTrk.clear();
92 NPTR=allCliques[cq].size();
93 for(i=0;i<NPTR;i++){ newvrt.selTrk.push_back(allCliques[cq][i]);}
94//==================================================
95 xAODwrk->tmpListTracks.clear();
96 TLorentzVector vpsum;
97 for(i=0;i<NPTR;i++) {
98 xAODwrk->tmpListTracks.push_back( xAODwrk->listSelTracks.at(newvrt.selTrk[i]) );
99 vpsum += xAODwrk->listSelTracks.at(newvrt.selTrk[i])->p4();
100 }
101 std::vector<double> iniVrtPos=estimVrtPos(nTracks,newvrt.selTrk,foundVrt2t);
102 m_fitSvc->setApproximateVertex(iniVrtPos[0], iniVrtPos[1], iniVrtPos[2], *state); /*Use as starting point*/
103 sc=m_fitSvc->VKalVrtFit(xAODwrk->tmpListTracks, neutralPartDummy,
104 newvrt.vertex, newvrt.vertexMom, newvrt.vertexCharge, newvrt.vertexCov,
105 newvrt.chi2PerTrk, newvrt.trkAtVrt, newvrt.chi2,
106 *state, false);
107 if( sc.isFailure() ) continue; /* Bad fit - goto next solution */
108 ATH_MSG_VERBOSE("Found IniVertex="<<newvrt.vertex[0]<<", "<<newvrt.vertex[1]<<", "<<newvrt.vertex[2]);
109 ATH_MSG_VERBOSE("with Chi2="<<newvrt.chi2<<" Ntrk="<<NPTR<<" trk1,2="<<newvrt.selTrk[0]<<", "<<newvrt.selTrk[1]);
110 if(NPTR==2 && newvrt.chi2>10.) continue; /* Bad 2track vertex */
111 if(newvrt.chi2PerTrk.size()==2) newvrt.chi2PerTrk[0]=newvrt.chi2PerTrk[1]=newvrt.chi2/2.;
112 newvrt.Good = true;
113 newvrt.projectedVrt=MomProjDist(newvrt.vertex, primVrt, newvrt.vertexMom); //3D SV-PV distance
114 wrkVrtSet->push_back(newvrt);
115 }
116 std::sort(wrkVrtSet->begin(),wrkVrtSet->end(),[](const WrkVrt& a, const WrkVrt& b){return a.selTrk.size()>b.selTrk.size();});
117//==================================================================================
118// boost::adjacency_list<boost::listS, boost::vecS, boost::undirectedS>::vertex_iterator vertexIt, vertexEnd;
119// boost::adjacency_list<boost::listS, boost::vecS, boost::undirectedS>::adjacency_iterator neighbourIt, neighbourEnd;
120// tie(vertexIt, vertexEnd) = vertices(compatibilityGraph);
121// for (; vertexIt != vertexEnd; ++vertexIt) { std::cout << *vertexIt << " is connected with ";
122// tie(neighbourIt, neighbourEnd) = adjacent_vertices(*vertexIt, compatibilityGraph);
123// for (; neighbourIt != neighbourEnd; ++neighbourIt) std::cout << *neighbourIt << " "; std::cout << "\n"; }
124//==================================================================================
125 if((*wrkVrtSet).empty())return finalVertices;
126 if(msgLvl(MSG::DEBUG))printWrkSet(wrkVrtSet.get(),"Initial Vertices");
127 //
128 //--Count track participation in different vertices
129 std::vector<int> trkNPairs(nTracks,0);
130 for(auto &vrt : (*wrkVrtSet)){
131 int ntInV=vrt.selTrk.size()-1;
132 for(auto &trk : vrt.selTrk)trkNPairs.at(trk) += ntInV;
133 }
134//
135//- Resolve all overlapped vertices
136//
137 state = m_fitSvc->makeState();
138 std::multimap<double,std::pair<int,int>> vrtWithCommonTrk;
139 while(true){
140 int nSoluI=(*wrkVrtSet).size();
141 vrtWithCommonTrk.clear();
142 unsigned int nTComMax=0;
143 for(int iv=0; iv<nSoluI-1; iv++ ){
144 if(!(*wrkVrtSet)[iv].Good) continue;
145 if( (*wrkVrtSet)[iv].selTrk.size()<nTComMax) continue; // Optimisation. Only biggest overlap matters
146 for(int jv=iv+1; jv<nSoluI; jv++){
147 if(!(*wrkVrtSet)[jv].Good) continue;
148 if( (*wrkVrtSet)[jv].selTrk.size()<nTComMax) continue; // Optimisation. Only biggest overlap matters
149 unsigned int nTCom=nTrkCommon( wrkVrtSet.get(), iv, jv);
150 if(!nTCom) continue;
151 if(nTCom<nTComMax) continue;
152 double sumChi2=(*wrkVrtSet)[iv].chi2+(*wrkVrtSet)[jv].chi2;
153 sumChi2=std::min(sumChi2,999.)*1.e-3;
154 vrtWithCommonTrk.emplace(nTCom+sumChi2,std::make_pair(iv,jv));
155 nTComMax=std::max(nTComMax,nTCom);
156 } }
157 if(vrtWithCommonTrk.empty())break;
158 //============================== DEBUG output
159 //for(auto ku : vrtWithCommonTrk)std::cout<<" nCom="<<ku.first<<" v1="<<ku.second.first<<" v2="<<ku.second.second<<'\n';
160 //if(msgLvl(MSG::DEBUG))printWrkSet(wrkVrtSet.get(),"Overlapped Vertex Cleaning");
161 //===========================================
162 for( const auto& ov : std::ranges::reverse_view(vrtWithCommonTrk)) {
163 WrkVrt & v1 = (*wrkVrtSet)[ov.second.first];
164 WrkVrt & v2 = (*wrkVrtSet)[ov.second.second];
165 if(!v1.Good)continue; //----One of the vertices is already preocessed
166 if(!v2.Good)continue;
167 //--Recheck amount of common tracks
168 unsigned int nTCom=nTrkCommon( wrkVrtSet.get(), ov.second.first, ov.second.second);
169 if(nTCom<nTComMax)continue; //----One of the vertices is already preocessed
170 //--First check if one vertex is fully contained in another
171 if( nTCom==v1.selTrk.size() || nTCom==v2.selTrk.size() ){
172 if(nTCom==v1.selTrk.size()){v1.Good = false; continue;}
173 if(nTCom==v2.selTrk.size()){v2.Good = false; continue;}
174 }
175 //--Then check if 2 vertices with common tracks can be simply merged
176 if( nTCom>1 && TMath::Prob( v1.chi2, 2*v1.selTrk.size()-3) > probVrtMergeLimit
177 && TMath::Prob( v2.chi2, 2*v2.selTrk.size()-3) > probVrtMergeLimit){
178 double prbV=mergeAndRefitVertices( v1, v2, newvrt, xAODwrk->listSelTracks, *state);
179 if(prbV>probVrtMergeLimit){
180 v1.Good = false; v2.Good = false;
181 newvrt.Good = true;
182 newvrt.projectedVrt=MomProjDist(newvrt.vertex, primVrt, newvrt.vertexMom); //3D SV-PV distance
183 std::swap(v1,newvrt); //Replace v1 by new vertex
184 continue;
185 } }
186 //--If not mergeable - refine them
187 refineVerticesWithCommonTracks( v1, v2, xAODwrk->listSelTracks, *state);
188 }
189 }
190 if(m_fillHist){
191 int cvgood=0; for(const auto& vrt:(*wrkVrtSet)) if(vrt.Good)cvgood++;
192 Hists& h = getHists();
193 h.m_hb_rawVrtN->Fill( (float)cvgood, m_w_1);
194 }
195//
196//-Clean duplicated 1track vertices if they exist
197//
199 for(auto &v1t : (*wrkVrtSet)){
200 if(v1t.selTrk.size()!=1 || !v1t.Good)continue;
201 int ind_t=v1t.selTrk[0];
202 if(trkNPairs[ind_t]<2){ v1t.Good=false; continue; } //Remove 1tr-vertex if track crosses only one other track
203 if( xAODwrk->listSelTracks[ind_t]->pt()<m_cutPt*2){ v1t.Good=false; continue; }; //Tighten track_pt cut for 1-track vertex
204 for(auto &vrt :(*wrkVrtSet)){ // Check if the track is present in another vertex, including other 1-track ones
205 if(!vrt.Good || &v1t==&vrt)continue;
206 if(std::find(vrt.selTrk.begin(),vrt.selTrk.end(),ind_t) != vrt.selTrk.end()){ v1t.Good=false; break; }
207 } }
208 }
209//
210//-Remove all bad vertices from the working set
211//
212 int tmpV=0; while( tmpV<(int)(*wrkVrtSet).size() )if( !(*wrkVrtSet)[tmpV].Good ) { (*wrkVrtSet).erase((*wrkVrtSet).begin()+tmpV);} else {tmpV++;}
213 if((*wrkVrtSet).empty())return finalVertices;
214 if(msgLvl(MSG::DEBUG))printWrkSet(wrkVrtSet.get(),"Intermediate Vertices");
215 for( auto &tmpV : (*wrkVrtSet) ) tmpV.projectedVrt=MomProjDist(tmpV.vertex, primVrt, tmpV.vertexMom ); //Setup ProjectedVrt
216//----------------------------------------------------------------------------
217
218//----------------------------------------------------------------------------
219//
220// Final check/merge for close vertices
221//
222 int foundV1=-1, foundV2=-1;
223 std::vector<double> checkedDst(0);
224 double minDistVV=minVrtVrtDist( wrkVrtSet.get(), foundV1, foundV2, checkedDst); //recalculate VV distances
225 if(m_fillHist){
226 Hists& h = getHists();
227 h.m_hb_distVV->Fill( minDistVV, m_w_1);
228 }
229 while ( minDistVV < m_vertexMergeCut) {
230 if(foundV1<foundV2) { int tmp=foundV1; foundV1=foundV2; foundV2=tmp;}
231 double probV=mergeAndRefitVertices( (*wrkVrtSet)[foundV1], (*wrkVrtSet)[foundV2], newvrt, xAODwrk->listSelTracks, *state, 0);
232 ATH_MSG_DEBUG( "Merged vertex prob=" << probV<<" Vrt1="<<foundV1<<" Vrt2="<<foundV2<<" dst="<<minDistVV);
233 if(probV<probVrtMergeLimit){ //--- If merged vertex is bad - try to remove the worst track
234 int pos=std::max_element(newvrt.chi2PerTrk.begin(),newvrt.chi2PerTrk.end())-newvrt.chi2PerTrk.begin();
235 newvrt.detachedTrack=newvrt.selTrk[pos];
236 newvrt.selTrk.erase(newvrt.selTrk.begin()+pos);
237 probV = refitVertex( newvrt, xAODwrk->listSelTracks, *state, false);
238 ATH_MSG_DEBUG( "Attempt to improve prob=" << probV);
239 }
240 if(probV>probVrtMergeLimit && newvrt.vertexMom.M()<m_vrtMassLimit){ // Good merged vertex found
241 newvrt.projectedVrt=MomProjDist(newvrt.vertex, primVrt, newvrt.vertexMom);
242 std::swap((*wrkVrtSet)[foundV1],newvrt);
243 (*wrkVrtSet)[foundV2].Good=false; //Drop vertex
244 (*wrkVrtSet)[foundV2].selTrk.clear(); //Clean dropped vertex
245 } else checkedDst.push_back(minDistVV);
246 minDistVV=minVrtVrtDist( wrkVrtSet.get(), foundV1, foundV2, checkedDst);
247 }
248//
249// Try to improve vertices with big Chi2 if something went wrong. Just precaution.
250 for(int iv=0; iv<(int)wrkVrtSet->size(); iv++) {
251 if(!(*wrkVrtSet)[iv].Good ) continue; //don't work on vertex which is already bad
252 if( (*wrkVrtSet)[iv].selTrk.size()<3 ) continue;
253 double tmpProb=TMath::Prob( (*wrkVrtSet)[iv].chi2, 2*(*wrkVrtSet)[iv].selTrk.size()-3 ); //Chi2 of the original vertex
254 if(tmpProb<m_globVrtProbCut){
255 ATH_MSG_DEBUG( "BAD vertex found prob=" << tmpProb);
256 tmpProb=improveVertexChi2( (*wrkVrtSet)[iv], xAODwrk->listSelTracks, *state, false);
257 (*wrkVrtSet)[iv].projectedVrt=MomProjDist((*wrkVrtSet)[iv].vertex, primVrt, (*wrkVrtSet)[iv].vertexMom);
258 }
259 }
260//
261//-Modify too heavy vertices
262 for(auto & iv : (*wrkVrtSet)){
263 if( iv.vertexMom.M()>m_vrtMassLimit ) {
264 ATH_MSG_DEBUG( "Heavy vertex found Mass=" << iv.vertexMom.M());
265 int it_bad=mostHeavyTrk(iv,xAODwrk->listSelTracks);
266 if(it_bad>-1){
267 iv.selTrk.erase( iv.selTrk.begin() + it_bad );
268 refitVertex(iv, xAODwrk->listSelTracks, *state, false);
269 iv.projectedVrt=MomProjDist(iv.vertex, primVrt, iv.vertexMom);
270 } } }
271//--------------------------------------------------------------------------------------------------------
272// Final vertex selection/cleaning
273//
274 double signif3D=0., signif2D=0.;
275
276//----- Vertices with >1 tracks
277 for(int iv=0; iv<(int)wrkVrtSet->size(); iv++) {
278 WrkVrt & curVrt=(*wrkVrtSet)[iv];
279 nth=(*wrkVrtSet)[iv].selTrk.size();
280 if(nth == 1) continue; // 1track vertices for later...
281 if(!curVrt.Good ) continue; //don't work on vertex which is already bad
282 (*wrkVrtSet)[iv].Good = false; /* Make all vertices bad for futher selection */
283 if(nth < 1) continue; /* Definitely bad vertices */
284 if((*wrkVrtSet)[iv].projectedVrt<0.) continue; /* Remove vertices behind primary one */
285 if( TMath::Prob( curVrt.chi2, 2*nth-3)<m_globVrtProbCut) continue; /* Bad Chi2 of refitted vertex */
286//-----------------------------------------------------------------------------------------
287 vrtVrtDist(primVrt,curVrt.vertex, curVrt.vertexCov, signif3D);
288 if(m_fillHist){
289 Hists& h = getHists();
290 if(nth==2 && curVrt.vertexCharge==0) h.m_hb_massPiPi1->Fill(curVrt.vertexMom.M(), m_w_1);
291 h.m_hb_sig3DTot->Fill( signif3D, m_w_1);
292 if(nth==2)h.m_hb_sig3D2tr->Fill( signif3D, m_w_1);
293 if(nth >2)h.m_hb_sig3DNtr->Fill( signif3D, m_w_1);
294 }
295//
296//--- Check V0s and conversions. Necessity must be checked in physics applications
297#if 0
298 if(nth==2 && curVrt.vertexCharge==0 && curVrt.detachedTrack<0){
299 double mass_PiPi = curVrt.vertexMom.M();
300 double mass_PPi = massV0(curVrt.trkAtVrt,m_massP,m_massPi);
301 double mass_EE = massV0(curVrt.trkAtVrt,m_massE,m_massE);
302 if(m_fillHist){
303 Hists& h = getHists();
304 h.m_hb_massPiPi->Fill( mass_PiPi, m_w_1);
305 h.m_hb_massPPi ->Fill( mass_PPi, m_w_1);
306 if( curVrt.vertex.perp()>20.)h.m_hb_massEE ->Fill( mass_EE, m_w_1);
307 }
308 if( std::abs(mass_PiPi-m_massK0) < 22.) continue;
309 if( std::abs(mass_PPi-m_massLam) < 8.) continue;
310 if( mass_EE < 60. && curVrt.vertex.perp() > 20.) continue;
311 }
312#endif
313//---
314 if(signif3D<m_selVrtSigCut) continue; //Main PV-SV distance quality cut
315 if(curVrt.vertex.perp() > m_maxSVRadiusCut) continue; // Too far from interaction point
316 curVrt.Good = true; /* Vertex is absolutely good */
317 }
318 if(msgLvl(MSG::DEBUG))printWrkSet(wrkVrtSet.get(),"Final Vertices");
319//-------------------------------------------
320// Debug ntuple filling and BDT application
321//
322 std::vector<double> impact,impactError;
323 for(int iv=0; iv<(int)wrkVrtSet->size(); iv++) {
324 WrkVrt & curVrt=(*wrkVrtSet)[iv];
325 nth=curVrt.selTrk.size();
326 if(!curVrt.Good || nth<2) continue;
327 double minPtT=1.e6, minSig3DT=1.e6, maxSig3DT=0.;
328 int ntrkBC=0,ntrkI=0,sumIBLHits=0,sumBLHits=0;
329 for(i=0;i<nth;i++) {
330 j=curVrt.selTrk[i]; //Track number
331 minPtT=std::min( minPtT, xAODwrk->listSelTracks[j]->pt());
332 m_fitSvc->VKalGetImpact(xAODwrk->listSelTracks[j], primVrt.position(), 1, impact, impactError);
333 double SigR2 = impact[0]*impact[0]/impactError[0];
334 double SigZ2 = impact[1]*impact[1]/impactError[2];
335 minSig3DT=std::min( minSig3DT, sqrt( SigR2 + SigZ2) );
336 maxSig3DT=std::max( maxSig3DT, sqrt( SigR2 + SigZ2) );
337 sumIBLHits += std::max(getIBLHit(xAODwrk->listSelTracks[j]),0);
338 sumBLHits += std::max(getBLHit(xAODwrk->listSelTracks[j]),0);
339 if(m_fillHist) {
340 ntrkBC += getIdHF(xAODwrk->listSelTracks[j]);
341 ntrkI += getG4Inter(xAODwrk->listSelTracks[j]);
342 }
343 }
344 float vProb=TMath::Prob(curVrt.chi2, 2*nth-3);
345 float cosSVPVM=projSV_PV(curVrt.vertex, primVrt, curVrt.vertexMom);
346 float vrtR=curVrt.vertex.perp();
347 TLorentzVector SVPV(curVrt.vertex.x()-primVrt.x(),curVrt.vertex.y()-primVrt.y(),curVrt.vertex.z()-primVrt.z(), 10.);
348 if(m_fillHist){
349 Hists& h = getHists();
350 if( nth>1 ){
351 vrtVrtDist(primVrt,curVrt.vertex, curVrt.vertexCov, signif3D);
352 float Dist2D=vrtVrtDist2D(primVrt,curVrt.vertex, curVrt.vertexCov, signif2D);
353 h.m_curTup->NVrtTrk [h.m_curTup->nNVrt] = nth;
354 h.m_curTup->NVrtTrkHF [h.m_curTup->nNVrt] = ntrkBC;
355 h.m_curTup->NVrtTrkI [h.m_curTup->nNVrt] = ntrkI;
356 h.m_curTup->NVrtProb [h.m_curTup->nNVrt] = vProb;
357 h.m_curTup->NVrtSig3D [h.m_curTup->nNVrt] = signif3D;
358 h.m_curTup->NVrtSig2D [h.m_curTup->nNVrt] = signif2D;
359 h.m_curTup->NVrtDist2D[h.m_curTup->nNVrt] = vrtR<20. ? Dist2D : vrtR;
360 h.m_curTup->NVrtM [h.m_curTup->nNVrt] = curVrt.vertexMom.M();
361 h.m_curTup->NVrtPt [h.m_curTup->nNVrt] = curVrt.vertexMom.Pt();
362 h.m_curTup->NVrtEta [h.m_curTup->nNVrt] = SVPV.Eta();
363 h.m_curTup->NVrtIBL [h.m_curTup->nNVrt] = sumIBLHits;
364 h.m_curTup->NVrtBL [h.m_curTup->nNVrt] = sumBLHits;
365 h.m_curTup->NVrtCosSPM[h.m_curTup->nNVrt] = cosSVPVM;
366 h.m_curTup->NVrtCh [h.m_curTup->nNVrt] = curVrt.vertexCharge;
367 h.m_curTup->NVMinPtT [h.m_curTup->nNVrt] = minPtT;
368 h.m_curTup->NVMinS3DT [h.m_curTup->nNVrt] = minSig3DT;
369 h.m_curTup->NVrtBDT [h.m_curTup->nNVrt] = 1.1;
370 h.m_curTup->NVrtHR1 [h.m_curTup->nNVrt] = xAODwrk->listSelTracks[curVrt.selTrk[0]]->radiusOfFirstHit();
371 h.m_curTup->NVrtHR2 [h.m_curTup->nNVrt] = xAODwrk->listSelTracks[curVrt.selTrk[1]]->radiusOfFirstHit();
372 if( h.m_curTup->nNVrt < DevTuple::maxNVrt-1 )h.m_curTup->nNVrt++;
373 }
374 }
375//--------Final selection of the remaining 2-track vertices
376// Should be tighter than initila 2-track vertex selection for compatibility graph construction
377//
378 if(nth==2){
379 float wgtSelect=-1.1;
380 xAOD::Vertex testV;
381 testV.makePrivateStore();
382 testV.setPosition(curVrt.vertex);
383 std::vector<float> testVcov(curVrt.vertexCov.begin(),curVrt.vertexCov.end());
384 testV.setCovariance(testVcov);
385 testV.setFitQuality(curVrt.chi2,1.);
386 bool acceptV=m_fin_v2trselector->isgood(std::make_pair(xAODwrk->listSelTracks[curVrt.selTrk[0]],
387 xAODwrk->listSelTracks[curVrt.selTrk[1]]),
388 testV,
389 std::make_pair(momAtVrt(curVrt.trkAtVrt[0]),momAtVrt(curVrt.trkAtVrt[1])), primVrt, wgtSelect);
390 curVrt.BDT=wgtSelect;
391 if(m_fillHist){
392 Hists& h = getHists();
393 h.m_hb_fakeSVBDT->Fill(wgtSelect,1.);
394 h.m_curTup->NVrtBDT[h.m_curTup->nNVrt-1] = wgtSelect;
395 }
396 if(!acceptV){
397 curVrt.Good = false; // Disable 2-track vertex with bad BDT score
398 if(m_multiWithOneTrkVrt){ // Check if linked 1-track vertex exists and disable it
399 for(auto it : curVrt.selTrk){
400 for(auto &vtmp : (*wrkVrtSet)){
401 if(vtmp.selTrk.size()!=1 || (!vtmp.Good)) continue;
402 if(it==vtmp.detachedTrack)vtmp.Good=false;
403 } } }
404 }
405 }
406 } //End vertex set loop
407//
408//-- Debug ntuple for 1track vertex is filled here
409//
411 Hists& h = getHists();
412 for(auto & vrt : (*wrkVrtSet)) {
413 if( !vrt.Good || vrt.selTrk.size() != 1 ) continue; // Good 1track vertices
414 const auto *xaodtp=xAODwrk->listSelTracks[vrt.selTrk[0]];
415 m_fitSvc->VKalGetImpact(xaodtp, primVrt.position(), 1, impact, impactError);
416 double SigR2 = std::abs(impact[0]*impact[0]/impactError[0]);
417 double SigZ2 = std::abs(impact[1]*impact[1]/impactError[2]);
418 float dist2D=vrtVrtDist2D(primVrt,vrt.vertex, vrt.vertexCov, signif2D);
419 h.m_curTup->NVrtTrk [h.m_curTup->nNVrt] = 1;
420 h.m_curTup->NVrtTrkHF [h.m_curTup->nNVrt] = getIdHF(xaodtp);
421 h.m_curTup->NVrtProb [h.m_curTup->nNVrt] = trkNPairs[vrt.selTrk[0]];
422 h.m_curTup->NVrtSig3D [h.m_curTup->nNVrt] = 0.;
423 h.m_curTup->NVrtSig2D [h.m_curTup->nNVrt] = signif2D;
424 h.m_curTup->NVrtDist2D[h.m_curTup->nNVrt] = dist2D;
425 h.m_curTup->NVrtM [h.m_curTup->nNVrt] = vrt.vertexMom.M();
426 h.m_curTup->NVrtPt [h.m_curTup->nNVrt] = vrt.vertexMom.Pt();
427 h.m_curTup->NVrtCosSPM[h.m_curTup->nNVrt] = 0.;
428 h.m_curTup->NVrtCh [h.m_curTup->nNVrt] = vrt.vertexCharge;
429 h.m_curTup->NVMinPtT [h.m_curTup->nNVrt] = xaodtp->pt();
430 h.m_curTup->NVMinS3DT [h.m_curTup->nNVrt] = sqrt(SigR2 + SigZ2);
431 h.m_curTup->NVrtBDT [h.m_curTup->nNVrt] = 1.1;
432 h.m_curTup->NVrtIBL [h.m_curTup->nNVrt] = std::max(getIBLHit(xaodtp),0);
433 h.m_curTup->NVrtBL [h.m_curTup->nNVrt] = std::max(getBLHit (xaodtp),0);
434 if( h.m_curTup->nNVrt < DevTuple::maxNVrt-1 )h.m_curTup->nNVrt++;
435 } }
436//-------------------------------------------
437//Sorting and check
438 std::multimap<double,WrkVrt,std::greater<double>> goodVertexMap;
439 int nNtrVrt=0;
440 for(auto & iv : (*wrkVrtSet) ) {
441 nth=iv.selTrk.size();
442 if(nth==1)iv.BDT=-2.; //To move 1-track vertices to the end of the vertex list later
443 double selector=iv.BDT;
444 if(nth==1) selector=iv.BDT+std::min(iv.vertexMom.Pt()*1.e-5,1.);
445 else if(nth>2) selector=iv.BDT+iv.vertexMom.M()*1.e-5;
446 if( iv.Good && nth>0 ) {
447 goodVertexMap.emplace(selector,iv); // add it and sort in the map
448 if(nth>1)nNtrVrt++;
449 }
450 }
451 if(nNtrVrt==0){ //--- No good vertices at all
452 if(m_fillHist) {
453 Hists& h = getHists();
454 h.m_curTup->nNVrt=0;
455 }
456 return finalVertices;
457 }
458//
459//-------------------------------------------
460// Final vertex refit for full covariance matrix and xAOD::Vertex creation
461//
462 static const SG::AuxElement::Decorator<float> wgtBDT("wgtBDT");
463 static const SG::AuxElement::Decorator<int> nTrksDec("nTracks");
464 static const SG::AuxElement::Decorator<int> vChrgTot("vCharge");
465 int n1trVrt=0; // Final number of good 1-track vertices
466 for(auto & iv : goodVertexMap){
467 WrkVrt & curVrt=iv.second;
468 nth=curVrt.selTrk.size();
469 xAODwrk->tmpListTracks.clear();
470 for(auto t : curVrt.selTrk)xAODwrk->tmpListTracks.push_back( xAODwrk->listSelTracks[t] );
471 if(m_fillHist){
472 Hists& h = getHists();
473 h.m_hb_totmass->Fill(curVrt.vertexMom.M(), m_w_1);
474 h.m_hb_r2d->Fill( curVrt.vertex.perp(), m_w_1);
475 }
476//--- Re-fit with full error matrix and xAOD::Vertex creation
477 xAOD::Vertex * tmpVertex=nullptr;
478 if(nth>1){ //-- Common case with full refit
479 tmpVertex=m_fitSvc->fit(xAODwrk->tmpListTracks,curVrt.vertex,*state);
480 } else if(nth==1){ //-- Special case for 1-track vertex
481 tmpVertex=new (std::nothrow) xAOD::Vertex();
482 if(!tmpVertex)continue;
483 tmpVertex->makePrivateStore();
484 tmpVertex->setPosition(curVrt.vertex);
485 std::vector<float> floatErrMtx(6);
486 for(int i=0; i<6; i++) floatErrMtx[i]=curVrt.vertexCov[i];
487 tmpVertex->setCovariance(floatErrMtx);
488 tmpVertex->setFitQuality(curVrt.chi2, (float)(nth*2-3));
489 std::vector<Trk::VxTrackAtVertex> & tmpVTAV=tmpVertex->vxTrackAtVertex(); tmpVTAV.clear();
490 AmgSymMatrix(5) CovMtxP;
491 CovMtxP.setIdentity();
492 Trk::Perigee * tmpMeasPer = new (std::nothrow) Trk::Perigee( 0.,0.,
493 curVrt.trkAtVrt[0][0],
494 curVrt.trkAtVrt[0][1],
495 curVrt.trkAtVrt[0][2],
497 std::move(CovMtxP) );
498 tmpVTAV.emplace_back( 1., tmpMeasPer );
500 const xAOD::TrackParticleContainer* cont = (const xAOD::TrackParticleContainer* ) (xAODwrk->tmpListTracks[0]->container() );
501 TEL.setStorableObject(*cont);
502 tmpVertex->addTrackAtVertex(TEL,1.);
503 n1trVrt++;
504 }
505 if(tmpVertex){
506 wgtBDT (*tmpVertex) =curVrt.BDT;
507 nTrksDec(*tmpVertex) =curVrt.selTrk.size();
508 vChrgTot(*tmpVertex) =curVrt.vertexCharge;
510 finalVertices.push_back(tmpVertex);
511 for (int ind=0; ind<nth; ind++) {
512 m_chi2_toSV(*xAODwrk->listSelTracks[curVrt.selTrk[ind]]) = curVrt.chi2PerTrk[ind] > FLT_MAX ? FLT_MAX : curVrt.chi2PerTrk[ind];
513 }
514 }
515 }
516 if(m_fillHist){
517 Hists& h = getHists();
518 h.m_hb_goodvrtN->Fill( finalVertices.size()+0.1, m_w_1);
519 h.m_hb_goodvrt1N->Fill( n1trVrt+0.1, m_w_1);
520 }
521//-----------------------------------------------------------------------------------
522// Saving of results
523//
524 return finalVertices;
525
526 }
527
528} //end namespace
529
530
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
static Double_t a
static Double_t sc
bool msgLvl(const MSG::Level lvl) const
Header file for AthHistogramAlgorithm.
double improveVertexChi2(WrkVrt &vertex, std::vector< const xAOD::TrackParticle * > &allTracks, Trk::IVKalState &istate, bool ifCovV0) const
void selGoodTrkParticle(workVectorArrxAOD *xAODwrk, const xAOD::Vertex &primVrt) const
double mergeAndRefitVertices(WrkVrt &v1, WrkVrt &v2, WrkVrt &newvrt, std::vector< const xAOD::TrackParticle * > &AllTrackList, Trk::IVKalState &istate, int robKey=0) const
Gaudi::Property< float > m_vrtMassLimit
double refitVertex(WrkVrt &Vrt, std::vector< const xAOD::TrackParticle * > &SelectedTracks, Trk::IVKalState &istate, bool ifCovV0) const
SG::AuxElement::Decorator< float > m_chi2_toSV
static double vrtVrtDist2D(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr, double &signif)
std::vector< xAOD::Vertex * > getVrtSecMulti(workVectorArrxAOD *inpParticlesxAOD, const xAOD::Vertex &primVrt, compatibilityGraph_t &compatibilityGraph) const
double refineVerticesWithCommonTracks(WrkVrt &v1, WrkVrt &v2, std::vector< const xAOD::TrackParticle * > &allTrackList, Trk::IVKalState &istate) const
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
Gaudi::Property< float > m_cutPt
Gaudi::Property< float > m_vertexMergeCut
Gaudi::Property< float > m_globVrtProbCut
static int nTrkCommon(std::vector< WrkVrt > *WrkVrtSet, int indexV1, int indexV2)
ToolHandle< Rec::ITwoTrackVertexSelector > m_fin_v2trselector
ToolHandle< Trk::TrkVKalVrtFitter > m_fitSvc
void select2TrVrt(std::vector< const xAOD::TrackParticle * > &SelectedTracks, const xAOD::Vertex &primVrt, std::map< long int, std::vector< double > > &vrt, compatibilityGraph_t &compatibilityGraph) const
static double massV0(const std::vector< std::vector< double > > &TrkAtVrt, double massP, double massPi)
Gaudi::Property< float > m_maxSVRadiusCut
ROOT::Math::PxPyPzEVector momAtVrt(const std::vector< double > &inpTrk) const
static double MomProjDist(const Amg::Vector3D &SV, const xAOD::Vertex &PV, const TLorentzVector &Direction)
void printWrkSet(const std::vector< WrkVrt > *WrkSet, const std::string &name) const
static int mostHeavyTrk(WrkVrt V, std::vector< const xAOD::TrackParticle * > AllTracks)
static double vrtVrtDist(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr, double &signif)
double minVrtVrtDist(std::vector< WrkVrt > *WrkVrtSet, int &indexV1, int &indexV2, std::vector< double > &check) const
Gaudi::Property< float > m_selVrtSigCut
Gaudi::Property< bool > m_fillHist
static std::vector< double > estimVrtPos(int nTrk, std::deque< long int > &selTrk, std::map< long int, std::vector< double > > &vrt)
Gaudi::Property< bool > m_multiWithOneTrkVrt
void makePrivateStore()
Create a new (empty) private store for this object.
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
Class describing the Line to which the Perigee refers to.
float z() const
Returns the z position.
void addTrackAtVertex(const ElementLink< TrackParticleContainer > &tr, float weight=1.0)
Add a new track to the vertex.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
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.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
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.
double chi2(TH1 *h0, TH1 *h1)
Gaudi Tools.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)
@ SecVtx
Secondary vertex.
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
std::vector< std::vector< double > > trkAtVrt
std::vector< const xAOD::TrackParticle * > listSelTracks
std::vector< const xAOD::TrackParticle * > inpTrk
std::vector< const xAOD::TrackParticle * > tmpListTracks