CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFDisplacedVertexFinder.cc
Go to the documentation of this file.
2 
4 
7 
11 
13 
16 
18 
19 #include "TMath.h"
20 
21 using namespace std;
22 using namespace reco;
23 
24 //for debug only
25 //#define PFLOW_DEBUG
26 
28  displacedVertexCandidates_( new PFDisplacedVertexCandidateCollection ),
29  displacedVertices_( new PFDisplacedVertexCollection ),
30  transvSize_(0.0),
31  longSize_(0.0),
32  primaryVertexCut_(0.0),
33  tobCut_(0.0),
34  tecCut_(0.0),
35  minAdaptWeight_(2.0),
36  debug_(false) {}
37 
38 
40 
41 void
43  const edm::Handle<reco::PFDisplacedVertexCandidateCollection>& displacedVertexCandidates) {
44 
45  if( displacedVertexCandidates_.get() ) {
47  }
48  else
50 
51 
52  if(displacedVertexCandidates.isValid()) {
53  for(unsigned i=0;i<displacedVertexCandidates->size(); i++) {
54  PFDisplacedVertexCandidateRef dvcref( displacedVertexCandidates, i);
55  displacedVertexCandidates_->push_back( (*dvcref));
56  }
57  }
58 
59 }
60 
61 
62 // -------- Main function which find vertices -------- //
63 
64 void
66 
67  if (debug_) cout << "========= Start Find Displaced Vertices =========" << endl;
68 
69  // The vertexCandidates have not been passed to the event
70  // So they need to be cleared is they are not empty
71  if( displacedVertices_.get() ) displacedVertices_->clear();
72  else
74 
75 
76  // Prepare the collections
77  PFDisplacedVertexSeedCollection tempDisplacedVertexSeeds;
78  tempDisplacedVertexSeeds.reserve(4*displacedVertexCandidates_->size());
79  PFDisplacedVertexCollection tempDisplacedVertices;
80  tempDisplacedVertices.reserve(4*displacedVertexCandidates_->size());
81 
82  if (debug_)
83  cout << "1) Parsing displacedVertexCandidates into displacedVertexSeeds" << endl;
84 
85  // 1) Parsing displacedVertexCandidates into displacedVertexSeeds which would
86  // be later used for vertex fitting
87 
88  int i = -1;
89 
90  for(IDVC idvc = displacedVertexCandidates_->begin();
91  idvc != displacedVertexCandidates_->end(); idvc++) {
92 
93  i++;
94  if (debug_) {
95  cout << "Analyse Vertex Candidate " << i << endl;
96  }
97 
98  findSeedsFromCandidate(*idvc, tempDisplacedVertexSeeds);
99 
100  }
101 
102  if (debug_) cout << "2) Merging Vertex Seeds" << endl;
103 
104  // 2) Some displacedVertexSeeds coming from different displacedVertexCandidates
105  // may be closed enough to be merged together. bLocked is an array which keeps the
106  // information about the seeds which are desabled.
107  vector<bool> bLockedSeeds;
108  bLockedSeeds.resize(tempDisplacedVertexSeeds.size());
109  mergeSeeds(tempDisplacedVertexSeeds, bLockedSeeds);
110 
111  if (debug_) cout << "3) Fitting Vertices From Seeds" << endl;
112 
113  // 3) Fit displacedVertices from displacedVertexSeeds
114  for(unsigned idv = 0; idv < tempDisplacedVertexSeeds.size(); idv++){
115 
116  if (!tempDisplacedVertexSeeds[idv].isEmpty() && !bLockedSeeds[idv]) {
117  PFDisplacedVertex displacedVertex;
118  bLockedSeeds[idv] = fitVertexFromSeed(tempDisplacedVertexSeeds[idv], displacedVertex);
119  if (!bLockedSeeds[idv]) tempDisplacedVertices.push_back(displacedVertex);
120  }
121  }
122 
123  if (debug_) cout << "4) Rejecting Bad Vertices and label them" << endl;
124 
125  // 4) Reject displaced vertices which may be considered as fakes
126  vector<bool> bLocked;
127  bLocked.resize(tempDisplacedVertices.size());
128  selectAndLabelVertices(tempDisplacedVertices, bLocked);
129 
130  if (debug_) cout << "5) Fill the Displaced Vertices" << endl;
131 
132  // 5) Fill the displacedVertex_ which would be transfered to the producer
133  displacedVertices_->reserve(tempDisplacedVertices.size());
134 
135  for(unsigned idv = 0; idv < tempDisplacedVertices.size(); idv++)
136  if (!bLocked[idv]) displacedVertices_->push_back(tempDisplacedVertices[idv]);
137 
138  if (debug_) cout << "========= End Find Displaced Vertices =========" << endl;
139 
140 
141 }
142 
143 // -------- Different steps of the finder algorithm -------- //
144 
145 
146 void
148 
149  const PFDisplacedVertexCandidate::DistMap r2Map = vertexCandidate.r2Map();
150  bool bNeedNewCandidate = false;
151 
152  tempDisplacedVertexSeeds.push_back( PFDisplacedVertexSeed() );
153 
154  IDVS idvc_current;
155 
156  for (PFDisplacedVertexCandidate::DistMap::const_iterator imap = r2Map.begin();
157  imap != r2Map.end(); imap++){
158 
159  unsigned ie1 = (*imap).second.first;
160  unsigned ie2 = (*imap).second.second;
161 
162  if (debug_) cout << "ie1 = " << ie1 << " ie2 = " << ie2 << " radius = " << sqrt((*imap).first) << endl;
163 
164  GlobalPoint dcaPoint = vertexCandidate.dcaPoint(ie1, ie2);
165  if (fabs(dcaPoint.x()) > 1e9) continue;
166 
167  bNeedNewCandidate = true;
168  for (idvc_current = tempDisplacedVertexSeeds.begin(); idvc_current != tempDisplacedVertexSeeds.end(); idvc_current++){
169  if ((*idvc_current).isEmpty()) {
170  bNeedNewCandidate = false;
171  break;
172  }
173  const GlobalPoint vertexPoint = (*idvc_current).seedPoint();
174  double Delta_Long = getLongDiff(vertexPoint, dcaPoint);
175  double Delta_Transv = getTransvDiff(vertexPoint, dcaPoint);
176  if (Delta_Long > longSize_) continue;
177  if (Delta_Transv > transvSize_) continue;
178  bNeedNewCandidate = false;
179  break;
180  }
181  if (bNeedNewCandidate) {
182  if (debug_) cout << "create new displaced vertex" << endl;
183  tempDisplacedVertexSeeds.push_back( PFDisplacedVertexSeed() );
184  idvc_current = tempDisplacedVertexSeeds.end();
185  idvc_current--;
186  bNeedNewCandidate = false;
187  }
188 
189 
190 
191  (*idvc_current).updateSeedPoint(dcaPoint, vertexCandidate.tref(ie1), vertexCandidate.tref(ie2));
192 
193 
194 
195  }
196 
197 
198 }
199 
200 
201 
202 
203 void
204 PFDisplacedVertexFinder::mergeSeeds(PFDisplacedVertexSeedCollection& tempDisplacedVertexSeeds, vector<bool>& bLocked){
205 
206  // loop over displaced vertex candidates
207  // and merge them if they are close to each other
208  for(unsigned idv_mother = 0;idv_mother < tempDisplacedVertexSeeds.size(); idv_mother++){
209  if (!bLocked[idv_mother]){
210 
211  for (unsigned idv_daughter = idv_mother+1;idv_daughter < tempDisplacedVertexSeeds.size(); idv_daughter++){
212 
213  if (!bLocked[idv_daughter]){
214  if (isCloseTo(tempDisplacedVertexSeeds[idv_mother], tempDisplacedVertexSeeds[idv_daughter])) {
215 
216  tempDisplacedVertexSeeds[idv_mother].mergeWith(tempDisplacedVertexSeeds[idv_daughter]);
217  bLocked[idv_daughter] = true;
218  if (debug_) cout << "Seeds " << idv_mother << " and " << idv_daughter << " merged" << endl;
219  }
220  }
221  }
222  }
223  }
224 
225 }
226 
227 
228 
229 
230 
231 
232 
233 
234 bool
236 
237 
238  if (debug_) cout << "== Start vertexing procedure ==" << endl;
239 
240 
241  // ---- Prepare transient track list ----
242 
243  set < TrackBaseRef, PFDisplacedVertexSeed::Compare > tracksToFit = displacedVertexSeed.elements();
244  GlobalPoint seedPoint = displacedVertexSeed.seedPoint();
245 
246  vector<TransientTrack> transTracks;
247  vector<TransientTrack> transTracksRaw;
248  vector<TrackBaseRef> transTracksRef;
249  vector<TrackBaseRef> transTracksRefRaw;
250 
251  transTracks.reserve(tracksToFit.size());
252  transTracksRaw.reserve(tracksToFit.size());
253  transTracksRef.reserve(tracksToFit.size());
254  transTracksRefRaw.reserve(tracksToFit.size());
255 
256 
257 
258  TransientVertex theVertexAdaptiveRaw;
259  TransientVertex theRecoVertex;
260 
261 
262  // ---- 1) Clean for potentially poor seeds ------- //
263  // --------------------------------------------- //
264 
265  if (tracksToFit.size() < 2) {
266  if (debug_) cout << "Only one to Fit Track" << endl;
267  return true;
268  }
269 
270  double rho = sqrt(seedPoint.x()*seedPoint.x()+seedPoint.y()*seedPoint.y());
271  double z = seedPoint.z();
272 
273  if (rho > tobCut_ || fabs(z) > tecCut_) {
274  if (debug_) cout << "Seed Point out of the tracker rho = " << rho << " z = "<< z << " nTracks = " << tracksToFit.size() << endl;
275  return true;
276  }
277 
278  if (debug_) displacedVertexSeed.Dump();
279 
280  int nStep45 = 0;
281  int nNotIterative = 0;
282 
283  // Fill vectors of TransientTracks and TrackRefs after applying preselection cuts.
284  for(IEset ie = tracksToFit.begin(); ie != tracksToFit.end(); ie++){
285  TransientTrack tmpTk( *((*ie).get()), magField_, globTkGeomHandle_);
286  transTracksRaw.push_back( tmpTk );
287  transTracksRefRaw.push_back( *ie );
288  switch((*ie)->algo()) {
291  case reco::TrackBase::rs:
293  nNotIterative++;
294  break;
299  break;
302  nStep45++;
303  break;
304  default:
305  nNotIterative++;
306  nStep45++; // why this should be increased for these cases?
307  break;
308  };
309 
310  }
311 
312  if (rho > 25 && nStep45 + nNotIterative < 1){
313  if (debug_) cout << "Seed point at rho > 25 cm but no step 4-5 tracks" << endl;
314  return true;
315  }
316 
317  // ----------------------------------------------- //
318  // ---- PRESELECT GOOD TRACKS FOR FINAL VERTEX --- //
319  // ----------------------------------------------- //
320 
321 
322 
323  // 1)If only two track are found do not prefit
324 
325 
326  if ( transTracksRaw.size() == 2 ){
327 
328  if (debug_) cout << "No raw fit done" << endl;
330  if (debug_)
331  cout << "Due to probably high pile-up conditions 2 track vertices switched off" << endl;
332  return true;
333 
334  }
335  GlobalError globalError;
336 
337  theVertexAdaptiveRaw = TransientVertex(seedPoint, globalError, transTracksRaw, 1.);
338 
339 
340 
341  } else {
342 
343 
344 
345  if (debug_) cout << "Raw fit done." << endl;
346 
352 
355 
356 
357  if ( transTracksRaw.size() < 1000 && transTracksRaw.size() > 3){
358 
359  if (debug_) cout << "First test with KFT" << endl;
360 
361  KalmanVertexFitter theKalmanFitter(true);
362  theVertexAdaptiveRaw = theKalmanFitter.vertex(transTracksRaw, seedPoint);
363 
364  if( !theVertexAdaptiveRaw.isValid() || theVertexAdaptiveRaw.totalChiSquared() < 0. ) {
365  if(debug_) cout << "Prefit failed : valid? " << theVertexAdaptiveRaw.isValid()
366  << " totalChi2 = " << theVertexAdaptiveRaw.totalChiSquared() << endl;
367  return true;
368  }
369 
370  if (debug_) cout << "We use KFT instead of seed point to set up a point for AVF "
371  << " x = " << theVertexAdaptiveRaw.position().x()
372  << " y = " << theVertexAdaptiveRaw.position().y()
373  << " z = " << theVertexAdaptiveRaw.position().z()
374  << endl;
375 
376  // To save time: reject the Displaced vertex if it is too close to the beam pipe.
377  // Frequently it is very big vertices, with some dosens of tracks
378 
379  Vertex vtx = theVertexAdaptiveRaw;
380  rho = vtx.position().rho();
381 
382  // cout << "primary vertex cut = " << primaryVertexCut_ << endl;
383 
384  if (rho < primaryVertexCut_ || rho > 100) {
385  if (debug_) cout << "KFT Vertex geometrically rejected with tracks #rho = " << rho << endl;
386  return true;
387  }
388 
389  // cout << "primary vertex cut = " << primaryVertexCut_ << " rho = " << rho << endl;
390 
391  theVertexAdaptiveRaw = theAdaptiveFitterRaw.vertex(transTracksRaw, theVertexAdaptiveRaw.position());
392 
393 
394  } else {
395 
396 
397  theVertexAdaptiveRaw = theAdaptiveFitterRaw.vertex(transTracksRaw, seedPoint);
398 
399  }
400 
401  if( !theVertexAdaptiveRaw.isValid() || theVertexAdaptiveRaw.totalChiSquared() < 0. ) {
402  if(debug_) cout << "Fit failed : valid? " << theVertexAdaptiveRaw.isValid()
403  << " totalChi2 = " << theVertexAdaptiveRaw.totalChiSquared() << endl;
404  return true;
405  }
406 
407  // To save time: reject the Displaced vertex if it is too close to the beam pipe.
408  // Frequently it is very big vertices, with some dosens of tracks
409 
410  Vertex vtx = theVertexAdaptiveRaw;
411  rho = vtx.position().rho();
412 
413  if (rho < primaryVertexCut_) {
414  if (debug_) cout << "Vertex " << " geometrically rejected with " << transTracksRaw.size() << " tracks #rho = " << rho << endl;
415  return true;
416  }
417 
418 
419  }
420 
421 
422 
423  // ---- Remove tracks with small weight or
424  // big first (last) hit_to_vertex distance
425  // and then refit ---- //
426 
427 
428 
429  for (unsigned i = 0; i < transTracksRaw.size(); i++) {
430 
431  if (debug_) cout << "Raw Weight track " << i << " = " << theVertexAdaptiveRaw.trackWeight(transTracksRaw[i]) << endl;
432 
433  if (theVertexAdaptiveRaw.trackWeight(transTracksRaw[i]) > minAdaptWeight_){
434 
435  PFTrackHitFullInfo pattern = hitPattern_.analyze(tkerGeomHandle_, transTracksRefRaw[i], theVertexAdaptiveRaw);
436 
437  PFDisplacedVertex::VertexTrackType vertexTrackType = getVertexTrackType(pattern);
438 
439  if (vertexTrackType != PFDisplacedVertex::T_NOT_FROM_VERTEX){
440 
441  bool bGoodTrack = helper_.isTrackSelected(transTracksRaw[i].track(), vertexTrackType);
442 
443  if (bGoodTrack){
444  transTracks.push_back(transTracksRaw[i]);
445  transTracksRef.push_back(transTracksRefRaw[i]);
446  } else {
447  if (debug_)
448  cout << "Track rejected nChi2 = " << transTracksRaw[i].track().normalizedChi2()
449  << " pt = " << transTracksRaw[i].track().pt()
450  << " dxy (wrt (0,0,0)) = " << transTracksRaw[i].track().dxy()
451  << " nHits = " << transTracksRaw[i].track().numberOfValidHits()
452  << " nOuterHits = " << transTracksRaw[i].track().hitPattern().numberOfHits(HitPattern::MISSING_OUTER_HITS) << endl;
453  }
454  } else {
455 
456  if (debug_){
457  cout << "Remove track because too far away from the vertex:" << endl;
458  }
459 
460  }
461 
462  }
463 
464  }
465 
466 
467 
468  if (debug_) cout << "All Tracks " << transTracksRaw.size()
469  << " with good weight " << transTracks.size() << endl;
470 
471 
472  // ---- Refit ---- //
473  FitterType vtxFitter = F_NOTDEFINED;
474 
475  if (transTracks.size() < 2) return true;
476  else if (transTracks.size() == 2){
477 
479  if (debug_)
480  cout << "Due to probably high pile-up conditions 2 track vertices switched off" << endl;
481  return true;
482  }
483  vtxFitter = F_KALMAN;
484  }
485  else if (transTracks.size() > 2 && transTracksRaw.size() > transTracks.size())
486  vtxFitter = F_ADAPTIVE;
487  else if (transTracks.size() > 2 && transTracksRaw.size() == transTracks.size())
488  vtxFitter = F_DONOTREFIT;
489  else return true;
490 
491  if (debug_) cout << "Vertex Fitter " << vtxFitter << endl;
492 
493  if(vtxFitter == F_KALMAN){
494 
495  KalmanVertexFitter theKalmanFitter(true);
496  theRecoVertex = theKalmanFitter.vertex(transTracks, seedPoint);
497 
498  } else if(vtxFitter == F_ADAPTIVE){
499 
500  AdaptiveVertexFitter theAdaptiveFitter(
506 
507  theRecoVertex = theAdaptiveFitter.vertex(transTracks, seedPoint);
508 
509  } else if (vtxFitter == F_DONOTREFIT) {
510  theRecoVertex = theVertexAdaptiveRaw;
511  } else {
512  return true;
513  }
514 
515 
516  // ---- Check if the fitted vertex is valid ---- //
517 
518  if( !theRecoVertex.isValid() || theRecoVertex.totalChiSquared() < 0. ) {
519  if (debug_) cout << "Refit failed : valid? " << theRecoVertex.isValid()
520  << " totalChi2 = " << theRecoVertex.totalChiSquared() << endl;
521  return true;
522  }
523 
524  // ---- Create vertex ---- //
525 
526  Vertex theRecoVtx = theRecoVertex;
527 
528  double chi2 = theRecoVtx.chi2();
529  double ndf = theRecoVtx.ndof();
530 
531 
532  if (chi2 > TMath::ChisquareQuantile(0.95, ndf)) {
533  if (debug_)
534  cout << "Rejected because of chi2 = " << chi2 << " ndf = " << ndf << " confid. level: " << TMath::ChisquareQuantile(0.95, ndf) << endl;
535  return true;
536  }
537 
538 
539 
540  // ---- Create and fill vector of refitted TransientTracks ---- //
541 
542  // -----------------------------------------------//
543  // ---- Prepare and Fill the Displaced Vertex ----//
544  // -----------------------------------------------//
545 
546 
547  displacedVertex = (PFDisplacedVertex) theRecoVtx;
548  displacedVertex.removeTracks();
549 
550  for(unsigned i = 0; i < transTracks.size();i++) {
551 
552  PFTrackHitFullInfo pattern = hitPattern_.analyze(tkerGeomHandle_, transTracksRef[i], theRecoVertex);
553 
554  PFDisplacedVertex::VertexTrackType vertexTrackType = getVertexTrackType(pattern);
555 
556  Track refittedTrack;
557  float weight = theRecoVertex.trackWeight(transTracks[i]);
558 
559 
560  if (weight < minAdaptWeight_) continue;
561 
562  try{
563  refittedTrack = theRecoVertex.refittedTrack(transTracks[i]).track();
564  }catch( cms::Exception& exception ){
565  continue;
566  }
567 
568  if (debug_){
569  cout << "Vertex Track Type = " << vertexTrackType << endl;
570 
571  cout << "nHitBeforeVertex = " << pattern.first.first
572  << " nHitAfterVertex = " << pattern.second.first
573  << " nMissHitBeforeVertex = " << pattern.first.second
574  << " nMissHitAfterVertex = " << pattern.second.second
575  << " Weight = " << weight << endl;
576  }
577 
578  displacedVertex.addElement(transTracksRef[i],
579  refittedTrack,
580  pattern, vertexTrackType, weight);
581 
582  }
583 
584  displacedVertex.setPrimaryDirection(helper_.primaryVertex());
585  displacedVertex.calcKinematics();
586 
587 
588 
589  if (debug_) cout << "== End vertexing procedure ==" << endl;
590 
591  return false;
592 
593 }
594 
595 
596 
597 
598 
599 
600 void
601 PFDisplacedVertexFinder::selectAndLabelVertices(PFDisplacedVertexCollection& tempDisplacedVertices, vector <bool>& bLocked){
602 
603  if (debug_) cout << " 4.1) Reject vertices " << endl;
604 
605  for(unsigned idv = 0; idv < tempDisplacedVertices.size(); idv++){
606 
607 
608  // ---- We accept a vertex only if it is not in TOB in the barrel
609  // and not in TEC in the end caps
610  // and not too much before the first pixel layer ---- //
611 
612  const float rho = tempDisplacedVertices[idv].position().rho();
613  const float z = tempDisplacedVertices[idv].position().z();
614 
615  if (rho > tobCut_ || rho < primaryVertexCut_ || fabs(z) > tecCut_) {
616  if (debug_) cout << "Vertex " << idv
617  << " geometrically rejected #rho = " << rho
618  << " z = " << z << endl;
619 
620  bLocked[idv] = true;
621 
622  continue;
623  }
624 
625  unsigned nPrimary = tempDisplacedVertices[idv].nPrimaryTracks();
626  unsigned nMerged = tempDisplacedVertices[idv].nMergedTracks();
627  unsigned nSecondary = tempDisplacedVertices[idv].nSecondaryTracks();
628 
629  if (nPrimary + nMerged > 1) {
630  bLocked[idv] = true;
631  if (debug_) cout << "Vertex " << idv
632  << " rejected because two primary or merged tracks" << endl;
633 
634 
635  }
636 
637  if (nPrimary + nMerged + nSecondary < 2){
638  bLocked[idv] = true;
639  if (debug_) cout << "Vertex " << idv
640  << " rejected because only one track related to the vertex" << endl;
641  }
642 
643 
644  }
645 
646 
647  if (debug_) cout << " 4.2) Check for common vertices" << endl;
648 
649  // ---- Among the remaining vertex we shall remove one
650  // of those which have two common tracks ---- //
651 
652  for(unsigned idv_mother = 0; idv_mother < tempDisplacedVertices.size(); idv_mother++){
653  for(unsigned idv_daughter = idv_mother+1;
654  idv_daughter < tempDisplacedVertices.size(); idv_daughter++){
655 
656  if(!bLocked[idv_daughter] && !bLocked[idv_mother]){
657 
658  const unsigned commonTrks = commonTracks(tempDisplacedVertices[idv_daughter], tempDisplacedVertices[idv_mother]);
659 
660  if (commonTrks > 1) {
661 
662  if (debug_) cout << "Vertices " << idv_daughter << " and " << idv_mother
663  << " has many common tracks" << endl;
664 
665  // We keep the vertex vertex which contains the most of the tracks
666 
667  const int mother_size = tempDisplacedVertices[idv_mother].nTracks();
668  const int daughter_size = tempDisplacedVertices[idv_daughter].nTracks();
669 
670  if (mother_size > daughter_size) bLocked[idv_daughter] = true;
671  else if (mother_size < daughter_size) bLocked[idv_mother] = true;
672  else {
673 
674  // If they have the same size we keep the vertex with the best normalised chi2
675 
676  const float mother_normChi2 = tempDisplacedVertices[idv_mother].normalizedChi2();
677  const float daughter_normChi2 = tempDisplacedVertices[idv_daughter].normalizedChi2();
678  if (mother_normChi2 < daughter_normChi2) bLocked[idv_daughter] = true;
679  else bLocked[idv_mother] = true;
680  }
681 
682  }
683  }
684  }
685  }
686 
687  for(unsigned idv = 0; idv < tempDisplacedVertices.size(); idv++)
688  if ( !bLocked[idv] ) bLocked[idv] = rejectAndLabelVertex(tempDisplacedVertices[idv]);
689 
690 
691 }
692 
693 bool
695 
697  dv.setVertexType(vertexType);
698 
699  return dv.isFake();
700 
701 }
702 
703 
705 
706 bool
708 
709  const GlobalPoint vP1 = dv1.seedPoint();
710  const GlobalPoint vP2 = dv2.seedPoint();
711 
712  double Delta_Long = getLongDiff(vP1, vP2);
713  if (Delta_Long > longSize_) return false;
714  double Delta_Transv = getTransvDiff(vP1, vP2);
715  if (Delta_Transv > transvSize_) return false;
716  // if (Delta_Long < longSize_ && Delta_Transv < transvSize_) isCloseTo = true;
717 
718  return true;
719 
720 }
721 
722 
723 double
725 
726  Basic3DVector<double>vRef(Ref);
727  Basic3DVector<double>vToProject(ToProject);
728  return fabs((vRef.dot(vToProject)-vRef.mag2())/vRef.mag());
729 
730 }
731 
732 double
734 
735  Basic3DVector<double>vRef(Ref);
736  Basic3DVector<double>vToProject(ToProject);
737  return (vRef.dot(vToProject))/vRef.mag();
738 
739 
740 }
741 
742 
743 double
745 
746  Basic3DVector<double>vRef(Ref);
747  Basic3DVector<double>vToProject(ToProject);
748  return fabs(vRef.cross(vToProject).mag()/vRef.mag());
749 
750 }
751 
752 
755 
756  unsigned int nHitBeforeVertex = pairTrackHitInfo.first.first;
757  unsigned int nHitAfterVertex = pairTrackHitInfo.second.first;
758 
759  unsigned int nMissHitBeforeVertex = pairTrackHitInfo.first.second;
760  unsigned int nMissHitAfterVertex = pairTrackHitInfo.second.second;
761 
762  // For the moment those definitions are given apriori a more general study would be useful to refine those criteria
763 
764  if (nHitBeforeVertex <= 1 && nHitAfterVertex >= 3 && nMissHitAfterVertex <= 1)
765  return PFDisplacedVertex::T_FROM_VERTEX;
766  else if (nHitBeforeVertex >= 3 && nHitAfterVertex <= 1 && nMissHitBeforeVertex <= 1)
767  return PFDisplacedVertex::T_TO_VERTEX;
768  else if ((nHitBeforeVertex >= 2 && nHitAfterVertex >= 3)
769  ||
770  (nHitBeforeVertex >= 3 && nHitAfterVertex >= 2))
771  return PFDisplacedVertex::T_MERGED;
772  else
773  return PFDisplacedVertex::T_NOT_FROM_VERTEX;
774 }
775 
776 
778 
779  vector<Track> vt1 = v1.refittedTracks();
780  vector<Track> vt2 = v2.refittedTracks();
781 
782  unsigned commonTracks = 0;
783 
784  for ( unsigned il1 = 0; il1 < vt1.size(); il1++){
785  unsigned il1_idx = v1.originalTrack(vt1[il1]).key();
786 
787  for ( unsigned il2 = 0; il2 < vt2.size(); il2++)
788  if (il1_idx == v2.originalTrack(vt2[il2]).key()) {commonTracks++; break;}
789 
790  }
791 
792  return commonTracks;
793 
794 }
795 
796 std::ostream& operator<<(std::ostream& out, const PFDisplacedVertexFinder& a) {
797 
798  if(! out) return out;
799  out << setprecision(3) << setw(5) << endl;
800  out << "" << endl;
801  out << " ====================================== " << endl;
802  out << " ====== Displaced Vertex Finder ======= " << endl;
803  out << " ====================================== " << endl;
804  out << " " << endl;
805 
806  a.helper_.Dump();
807  out << "" << endl
808  << " Adaptive Vertex Fitter parameters are :"<< endl
809  << " sigmacut = " << a.sigmacut_ << " T_ini = "
810  << a.t_ini_ << " ratio = " << a.ratio_ << endl << endl;
811 
812  const std::auto_ptr< reco::PFDisplacedVertexCollection >& displacedVertices_
813  = a.displacedVertices();
814 
815 
816  if(!displacedVertices_.get() ) {
817  out<<"displacedVertex already transfered"<<endl;
818  }
819  else{
820 
821  out<<"Number of displacedVertices found : "<< displacedVertices_->size()<<endl<<endl;
822 
823  int i = -1;
824 
825  for(PFDisplacedVertexFinder::IDV idv = displacedVertices_->begin();
826  idv != displacedVertices_->end(); idv++){
827  i++;
828  out << i << " "; idv->Dump(); out << "" << endl;
829  }
830  }
831 
832  return out;
833 }
double getTransvDiff(const GlobalPoint &, const GlobalPoint &) const
int i
Definition: DBlmapReader.cc:9
std::vector< PFDisplacedVertex > PFDisplacedVertexCollection
collection of PFDisplacedVertex objects
std::vector< PFDisplacedVertexCandidate > PFDisplacedVertexCandidateCollection
collection of PFDisplacedVertexCandidate objects
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
A block of tracks linked together.
void setPrimaryDirection(const math::XYZPoint &pvtx)
reco::PFDisplacedVertexSeedCollection::iterator IDVS
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:86
std::auto_ptr< reco::PFDisplacedVertexCollection > displacedVertices_
edm::ESHandle< GlobalTrackingGeometry > globTkGeomHandle_
Tracker geometry for discerning hit positions.
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
float totalChiSquared() const
void Dump(std::ostream &out=std::cout) const
Basic3DVector cross(const Basic3DVector &lh) const
Vector product, or &quot;cross&quot; product, with a vector of same type.
T y() const
Definition: PV3DBase.h:63
bool rejectAndLabelVertex(reco::PFDisplacedVertex &dv)
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
std::map< float, std::pair< int, int > > DistMap
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:149
void findDisplacedVertices()
-----— Main function which find vertices -----— ///
bool debug_
If true, debug printouts activated.
const Point & position() const
position
Definition: Vertex.h:106
const std::set< TrackBaseRef, Compare > & elements() const
std::pair< PFTrackHitInfo, PFTrackHitInfo > PFTrackHitFullInfo
bool fitVertexFromSeed(reco::PFDisplacedVertexSeed &, reco::PFDisplacedVertex &)
Fit one by one the vertex points with associated tracks to get displaced vertices.
unsigned commonTracks(const reco::PFDisplacedVertex &, const reco::PFDisplacedVertex &) const
void mergeSeeds(reco::PFDisplacedVertexSeedCollection &, std::vector< bool > &bLocked)
Sometimes two vertex candidates can be quite close and coming from the same vertex.
std::set< reco::TrackBaseRef >::iterator IEset
-----— Useful Types -----— ///
PFTrackHitFullInfo analyze(edm::ESHandle< TrackerGeometry >, const reco::TrackBaseRef track, const TransientVertex &vert)
reco::TransientTrack refittedTrack(const reco::TransientTrack &track) const
void Dump(std::ostream &out=std::cout) const
cout function
void addElement(const TrackBaseRef &r, const Track &refTrack, const PFTrackHitFullInfo &hitInfo, VertexTrackType trackType=T_NOT_FROM_VERTEX, float w=1.0)
Add a new track to the vertex.
bool isTrackSelected(const reco::Track &trk, const reco::PFDisplacedVertex::VertexTrackType vertexTrackType) const
Select tracks tool.
edm::ESHandle< TrackerGeometry > tkerGeomHandle_
doc?
const TrackBaseRef & tref(unsigned ie) const
void setVertexType(VertexType vertexType)
Set the type of this vertex.
const GlobalPoint dcaPoint(unsigned ie1, unsigned ie2) const
std::auto_ptr< reco::PFDisplacedVertexCandidateCollection > displacedVertexCandidates_
-----— Members -----— ///
std::vector< PFDisplacedVertexSeed > PFDisplacedVertexSeedCollection
collection of PFDisplacedVertexSeed objects
size_t key() const
Definition: RefToBase.h:241
T sqrt(T t)
Definition: SSEVec.h:48
reco::PFDisplacedVertex::VertexType identifyVertex(const reco::PFDisplacedVertex &v) const
Vertex identification tool.
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
void selectAndLabelVertices(reco::PFDisplacedVertexCollection &, std::vector< bool > &)
Remove potentially fakes displaced vertices.
const MagneticField * magField_
to be able to extrapolate tracks f
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
double chi2() const
chi-squares
Definition: Vertex.h:95
double getLongProj(const GlobalPoint &, const GlobalVector &) const
DistMap r2Map() const
--—— Provide useful information --—— ///
reco::PFDisplacedVertexCandidateCollection::iterator IDVC
bool isValid() const
Definition: HandleBase.h:75
void setInput(const edm::Handle< reco::PFDisplacedVertexCandidateCollection > &)
Set input collections of tracks.
float trackWeight(const reco::TransientTrack &track) const
double ndof() const
Definition: Vertex.h:102
const std::auto_ptr< reco::PFDisplacedVertexCollection > & displacedVertices() const
tuple out
Definition: dbtoconf.py:99
void findSeedsFromCandidate(reco::PFDisplacedVertexCandidate &, reco::PFDisplacedVertexSeedCollection &)
--—— Different steps of the finder algorithm --—— ///
const GlobalPoint & seedPoint() const
double transvSize_
--—— Parameters --—— ///
Block of elements.
const Track & track() const
double sigmacut_
Adaptive Vertex Fitter parameters.
reco::PFDisplacedVertexCollection::iterator IDV
math::XYZPoint primaryVertex() const
Set Vertex direction using the primary vertex.
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
double getLongDiff(const GlobalPoint &, const GlobalPoint &) const
void removeTracks()
Definition: Vertex.cc:64
volatile std::atomic< bool > shutdown_flag false
bool isCloseTo(const reco::PFDisplacedVertexSeed &, const reco::PFDisplacedVertexSeed &) const
-----— Tools -----— ///
T x() const
Definition: PV3DBase.h:62
bool isValid() const
reco::PFDisplacedVertex::VertexTrackType getVertexTrackType(PFTrackHitFullInfo &) const
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
PFDisplacedVertexHelper helper_
T dot(const Basic3DVector &rh) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.