CMS 3D CMS Logo

TrackAndECALLinker.cc
Go to the documentation of this file.
6 
8 public:
10  : BlockElementLinkerBase(conf),
11  useKDTree_(conf.getParameter<bool>("useKDTree")),
12  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
13 
14  bool linkPrefilter(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
15 
16  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
17 
18 private:
19  const bool useKDTree_, debug_;
20 };
21 
23 
25  bool result = false;
26  // Track-ECAL KDTree multilinks are stored to ecal's elem
27  switch (elem1->type()) {
29  result = (elem2->isMultilinksValide(elem1->type()) && !elem2->getMultilinks(elem1->type()).empty() &&
30  elem1->isMultilinksValide(elem2->type()));
31  break;
33  result = (elem1->isMultilinksValide(elem2->type()) && !elem1->getMultilinks(elem2->type()).empty() &&
34  elem2->isMultilinksValide(elem1->type()));
35  default:
36  break;
37  }
38  return (useKDTree_ ? result : true);
39 }
40 
43  const reco::PFBlockElementCluster* ecalelem(nullptr);
44  const reco::PFBlockElementTrack* tkelem(nullptr);
45  double dist(-1.0);
46  if (elem1->type() < elem2->type()) {
47  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
48  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
49  } else {
50  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
51  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
52  }
53  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
54  const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
55  const reco::PFCluster::REPPoint& ecalreppos = clusterref->positionREP();
56  const reco::PFTrajectoryPoint& tkAtECAL = trackref->extrapolatedPoint(ECALShowerMax);
57 
58  // Check if the linking has been done using the KDTree algo
59  // Glowinski & Gouzevitch
60  if (useKDTree_ && ecalelem->isMultilinksValide(tkelem->type())) { //KDTree Algo
61  const reco::PFMultilinksType& multilinks = ecalelem->getMultilinks(tkelem->type());
62  const double tracketa = tkAtECAL.positionREP().Eta();
63  const double trackphi = tkAtECAL.positionREP().Phi();
64  // Check if the link Track/Ecal exist
65  reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
66  for (; mlit != multilinks.end(); ++mlit)
67  if (mlit->trackRef == trackref)
68  break;
69 
70  // If the link exist, we fill dist and linktest.
71  if (mlit != multilinks.end()) {
72  dist = LinkByRecHit::computeDist(ecalreppos.Eta(), ecalreppos.Phi(), tracketa, trackphi);
73  }
74 
75  } else { // Old algorithm
76  if (tkAtECAL.isValid())
77  dist = LinkByRecHit::testTrackAndClusterByRecHit(*trackref, *clusterref, false, debug_);
78  }
79 
80  if (debug_) {
81  if (dist > 0.) {
82  std::cout << " Here a link has been established"
83  << " between a track an Ecal with dist " << dist << std::endl;
84  } else
85  std::cout << " No link found " << std::endl;
86  }
87  return dist;
88 }
Abstract base class for a PFBlock element (track, cluster...)
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
const PFMultilinksType & getMultilinks(Type type) const
bool isValid() const
is this point valid ?
std::vector< PFMultilink > PFMultilinksType
collection of PFSuperCluster objects
bool isMultilinksValide(Type type) const
const PFClusterRef & clusterRef() const override
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
const PFRecTrackRef & trackRefPF() const override
LayerType
Define the different layers where the track can be propagated.
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:18
TrackAndECALLinker(const edm::ParameterSet &conf)