CMS 3D CMS Logo

TrackAndTrackLinker.cc
Go to the documentation of this file.
6 
8 public:
11  _useKDTree(conf.getParameter<bool>("useKDTree")),
12  _debug(conf.getUntrackedParameter<bool>("debug",false)) {}
13 
15  const reco::PFBlockElement* ) const override;
16 
17  double testLink( const reco::PFBlockElement*,
18  const reco::PFBlockElement* ) const override;
19 
20 private:
22 };
23 
26  "TrackAndTrackLinker");
27 
30  const reco::PFBlockElement* e2 ) const {
31  return ( e1->isLinkedToDisplacedVertex() ||
33 }
34 
37  const reco::PFBlockElement* elem2) const {
42  double dist = -1.0;
43 
44  const reco::PFDisplacedTrackerVertexRef& ni1_TO_DISP =
45  elem1->displacedVertexRef(T_TO_DISP);
46  const reco::PFDisplacedTrackerVertexRef& ni2_TO_DISP =
47  elem2->displacedVertexRef(T_TO_DISP);
48  const reco::PFDisplacedTrackerVertexRef& ni1_FROM_DISP =
49  elem1->displacedVertexRef(T_FROM_DISP);
50  const reco::PFDisplacedTrackerVertexRef& ni2_FROM_DISP =
51  elem2->displacedVertexRef(T_FROM_DISP);
52 
53  if( ni1_TO_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
54  if( ni1_TO_DISP == ni2_FROM_DISP ) { dist = 1.0; }
55 
56  if( ni1_FROM_DISP.isNonnull() && ni2_TO_DISP.isNonnull())
57  if( ni1_FROM_DISP == ni2_TO_DISP ) { dist = 1.0; }
58 
59  if( ni1_FROM_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
60  if( ni1_FROM_DISP == ni2_FROM_DISP ) { dist = 1.0; }
61 
64  for( const auto& conv1 : elem1->convRefs() ) {
65  for( const auto& conv2 : elem2->convRefs() ) {
66  if( conv1.isNonnull() && conv2.isNonnull() &&
67  conv1 == conv2 ) {
68  dist=1.0;
69  break;
70  }
71  }
72  }
73  }
74 
77  if ( elem1->V0Ref().isNonnull() && elem2->V0Ref().isNonnull() ) {
78  if ( elem1->V0Ref() == elem2->V0Ref() ) {
79  dist=1.0;
80  }
81  }
82  }
83 
84  return dist;
85 }
Abstract base class for a PFBlock element (track, cluster...)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
TrackAndTrackLinker(const edm::ParameterSet &conf)
#define constexpr
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
virtual const ConversionRefVector & convRefs() const
Float e1
Definition: deltaR.h:20
virtual const PFDisplacedTrackerVertexRef & displacedVertexRef(TrackType trType) const
virtual bool trackType(TrackType trType) const
Float e2
Definition: deltaR.h:21
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
#define DEFINE_EDM_PLUGIN(factory, type, name)
virtual bool isLinkedToDisplacedVertex() const
virtual const VertexCompositeCandidateRef & V0Ref() const