CMS 3D CMS Logo

PreshowerAndECALLinker.cc
Go to the documentation of this file.
5 
7 public:
10  useKDTree_(conf.getParameter<bool>("useKDTree")),
11  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
12 
13  bool linkPrefilter(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
14 
15  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
16 
17 private:
19 };
20 
22 
24  bool result = false;
25  switch (elem1->type()) {
28  result = (elem1->isMultilinksValide() && !elem1->getMultilinks().empty());
29  break;
31  result = (elem2->isMultilinksValide() && !elem2->getMultilinks().empty());
32  break;
33  default:
34  break;
35  }
36  return (useKDTree_ ? result : true);
37 }
38 
40  const reco::PFBlockElementCluster *pselem(nullptr), *ecalelem(nullptr);
41  double dist(-1.0);
42  if (elem1->type() < elem2->type()) {
43  pselem = static_cast<const reco::PFBlockElementCluster*>(elem1);
44  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
45  } else {
46  pselem = static_cast<const reco::PFBlockElementCluster*>(elem2);
47  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
48  }
49  const reco::PFClusterRef& psref = pselem->clusterRef();
50  const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
51  if (psref.isNull() || ecalref.isNull()) {
52  throw cms::Exception("BadClusterRefs") << "PFBlockElementCluster's refs are null!";
53  }
54  // Check if the linking has been done using the KDTree algo
55  // Glowinski & Gouzevitch
56  if (useKDTree_ && pselem->isMultilinksValide()) { // KDTree algo
57  const reco::PFMultilinksType& multilinks = pselem->getMultilinks();
58  const reco::PFCluster::REPPoint& ecalreppos = ecalref->positionREP();
59  const math::XYZPoint& ecalxyzpos = ecalref->position();
60  const math::XYZPoint& psxyzpos = psref->position();
61  const double ecalPhi = ecalreppos.Phi();
62  const double ecalEta = ecalreppos.Eta();
63 
64  // Check if the link PS/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()) {
73  ecalxyzpos.X() / 1000., ecalxyzpos.Y() / 1000., psxyzpos.X() / 1000., psxyzpos.Y() / 1000., false);
74  }
75  } else { //Old algorithm
76  dist = LinkByRecHit::testECALAndPSByRecHit(*ecalref, *psref, debug_);
77  }
78  return dist;
79 }
std::vector< std::pair< double, double > > PFMultilinksType
Abstract This class is used by the KDTree Track / Ecal Cluster linker to store all found links...
Abstract base class for a PFBlock element (track, cluster...)
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
const PFClusterRef & clusterRef() const override
static double testECALAndPSByRecHit(const reco::PFCluster &clusterECAL, const reco::PFCluster &clusterPS, bool debug=false)
Type type() const
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
bool isNull() const
Checks for null.
Definition: Ref.h:235
PreshowerAndECALLinker(const edm::ParameterSet &conf)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isMultilinksValide() const
double ecalEta(const math::XYZVector &momentum, const math::XYZPoint &vertex)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:52
#define DEFINE_EDM_PLUGIN(factory, type, name)
double ecalPhi(const MagneticField &magField, const math::XYZVector &momentum, const math::XYZPoint &vertex, const int charge)
const PFMultilinksType & getMultilinks() const