CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackAndGSFLinker.cc
Go to the documentation of this file.
6 
8 public:
10  : BlockElementLinkerBase(conf),
11  useKDTree_(conf.getParameter<bool>("useKDTree")),
12  useConvertedBrems_(conf.getParameter<bool>("useConvertedBrems")),
13  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
14 
15  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
16 
17 private:
19 };
20 
22 
23 double TrackAndGSFLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
25  double dist = -1.0;
26  const reco::PFBlockElementGsfTrack* gsfelem(nullptr);
27  const reco::PFBlockElementTrack* tkelem(nullptr);
28  if (elem1->type() < elem2->type()) {
29  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
30  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem2);
31  } else {
32  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
33  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem1);
34  }
35 
36  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
37  const reco::GsfPFRecTrackRef& gsfref = gsfelem->GsftrackRefPF();
38  const reco::TrackRef& kftrackref = trackref->trackRef();
39  const reco::TrackBaseRef kftrackrefbase(kftrackref);
40  const reco::PFRecTrackRef& refkf = gsfref->kfPFRecTrackRef();
41  if (refkf.isNonnull()) {
42  const reco::TrackRef& gsftrackref = refkf->trackRef();
43  if (gsftrackref.isNonnull() && kftrackref.isNonnull() && kftrackref == gsftrackref) {
44  dist = 0.001;
45  }
46  }
47 
48  //override for converted brems
49  if (useConvertedBrems_) {
50  if (tkelem->isLinkedToDisplacedVertex()) {
51  const std::vector<reco::PFRecTrackRef>& convbrems = gsfref->convBremPFRecTrackRef();
52  for (const auto& convbrem : convbrems) {
53  if (tkelem->trackType(T_FROM_GAMMACONV) && kftrackref == convbrem->trackRef()) {
54  dist = 0.001;
55  } else { // check the base ref as well (for dedicated conversions?)
56  const reco::TrackBaseRef convbrembase(convbrem->trackRef());
57  if (convbrembase == kftrackrefbase) {
58  dist = 0.001;
59  }
60  }
61  }
62  }
63  }
64  return dist;
65 }
Abstract base class for a PFBlock element (track, cluster...)
TrackAndGSFLinker(const edm::ParameterSet &conf)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
Type type() const
const GsfPFRecTrackRef & GsftrackRefPF() const
bool trackType(TrackType trType) const override
bool isLinkedToDisplacedVertex() const override
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
#define DEFINE_EDM_PLUGIN(factory, type, name)
const PFRecTrackRef & trackRefPF() const override