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  // PS-ECAL KDTree multilinks are stored to PS's elem
26  switch (elem1->type()) {
29  result = (elem1->isMultilinksValide(elem2->type()) && !elem1->getMultilinks(elem2->type()).empty() &&
30  elem2->isMultilinksValide(elem1->type()));
31  break;
33  result = (elem2->isMultilinksValide(elem1->type()) && !elem2->getMultilinks(elem1->type()).empty() &&
34  elem1->isMultilinksValide(elem2->type()));
35  break;
36  default:
37  break;
38  }
39  return (useKDTree_ ? result : true);
40 }
41 
43  const reco::PFBlockElementCluster *pselem(nullptr), *ecalelem(nullptr);
44  double dist(-1.0);
45  if (elem1->type() < elem2->type()) {
46  pselem = static_cast<const reco::PFBlockElementCluster*>(elem1);
47  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
48  } else {
49  pselem = static_cast<const reco::PFBlockElementCluster*>(elem2);
50  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
51  }
52  const reco::PFClusterRef& psref = pselem->clusterRef();
53  const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
54  if (psref.isNull() || ecalref.isNull()) {
55  throw cms::Exception("BadClusterRefs") << "PFBlockElementCluster's refs are null!";
56  }
57  // Check if the linking has been done using the KDTree algo
58  // Glowinski & Gouzevitch
59  if (useKDTree_ && pselem->isMultilinksValide(ecalelem->type())) { // KDTree algo
60  const reco::PFMultilinksType& multilinks = pselem->getMultilinks(ecalelem->type());
61  const math::XYZPoint& ecalxyzpos = ecalref->position();
62  const math::XYZPoint& psxyzpos = psref->position();
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->clusterRef == ecalref)
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 }
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
static double testECALAndPSByRecHit(const reco::PFCluster &clusterECAL, const reco::PFCluster &clusterPS, bool debug=false)
const PFMultilinksType & getMultilinks(Type type) const
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
std::vector< PFMultilink > PFMultilinksType
collection of PFSuperCluster objects
bool isNull() const
Checks for null.
Definition: Ref.h:235
bool isMultilinksValide(Type type) const
const PFClusterRef & clusterRef() const override
PreshowerAndECALLinker(const edm::ParameterSet &conf)
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
#define DEFINE_EDM_PLUGIN(factory, type, name)