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 track's elem
27  switch (elem1->type()) {
29  result = (elem1->isMultilinksValide() && !elem1->getMultilinks().empty());
30  break;
32  result = (elem2->isMultilinksValide() && !elem2->getMultilinks().empty());
33  default:
34  break;
35  }
36  return (useKDTree_ ? result : true);
37 }
38 
41  const reco::PFBlockElementCluster* ecalelem(nullptr);
42  const reco::PFBlockElementTrack* tkelem(nullptr);
43  double dist(-1.0);
44  if (elem1->type() < elem2->type()) {
45  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
46  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
47  } else {
48  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
49  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
50  }
51  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
52  const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
53  const reco::PFCluster::REPPoint& ecalreppos = clusterref->positionREP();
54  const reco::PFTrajectoryPoint& tkAtECAL = trackref->extrapolatedPoint(ECALShowerMax);
55  const reco::PFCluster::REPPoint& tkreppos = tkAtECAL.positionREP();
56 
57  // Check if the linking has been done using the KDTree algo
58  // Glowinski & Gouzevitch
59  if (useKDTree_ && tkelem->isMultilinksValide()) { //KDTree Algo
60  const reco::PFMultilinksType& multilinks = tkelem->getMultilinks();
61  const double ecalphi = ecalreppos.Phi();
62  const double ecaleta = ecalreppos.Eta();
63 
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->first == ecalphi) && (mlit->second == ecaleta))
68  break;
69 
70  // If the link exist, we fill dist and linktest.
71  if (mlit != multilinks.end()) {
72  dist = LinkByRecHit::computeDist(ecaleta, ecalphi, tkreppos.Eta(), tkreppos.Phi());
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 }
electrons_cff.bool
bool
Definition: electrons_cff.py:372
funct::false
false
Definition: Factorize.h:34
reco::PFMultilinksType
std::vector< std::pair< double, double > > PFMultilinksType
Abstract This class is used by the KDTree Track / Ecal Cluster linker to store all found links.
Definition: PFMultilinksTC.h:13
PFBlockElementCluster.h
TrackAndECALLinker::testLink
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
Definition: TrackAndECALLinker.cc:39
TrackAndECALLinker::useKDTree_
const bool useKDTree_
Definition: TrackAndECALLinker.cc:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
LinkByRecHit::testTrackAndClusterByRecHit
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:16
reco::PFBlockElement::getMultilinks
const PFMultilinksType & getMultilinks() const
Definition: PFBlockElement.h:119
BlockElementLinkerBase
Definition: BlockElementLinkerBase.h:10
edm::Ref< PFRecTrackCollection >
reco::PFTrajectoryPoint::positionREP
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
Definition: PFTrajectoryPoint.h:103
TrackAndECALLinker
Definition: TrackAndECALLinker.cc:7
reco::PFBlockElementTrack::trackRefPF
const PFRecTrackRef & trackRefPF() const override
Definition: PFBlockElementTrack.h:46
reco::PFBlockElement::TRACK
Definition: PFBlockElement.h:32
reco::PFTrajectoryPoint::ECALShowerMax
Definition: PFTrajectoryPoint.h:46
PFCluster.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
edm::ParameterSet
Definition: ParameterSet.h:36
edmplugin::PluginFactory
Definition: PluginFactory.h:34
LinkByRecHit.h
reco::PFBlockElement::ECAL
Definition: PFBlockElement.h:35
TrackAndECALLinker::TrackAndECALLinker
TrackAndECALLinker(const edm::ParameterSet &conf)
Definition: TrackAndECALLinker.cc:9
reco::PFTrajectoryPoint::LayerType
LayerType
Define the different layers where the track can be propagated.
Definition: PFTrajectoryPoint.h:34
reco::PFTrajectoryPoint::isValid
bool isValid() const
is this point valid ?
Definition: PFTrajectoryPoint.h:84
reco::PFBlockElement
Abstract base class for a PFBlock element (track, cluster...)
Definition: PFBlockElement.h:26
reco::PFTrajectoryPoint
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
Definition: PFTrajectoryPoint.h:26
reco::PFBlockElementCluster
Cluster Element.
Definition: PFBlockElementCluster.h:16
reco::PFBlockElementTrack
Track Element.
Definition: PFBlockElementTrack.h:17
BlockElementLinkerBase.h
reco::PFBlockElementCluster::clusterRef
const PFClusterRef & clusterRef() const override
Definition: PFBlockElementCluster.h:29
reco::PFBlockElement::type
Type type() const
Definition: PFBlockElement.h:69
mps_fire.result
result
Definition: mps_fire.py:303
reco::PFCluster::REPPoint
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:52
TrackAndECALLinker::linkPrefilter
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
Definition: TrackAndECALLinker.cc:24
LinkByRecHit::computeDist
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
Definition: LinkByRecHit.cc:517
PFBlockElementTrack.h
TrackAndECALLinker::debug_
const bool debug_
Definition: TrackAndECALLinker.cc:19
reco::PFBlockElement::isMultilinksValide
bool isMultilinksValide() const
Definition: PFBlockElement.h:118