CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PF_PU_AssoMapAlgos.cc
Go to the documentation of this file.
1 
3 
5 
7 
12 
14 
17 
18 using namespace edm;
19 using namespace std;
20 using namespace reco;
21 
22 /*************************************************************************************/
23 /* dedicated constructor for the algorithms */
24 /*************************************************************************************/
25 
27  : maxNumWarnings_(3),
28  numWarnings_(0)
29 {
30 
31  input_MaxNumAssociations_ = iConfig.getParameter<int>("MaxNumberOfAssociations");
32 
33  token_VertexCollection_= iC.consumes<VertexCollection>(iConfig.getParameter<InputTag>("VertexCollection"));
34 
35  token_BeamSpot_= iC.consumes<BeamSpot>(iConfig.getParameter<InputTag>("BeamSpot"));
36 
37  input_doReassociation_= iConfig.getParameter<bool>("doReassociation");
38  cleanedColls_ = iConfig.getParameter<bool>("GetCleanedCollections");
39 
40  ConversionsCollectionToken_= iC.consumes<ConversionCollection>(iConfig.getParameter<InputTag>("ConversionsCollection"));
41 
44 
46 
47  input_FinalAssociation_= iConfig.getUntrackedParameter<int>("FinalAssociation", 0);
48 
49  ignoremissingpfcollection_ = iConfig.getParameter<bool>("ignoreMissingCollection");
50 
51  input_nTrack_ = iConfig.getParameter<double>("nTrackWeight");
52 
53 }
54 
55 /*************************************************************************************/
56 /* get all needed collections at the beginning */
57 /*************************************************************************************/
58 
59 void
61 {
62 
63  //get the offline beam spot
65 
66  //get the conversion collection for the gamma conversions
69 
70  //get the vertex composite candidate collection for the Kshort's
73 
74  //get the vertex composite candidate collection for the Lambda's
77 
78  //get the displaced vertex collection for nuclear interactions
79  //create a new bool, true if no displaced vertex collection is in the event, mostly for AOD
80  missingColls = false;
83 
84  missingColls = true;
85 
86  if ( numWarnings_ < maxNumWarnings_ ) {
87  LogWarning("PF_PU_AssoMapAlgos::GetInputCollections")
88  << "No Extra objects available in input file --> skipping reconstruction of displaced vertices !!" << endl;
89  ++numWarnings_;
90  }
91 
92  }
93  } else {
94 
96 
97  }
98 
99  //get the input vertex collection
101 
102  iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
103 
104 }
105 
106 /*************************************************************************************/
107 /* create the track to vertex association map */
108 /*************************************************************************************/
109 std::auto_ptr<TrackToVertexAssMap>
111 {
112 
113  auto_ptr<TrackToVertexAssMap> track2vertex(new TrackToVertexAssMap(vtxcollH, trkcollH));
114 
115  int num_vertices = vtxcollH->size();
116  if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
117 
118  //loop over all tracks of the track collection
119  for ( size_t idxTrack = 0; idxTrack < trkcollH->size(); ++idxTrack ) {
120 
121  TrackRef trackref = TrackRef(trkcollH, idxTrack);
122 
123  TransientTrack transtrk(trackref, &(*bFieldH) );
124  transtrk.setBeamSpot(*beamspotH);
125  transtrk.setES(iSetup);
126 
127  vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
128 
129  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
130 
131  VertexStepPair assocVtx = FindAssociation(trackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
132  int step = assocVtx.second;
133  double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
134 
135  int quality = DefineQuality(assoc_ite, step, distance);
136 
137  //std::cout << "associating track: Pt = " << trackref->pt() << ","
138  // << " eta = " << trackref->eta() << ", phi = " << trackref->phi()
139  // << " to vertex: z = " << associatedVertex.first->position().z() << " with quality q = " << quality << std::endl;
140 
141 
142  // Insert the best vertex and the pair of track and the quality of this association in the map
143  track2vertex->insert( assocVtx.first, make_pair(trackref, quality) );
144 
145  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
146 
147  }
148 
149  delete vtxColl_help;
150 
151  }
152 
153  return track2vertex;
154 
155 }
156 
157 /*************************************************************************************/
158 /* create the vertex to track association map */
159 /*************************************************************************************/
160 
161 std::auto_ptr<VertexToTrackAssMap>
163 {
164 
165  auto_ptr<VertexToTrackAssMap> vertex2track(new VertexToTrackAssMap(trkcollH, vtxcollH));
166 
167  int num_vertices = vtxcollH->size();
168  if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
169 
170  //loop over all tracks of the track collection
171  for ( size_t idxTrack = 0; idxTrack < trkcollH->size(); ++idxTrack ) {
172 
173  TrackRef trackref = TrackRef(trkcollH, idxTrack);
174 
175  TransientTrack transtrk(trackref, &(*bFieldH) );
176  transtrk.setBeamSpot(*beamspotH);
177  transtrk.setES(iSetup);
178 
179  vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
180 
181  for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
182 
183  VertexStepPair assocVtx = FindAssociation(trackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
184  int step = assocVtx.second;
185  double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
186 
187  int quality = DefineQuality(assoc_ite, step, distance);
188 
189  //std::cout << "associating track: Pt = " << trackref->pt() << ","
190  // << " eta = " << trackref->eta() << ", phi = " << trackref->phi()
191  // << " to vertex: z = " << associatedVertex.first->position().z() << " with quality q = " << quality << std::endl;
192 
193  // Insert the best vertex and the pair of track and the quality of this association in the map
194  vertex2track->insert( trackref, make_pair(assocVtx.first, quality) );
195 
196  PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
197 
198  }
199 
200  delete vtxColl_help;
201 
202  }
203 
204  return vertex2track;
205 
206 }
207 
208 /*****************************************************************************************/
209 /* function to sort the vertices in the AssociationMap by the sum of (pT - pT_Error)**2 */
210 /*****************************************************************************************/
211 
212 auto_ptr<TrackToVertexAssMap>
214 {
215  //create a new TrackVertexAssMap for the Output which will be sorted
216  auto_ptr<TrackToVertexAssMap> trackvertexassOutput(new TrackToVertexAssMap() );
217 
218  //Create and fill a vector of pairs of vertex and the summed (pT-pT_Error)**2 of the tracks associated to the vertex
219  VertexPtsumVector vertexptsumvector;
220 
221  //loop over all vertices in the association map
222  for(TrackToVertexAssMap::const_iterator assomap_ite=trackvertexassInput->begin(); assomap_ite!=trackvertexassInput->end(); assomap_ite++){
223 
224  const VertexRef assomap_vertexref = assomap_ite->key;
225  const TrackQualityPairVector trckcoll = assomap_ite->val;
226 
227  float ptsum = 0;
228 
229  TrackRef trackref;
230 
231  //get the tracks associated to the vertex and calculate the manipulated pT**2
232  for(unsigned int trckcoll_ite=0; trckcoll_ite<trckcoll.size(); trckcoll_ite++){
233 
234  trackref = trckcoll[trckcoll_ite].first;
235  int quality = trckcoll[trckcoll_ite].second;
236 
237  if ( quality<=2 ) continue;
238 
239  double man_pT = trackref->pt() - trackref->ptError();
240  if(man_pT>0.) ptsum+=man_pT*man_pT;
241 
242  }
243 
244  vertexptsumvector.push_back(make_pair(assomap_vertexref,ptsum));
245 
246  }
247 
248  while (vertexptsumvector.size()!=0){
249 
250  VertexRef vertexref_highestpT;
251  float highestpT = 0.;
252  int highestpT_index = 0;
253 
254  for(unsigned int vtxptsumvec_ite=0; vtxptsumvec_ite<vertexptsumvector.size(); vtxptsumvec_ite++){
255 
256  if(vertexptsumvector[vtxptsumvec_ite].second>highestpT){
257 
258  vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
259  highestpT = vertexptsumvector[vtxptsumvec_ite].second;
260  highestpT_index = vtxptsumvec_ite;
261 
262  }
263 
264  }
265 
266  //loop over all vertices in the association map
267  for(TrackToVertexAssMap::const_iterator assomap_ite=trackvertexassInput->begin(); assomap_ite!=trackvertexassInput->end(); assomap_ite++){
268 
269  const VertexRef assomap_vertexref = assomap_ite->key;
270  const TrackQualityPairVector trckcoll = assomap_ite->val;
271 
272  //if the vertex from the association map the vertex with the highest manipulated pT
273  //insert all associated tracks in the output Association Map
274  if(assomap_vertexref==vertexref_highestpT)
275  for(unsigned int trckcoll_ite=0; trckcoll_ite<trckcoll.size(); trckcoll_ite++)
276  trackvertexassOutput->insert(assomap_vertexref,trckcoll[trckcoll_ite]);
277 
278  }
279 
280  vertexptsumvector.erase(vertexptsumvector.begin()+highestpT_index);
281 
282  }
283 
284  return trackvertexassOutput;
285 
286 }
287 
288 /********************/
289 /* */
290 /* Member Functions */
291 /* */
292 /********************/
293 
294 /*************************************************************************************/
295 /* create helping vertex vector to remove associated vertices */
296 /*************************************************************************************/
297 
298 std::vector<reco::VertexRef>*
300 {
301 
302  vector<VertexRef>* output = new vector<VertexRef>();
303 
304  for(unsigned int index_vtx=0; index_vtx<vtxcollH->size(); ++index_vtx){
305 
306  VertexRef vertexref(vtxcollH,index_vtx);
307 
308  output->push_back(vertexref);
309 
310  }
311 
312  return output;
313 
314 }
315 
316 /****************************************************************************/
317 /* erase one vertex from the vertex vector */
318 /****************************************************************************/
319 
320 void
321 PF_PU_AssoMapAlgos::EraseVertex(std::vector<reco::VertexRef>* vtxcollV, reco::VertexRef toErase)
322 {
323 
324  for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
325 
326  VertexRef vertexref = vtxcollV->at(index_vtx);
327 
328  if ( vertexref == toErase ){
329  vtxcollV->erase(vtxcollV->begin()+index_vtx);
330  break;
331  }
332 
333  }
334 
335 }
336 
337 
338 /*************************************************************************************/
339 /* function to find the closest vertex in 3D for a certain track */
340 /*************************************************************************************/
341 
342 VertexRef
343 PF_PU_AssoMapAlgos::FindClosestZ(const reco::TrackRef trkref, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
344 {
345 
346  double ztrack = trkref->vertex().z();
347 
348  VertexRef foundVertexRef = vtxcollV->at(0);
349 
350  double dzmin = 1e5;
351 
352  //loop over all vertices with a good quality in the vertex collection
353  for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
354 
355  VertexRef vertexref = vtxcollV->at(index_vtx);
356 
357  double nTracks = sqrt(vertexref->tracksSize());
358 
359  double z_distance = fabs(ztrack - vertexref->z());
360 
361  double weightedDistance = z_distance-tWeight*nTracks;
362 
363  if(weightedDistance<dzmin) {
364  dzmin = weightedDistance;
365  foundVertexRef = vertexref;
366  }
367 
368  }
369 
370  return foundVertexRef;
371 }
372 
373 
374 /*************************************************************************************/
375 /* function to find the closest vertex in 3D for a certain track */
376 /*************************************************************************************/
377 
378 VertexRef
379 PF_PU_AssoMapAlgos::FindClosest3D(TransientTrack transtrk, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
380 {
381 
382  VertexRef foundVertexRef = vtxcollV->at(0);
383 
384  double d3min = 1e5;
385 
386  //loop over all vertices with a good quality in the vertex collection
387  for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
388 
389  VertexRef vertexref = vtxcollV->at(index_vtx);
390 
391  double nTracks = sqrt(vertexref->tracksSize());
392 
393  double distance = 1e5;
394  pair<bool,Measurement1D> IpPair = IPTools::absoluteImpactParameter3D(transtrk, *vertexref);
395 
396  if(IpPair.first)
397  distance = IpPair.second.value();
398 
399  double weightedDistance = distance-tWeight*nTracks;
400 
401  if(weightedDistance<d3min) {
402  d3min = weightedDistance;
403  foundVertexRef = vertexref;
404  }
405 
406  }
407 
408  return foundVertexRef;
409 }
410 
411 /****************************************************************************************/
412 /* function to calculate the deltaR between a vector and a vector connecting two points */
413 /****************************************************************************************/
414 
415 double
417 {
418 
419  double bs_x = bsH->x0();
420  double bs_y = bsH->y0();
421  double bs_z = bsH->z0();
422 
423  double connVec_x = vtx_pos.x() - bs_x;
424  double connVec_y = vtx_pos.y() - bs_y;
425  double connVec_z = vtx_pos.z() - bs_z;
426 
427  double connVec_r = sqrt(connVec_x*connVec_x + connVec_y*connVec_y + connVec_z*connVec_z);
428  double connVec_theta = acos(connVec_z*1./connVec_r);
429 
430  double connVec_eta = -1.*log(tan(connVec_theta*1./2.));
431  double connVec_phi = atan2(connVec_y,connVec_x);
432 
433  return deltaR(vtx_mom.eta(),vtx_mom.phi(),connVec_eta,connVec_phi);
434 
435 }
436 
437 /*************************************************************************************/
438 /* function to filter the conversion collection */
439 /*************************************************************************************/
440 
441 auto_ptr<ConversionCollection>
443 {
444  auto_ptr<ConversionCollection> cleanedConvColl(new ConversionCollection() );
445 
446  for (unsigned int convcoll_idx=0; convcoll_idx<convCollH->size(); convcoll_idx++){
447 
448  ConversionRef convref(convCollH,convcoll_idx);
449 
450  if(!cleanedColl){
451  cleanedConvColl->push_back(*convref);
452  continue;
453  }
454 
455  if( (convref->nTracks()==2) &&
456  (fabs(convref->pairInvariantMass())<=0.1) ){
457 
458  cleanedConvColl->push_back(*convref);
459 
460  }
461 
462  }
463 
464  return cleanedConvColl;
465 
466 }
467 
468 
469 /*************************************************************************************/
470 /* function to find out if the track comes from a gamma conversion */
471 /*************************************************************************************/
472 
473 bool
475 {
476 
477  for(unsigned int convcoll_ite=0; convcoll_ite<cleanedConvColl.size(); convcoll_ite++){
478 
479  if(ConversionTools::matchesConversion(trackref,cleanedConvColl.at(convcoll_ite))){
480 
481  *gamma = cleanedConvColl.at(convcoll_ite);
482  return true;
483 
484  }
485 
486  }
487 
488  return false;
489 }
490 
491 
492 /********************************************************************************/
493 /* function to find the closest vertex for a track from a conversion */
494 /********************************************************************************/
495 
496 VertexRef
497 PF_PU_AssoMapAlgos::FindConversionVertex(const reco::TrackRef trackref, const reco::Conversion& gamma, ESHandle<MagneticField> bfH, const EventSetup& iSetup, edm::Handle<reco::BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
498 {
499 
500  math::XYZPoint conv_pos = gamma.conversionVertex().position();
501 
502  math::XYZVector conv_mom(gamma.refittedPair4Momentum().x(),
503  gamma.refittedPair4Momentum().y(),
504  gamma.refittedPair4Momentum().z());
505 
506  Track photon(trackref->chi2(), trackref->ndof(), conv_pos, conv_mom, 0, trackref->covariance());
507 
508  TransientTrack transpho(photon, &(*bfH) );
509  transpho.setBeamSpot(*bsH);
510  transpho.setES(iSetup);
511 
512  return FindClosest3D(transpho, vtxcollV, tWeight);
513 
514 }
515 
516 /*************************************************************************************/
517 /* function to filter the Kshort collection */
518 /*************************************************************************************/
519 
520 auto_ptr<VertexCompositeCandidateCollection>
522 {
523 
524  auto_ptr<VertexCompositeCandidateCollection> cleanedKaonColl(new VertexCompositeCandidateCollection() );
525 
526  for (unsigned int kscoll_idx=0; kscoll_idx<KshortsH->size(); kscoll_idx++){
527 
528  VertexCompositeCandidateRef ksref(KshortsH,kscoll_idx);
529 
530  if(!cleanedColl){
531  cleanedKaonColl->push_back(*ksref);
532  continue;
533  }
534 
535  VertexDistance3D distanceComputer;
536 
537  GlobalPoint dec_pos = RecoVertex::convertPos(ksref->vertex());
538 
539  GlobalError decayVertexError = GlobalError(ksref->vertexCovariance(0,0), ksref->vertexCovariance(0,1), ksref->vertexCovariance(1,1), ksref->vertexCovariance(0,2), ksref->vertexCovariance(1,2), ksref->vertexCovariance(2,2));
540 
541  math::XYZVector dec_mom(ksref->momentum().x(),
542  ksref->momentum().y(),
543  ksref->momentum().z());
544 
545  GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
546  GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
547 
548  double kaon_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(dec_pos, decayVertexError))).significance();
549 
550  if ((ksref->vertex().rho()>=3.) &&
551  (ksref->vertexNormalizedChi2()<=3.) &&
552  (fabs(ksref->mass() - kMass)<=0.01) &&
553  (kaon_significance>15.) &&
554  (PF_PU_AssoMapAlgos::dR(ksref->vertex(),dec_mom,bsH)<=0.3) ){
555 
556  cleanedKaonColl->push_back(*ksref);
557 
558  }
559 
560  }
561 
562  return cleanedKaonColl;
563 }
564 
565 /*************************************************************************************/
566 /* function to filter the Lambda collection */
567 /*************************************************************************************/
568 
569 auto_ptr<VertexCompositeCandidateCollection>
571 {
572 
573  auto_ptr<VertexCompositeCandidateCollection> cleanedLambdaColl(new VertexCompositeCandidateCollection() );
574 
575  for (unsigned int lambdacoll_idx=0; lambdacoll_idx<LambdasH->size(); lambdacoll_idx++){
576 
577  VertexCompositeCandidateRef lambdaref(LambdasH,lambdacoll_idx);
578 
579  if(!cleanedColl){
580  cleanedLambdaColl->push_back(*lambdaref);
581  continue;
582  }
583 
584  VertexDistance3D distanceComputer;
585 
586  GlobalPoint dec_pos = RecoVertex::convertPos(lambdaref->vertex());
587 
588  GlobalError decayVertexError = GlobalError(lambdaref->vertexCovariance(0,0), lambdaref->vertexCovariance(0,1), lambdaref->vertexCovariance(1,1), lambdaref->vertexCovariance(0,2), lambdaref->vertexCovariance(1,2), lambdaref->vertexCovariance(2,2));
589 
590  math::XYZVector dec_mom(lambdaref->momentum().x(),
591  lambdaref->momentum().y(),
592  lambdaref->momentum().z());
593 
594  GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
595  GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
596 
597  double lambda_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(dec_pos, decayVertexError))).significance();
598 
599  if ((lambdaref->vertex().rho()>=3.) &&
600  (lambdaref->vertexNormalizedChi2()<=3.) &&
601  (fabs(lambdaref->mass() - lamMass)<=0.005) &&
602  (lambda_significance>15.) &&
603  (PF_PU_AssoMapAlgos::dR(lambdaref->vertex(),dec_mom,bsH)<=0.3) ){
604 
605  cleanedLambdaColl->push_back(*lambdaref);
606 
607  }
608 
609  }
610 
611  return cleanedLambdaColl;
612 }
613 
614 /*************************************************************************************/
615 /* function to find out if the track comes from a V0 decay */
616 /*************************************************************************************/
617 
618 bool
620 {
621 
622  //the part for the reassociation of particles from Kshort decays
623  for(VertexCompositeCandidateCollection::const_iterator iKS=cleanedKshort.begin(); iKS!=cleanedKshort.end(); iKS++){
624 
625  const RecoChargedCandidate *dauCand1 = dynamic_cast<const RecoChargedCandidate*>(iKS->daughter(0));
626  TrackRef dauTk1 = dauCand1->track();
627  const RecoChargedCandidate *dauCand2 = dynamic_cast<const RecoChargedCandidate*>(iKS->daughter(1));
628  TrackRef dauTk2 = dauCand2->track();
629 
630  if((trackref==dauTk1) || (trackref==dauTk2)){
631 
632  *V0 = *iKS;
633  return true;
634 
635  }
636 
637  }
638 
639  //the part for the reassociation of particles from Lambda decays
640  for(VertexCompositeCandidateCollection::const_iterator iLambda=cleanedLambda.begin(); iLambda!=cleanedLambda.end(); iLambda++){
641 
642  const RecoChargedCandidate *dauCand1 = dynamic_cast<const RecoChargedCandidate*>(iLambda->daughter(0));
643  TrackRef dauTk1 = dauCand1->track();
644  const RecoChargedCandidate *dauCand2 = dynamic_cast<const RecoChargedCandidate*>(iLambda->daughter(1));
645  TrackRef dauTk2 = dauCand2->track();
646 
647  if((trackref==dauTk1) || (trackref==dauTk2)){
648 
649  *V0 = *iLambda;
650  return true;
651 
652  }
653 
654  }
655 
656  return false;
657 }
658 
659 
660 /*************************************************************************************/
661 /* function to find the closest vertex in z for a track from a V0 */
662 /*************************************************************************************/
663 
664 VertexRef
665 PF_PU_AssoMapAlgos::FindV0Vertex(const TrackRef trackref, const VertexCompositeCandidate& V0_vtx, ESHandle<MagneticField> bFieldH, const EventSetup& iSetup, Handle<BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
666 {
667 
668  math::XYZPoint dec_pos = V0_vtx.vertex();
669 
670  math::XYZVector dec_mom(V0_vtx.momentum().x(),
671  V0_vtx.momentum().y(),
672  V0_vtx.momentum().z());
673 
674  Track V0(trackref->chi2(), trackref->ndof(), dec_pos, dec_mom, 0, trackref->covariance());
675 
676  TransientTrack transV0(V0, &(*bFieldH) );
677  transV0.setBeamSpot(*bsH);
678  transV0.setES(iSetup);
679 
680  return FindClosest3D(transV0, vtxcollV, tWeight);
681 
682 }
683 
684 
685 /*************************************************************************************/
686 /* function to filter the nuclear interaction collection */
687 /*************************************************************************************/
688 
689 auto_ptr<PFDisplacedVertexCollection>
691 {
692 
693  auto_ptr<PFDisplacedVertexCollection> cleanedNIColl(new PFDisplacedVertexCollection() );
694 
695  for (PFDisplacedVertexCollection::const_iterator niref=NuclIntH->begin(); niref!=NuclIntH->end(); niref++){
696 
697 
698  if( (niref->isFake()) || !(niref->isNucl()) ) continue;
699 
700  if(!cleanedColl){
701  cleanedNIColl->push_back(*niref);
702  continue;
703  }
704 
705  VertexDistance3D distanceComputer;
706 
707  GlobalPoint ni_pos = RecoVertex::convertPos(niref->position());
708  GlobalError interactionVertexError = RecoVertex::convertError(niref->error());
709 
710  math::XYZVector ni_mom(niref->primaryMomentum().x(),
711  niref->primaryMomentum().y(),
712  niref->primaryMomentum().z());
713 
714  GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
715  GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
716 
717  double nuclint_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(ni_pos, interactionVertexError))).significance();
718 
719  if ((niref->position().rho()>=3.) &&
720  (nuclint_significance>15.) &&
721  (PF_PU_AssoMapAlgos::dR(niref->position(),ni_mom,bsH)<=0.3) ){
722 
723  cleanedNIColl->push_back(*niref);
724 
725  }
726 
727  }
728 
729  return cleanedNIColl;
730 }
731 
732 /*************************************************************************************/
733 /* function to find out if the track comes from a nuclear interaction */
734 /*************************************************************************************/
735 
736 bool
738 {
739 
740  //the part for the reassociation of particles from nuclear interactions
741  for(PFDisplacedVertexCollection::const_iterator iDisplV=cleanedNI.begin(); iDisplV!=cleanedNI.end(); iDisplV++){
742 
743  if(iDisplV->trackWeight(trackref)>1.e-5){
744 
745  *displVtx = *iDisplV;
746  return true;
747 
748  }
749 
750  }
751 
752  return false;
753 }
754 
755 
756 /*************************************************************************************/
757 /* function to find the closest vertex in z for a track from a nuclear interaction */
758 /*************************************************************************************/
759 
760 VertexRef
761 PF_PU_AssoMapAlgos::FindNIVertex(const TrackRef trackref, const PFDisplacedVertex& displVtx, ESHandle<MagneticField> bFieldH, const EventSetup& iSetup, Handle<BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
762 {
763 
765 
766  if((displVtx.isTherePrimaryTracks()) || (displVtx.isThereMergedTracks())){
767 
768  for(TrackCollection::const_iterator trkcoll_ite=refittedTracks.begin(); trkcoll_ite!=refittedTracks.end(); trkcoll_ite++){
769 
770  const TrackBaseRef retrackbaseref = displVtx.originalTrack(*trkcoll_ite);
771 
772  if(displVtx.isIncomingTrack(retrackbaseref)){
773 
774  VertexRef VOAssociation = TrackWeightAssociation(retrackbaseref, vtxcollV);
775 
776  if( VOAssociation->trackWeight(retrackbaseref) >= 1.e-5 ){
777  return VOAssociation;
778  }
779 
780  TransientTrack transIncom(*retrackbaseref, &(*bFieldH) );
781  transIncom.setBeamSpot(*bsH);
782  transIncom.setES(iSetup);
783 
784  return FindClosest3D(transIncom, vtxcollV, tWeight);
785 
786  }
787 
788  }
789 
790  }
791 
792  math::XYZPoint ni_pos = displVtx.position();
793 
794  math::XYZVector ni_mom(displVtx.primaryMomentum().x(),
795  displVtx.primaryMomentum().y(),
796  displVtx.primaryMomentum().z());
797 
798  Track incom(trackref->chi2(), trackref->ndof(), ni_pos, ni_mom, 0, trackref->covariance());
799 
800  TransientTrack transIncom(incom, &(*bFieldH) );
801  transIncom.setBeamSpot(*bsH);
802  transIncom.setES(iSetup);
803 
804  return FindClosest3D(transIncom, vtxcollV, tWeight);
805 
806 }
807 
808 /*************************************************************************************/
809 /* function to find the vertex with the highest TrackWeight for a certain track */
810 /*************************************************************************************/
811 
812 VertexRef
813 PF_PU_AssoMapAlgos::TrackWeightAssociation(const TrackBaseRef& trackbaseRef, std::vector<reco::VertexRef>* vtxcollV)
814 {
815 
816  VertexRef bestvertexref = vtxcollV->at(0);
817  float bestweight = 0.;
818 
819  //loop over all vertices in the vertex collection
820  for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
821 
822  VertexRef vertexref = vtxcollV->at(index_vtx);
823 
824  //get the most probable vertex for the track
825  float weight = vertexref->trackWeight(trackbaseRef);
826  if(weight>bestweight){
827  bestweight = weight;
828  bestvertexref = vertexref;
829  }
830 
831  }
832 
833  return bestvertexref;
834 
835 }
836 
837 /*************************************************************************************/
838 /* find an association for a certain track */
839 /*************************************************************************************/
840 
842 PF_PU_AssoMapAlgos::FindAssociation(const reco::TrackRef& trackref, std::vector<reco::VertexRef>* vtxColl, edm::ESHandle<MagneticField> bfH, const edm::EventSetup& iSetup, edm::Handle<reco::BeamSpot> bsH, int assocNum)
843 {
844 
845  const TrackBaseRef& trackbaseRef = TrackBaseRef(trackref);
846 
847  VertexRef foundVertex;
848 
849  //if it is not the first try of an association jump to the final association
850  //to avoid multiple (secondary) associations and/or unphysical (primary and secondary) associations
851  if ( assocNum>0 ) goto finalStep;
852 
853  // Step 1: First round of association:
854  // Find the vertex with the highest track-to-vertex association weight
855  foundVertex = TrackWeightAssociation(trackbaseRef, vtxColl);
856 
857  if ( foundVertex->trackWeight(trackbaseRef) >= 1.e-5 ){
858  return make_pair( foundVertex, 0. );
859  }
860 
861  // Step 2: Reassociation
862  // Second round of association:
863  // In case no vertex with track-to-vertex association weight > 1.e-5 is found,
864  // check the track originates from a neutral hadron decay, photon conversion or nuclear interaction
865 
866  if ( input_doReassociation_ ) {
867 
868  // Test if the track comes from a photon conversion:
869  // If so, try to find the vertex of the mother particle
870  Conversion gamma;
871  if ( ComesFromConversion(trackref, *cleanedConvCollP, &gamma) ){
872  foundVertex = FindConversionVertex(trackref, gamma, bfH, iSetup, bsH, vtxColl, input_nTrack_);
873  return make_pair( foundVertex, 1. );
874  }
875 
876  // Test if the track comes from a Kshort or Lambda decay:
877  // If so, reassociate the track to the vertex of the V0
879  if ( ComesFromV0Decay(trackref, *cleanedKshortCollP, *cleanedLambdaCollP, &V0) ) {
880  foundVertex = FindV0Vertex(trackref, V0, bfH, iSetup, bsH, vtxColl, input_nTrack_);
881  return make_pair( foundVertex, 1. );
882  }
883 
884  if ( !missingColls ) {
885 
886  // Test if the track comes from a nuclear interaction:
887  // If so, reassociate the track to the vertex of the incoming particle
888  PFDisplacedVertex displVtx;
889  if ( ComesFromNI(trackref, *cleanedNICollP, &displVtx) ){
890  foundVertex = FindNIVertex(trackref, displVtx, bfH, iSetup, bsH, vtxColl, input_nTrack_);
891  return make_pair( foundVertex, 1. );
892  }
893 
894  }
895 
896  }
897 
898  // Step 3: Final association
899  // If no vertex is found with track-to-vertex association weight > 1.e-5
900  // and no reassociation was done do the final association
901  // look for the closest vertex in 3D or in z/longitudinal distance
902  // or associate the track always to the first vertex (default)
903 
904  finalStep:
905 
906  switch (input_FinalAssociation_) {
907 
908  case 1:{
909 
910  // closest in z
911  foundVertex = FindClosestZ(trackref,vtxColl,input_nTrack_);
912  break;
913 
914 
915  }
916 
917  case 2:{
918 
919  // closest in 3D
920  TransientTrack transtrk(trackref, &(*bfH) );
921  transtrk.setBeamSpot(*bsH);
922  transtrk.setES(iSetup);
923 
924  foundVertex = FindClosest3D(transtrk,vtxColl,input_nTrack_);
925  break;
926 
927  }
928 
929  default:{
930 
931  // allways first vertex
932  foundVertex = vtxColl->at(0);
933  break;
934 
935  }
936 
937  }
938 
939  return make_pair( foundVertex, 2. );
940 
941 }
942 
943 /*************************************************************************************/
944 /* get the quality for a certain association */
945 /*************************************************************************************/
946 
947 int
949 {
950 
951  int quality = 0;
952 
953  switch (step) {
954 
955  case 0:{
956 
957  //TrackWeight association
958  if ( distance <= tw_90 ) {
959  quality = 5;
960  } else {
961  if ( distance <= tw_70 ) {
962  quality = 4;
963  } else {
964  if ( distance <= tw_50 ) {
965  quality = 3;
966  } else {
967  quality = 2;
968  }
969  }
970  }
971  break;
972 
973  }
974 
975  case 1:{
976 
977  //Secondary association
978  if ( distance <= sec_70 ) {
979  quality = 4;
980  } else {
981  if ( distance <= sec_50 ) {
982  quality = 3;
983  } else {
984  quality = 2;
985  }
986  }
987  break;
988 
989  }
990 
991  case 2:{
992 
993  //Final association
994  if ( assoc_ite == 1 ) {
995  quality = 1;
996  } else {
997  if ( assoc_ite >= 2 ) {
998  quality = 0;
999  } else {
1000  if ( distance <= fin_70 ) {
1001  quality = 4;
1002  } else {
1003  if ( distance <= fin_50 ) {
1004  quality = 3;
1005  } else {
1006  quality = 2;
1007  }
1008  }
1009  }
1010  }
1011  break;
1012 
1013  }
1014 
1015  default:{
1016 
1017  quality = -1;
1018  break;
1019  }
1020 
1021  }
1022 
1023  return quality;
1024 
1025 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > LambdaCollectionToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
virtual Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const
VertexStepPair FindAssociation(const reco::TrackRef &, std::vector< reco::VertexRef > *, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, int)
reco::Vertex::Point convertPos(const GlobalPoint &p)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< PFDisplacedVertex > PFDisplacedVertexCollection
collection of PFDisplacedVertex objects
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
std::vector< TrackQualityPair > TrackQualityPairVector
edm::ESHandle< MagneticField > bFieldH
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot_
const unsigned int nTracks(const reco::Vertex &sv)
const double lamMass
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
static reco::VertexRef FindNIVertex(const reco::TrackRef, const reco::PFDisplacedVertex &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
const bool isTherePrimaryTracks() const
--—— Provide useful information --—— ///
const_iterator end() const
last iterator over the map (read only)
void setBeamSpot(const reco::BeamSpot &beamSpot)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
const double kMass
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:71
static bool ComesFromNI(const reco::TrackRef, const reco::PFDisplacedVertexCollection &, reco::PFDisplacedVertex *)
static reco::VertexRef FindConversionVertex(const reco::TrackRef, const reco::Conversion &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
virtual const Point & vertex() const
vertex position (overwritten by PF...)
int DefineQuality(int, int, double)
static std::auto_ptr< reco::VertexCompositeCandidateCollection > GetCleanedKshort(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
void EraseVertex(std::vector< reco::VertexRef > *, reco::VertexRef)
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
static reco::VertexRef FindClosest3D(reco::TransientTrack, std::vector< reco::VertexRef > *, double tWeight=0.)
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
const double tw_70
const double fin_70
const bool isIncomingTrack(const reco::TrackBaseRef &originalTrack) const
Is primary or merged track.
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::auto_ptr< TrackToVertexAssMap > CreateTrackToVertexMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
key_type key() const
Accessor for product key.
Definition: Ref.h:264
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:142
const Point & position() const
position
Definition: Vertex.h:99
std::pair< reco::VertexRef, int > VertexStepPair
static std::auto_ptr< reco::VertexCompositeCandidateCollection > GetCleanedLambda(edm::Handle< reco::VertexCompositeCandidateCollection >, edm::Handle< reco::BeamSpot >, bool)
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
static std::auto_ptr< reco::PFDisplacedVertexCollection > GetCleanedNI(edm::Handle< reco::PFDisplacedVertexCollection >, edm::Handle< reco::BeamSpot >, bool)
U second(std::pair< T, U > const &p)
static reco::VertexRef FindV0Vertex(const reco::TrackRef, const reco::VertexCompositeCandidate &, edm::ESHandle< MagneticField >, const edm::EventSetup &, edm::Handle< reco::BeamSpot >, std::vector< reco::VertexRef > *, double)
GlobalErrorBase< double, ErrorMatrixTag > GlobalError
Definition: GlobalError.h:12
std::vector< VertexPtsumPair > VertexPtsumVector
int iEvent
Definition: GenABIO.cc:230
static reco::VertexRef FindClosestZ(const reco::TrackRef, std::vector< reco::VertexRef > *, double tWeight=0.)
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:32
T sqrt(T t)
Definition: SSEVec.h:18
PF_PU_AssoMapAlgos(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
edm::Handle< reco::PFDisplacedVertexCollection > displVertexCollH
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
virtual void GetInputCollections(edm::Event &, const edm::EventSetup &)
std::auto_ptr< reco::VertexCompositeCandidateCollection > cleanedKshortCollP
static reco::VertexRef TrackWeightAssociation(const reco::TrackBaseRef &, std::vector< reco::VertexRef > *)
static bool ComesFromV0Decay(const reco::TrackRef, const reco::VertexCompositeCandidateCollection &, const reco::VertexCompositeCandidateCollection &, reco::VertexCompositeCandidate *)
edm::AssociationMap< edm::OneToManyWithQuality< reco::VertexCollection, reco::TrackCollection, int > > TrackToVertexAssMap
std::auto_ptr< reco::VertexCompositeCandidateCollection > cleanedLambdaCollP
static double dR(const math::XYZPoint &, const math::XYZVector &, edm::Handle< reco::BeamSpot >)
edm::EDGetTokenT< reco::PFDisplacedVertexCollection > NIVertexCollectionToken_
void setES(const edm::EventSetup &es)
edm::Handle< reco::VertexCollection > vtxcollH
std::auto_ptr< TrackToVertexAssMap > SortAssociationMap(TrackToVertexAssMap *)
std::auto_ptr< reco::PFDisplacedVertexCollection > cleanedNICollP
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
const double sec_50
Block of elements.
virtual Vector momentum() const final
spatial momentum vector
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > KshortCollectionToken_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const math::XYZTLorentzVector primaryMomentum(std::string massHypo="PI", bool useRefitted=true, double mass=0.0) const
Momentum of primary or merged track calculated with a mass hypothesis.
const T & get() const
Definition: EventSetup.h:56
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
edm::EDGetTokenT< reco::ConversionCollection > ConversionsCollectionToken_
edm::Handle< reco::VertexCompositeCandidateCollection > vertCompCandCollLambdaH
edm::EDGetTokenT< reco::VertexCollection > token_VertexCollection_
virtual reco::TrackRef track() const
reference to a track
std::vector< reco::VertexRef > * CreateVertexVector(edm::Handle< reco::VertexCollection >)
const double tw_50
const double sec_70
static std::auto_ptr< reco::ConversionCollection > GetCleanedConversions(edm::Handle< reco::ConversionCollection >, edm::Handle< reco::BeamSpot >, bool)
static bool ComesFromConversion(const reco::TrackRef, const reco::ConversionCollection &, reco::Conversion *)
edm::Handle< reco::BeamSpot > beamspotH
const_iterator begin() const
first iterator over the map (read only)
edm::Handle< reco::ConversionCollection > convCollH
int weight
Definition: histoStyle.py:50
const bool isThereMergedTracks() const
If a merged track was identified.
std::auto_ptr< VertexToTrackAssMap > CreateVertexToTrackMap(edm::Handle< reco::TrackCollection >, const edm::EventSetup &)
edm::AssociationMap< edm::OneToManyWithQuality< reco::TrackCollection, reco::VertexCollection, int > > VertexToTrackAssMap
const double tw_90
edm::Handle< reco::VertexCompositeCandidateCollection > vertCompCandCollKshortH
std::auto_ptr< reco::ConversionCollection > cleanedConvCollP
const double fin_50