CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SCAndECALLinker.cc
Go to the documentation of this file.
7 
9 public:
11  : BlockElementLinkerBase(conf),
12  useKDTree_(conf.getParameter<bool>("useKDTree")),
13  debug_(conf.getUntrackedParameter<bool>("debug", false)),
14  superClusterMatchByRef_(conf.getParameter<bool>("SuperClusterMatchByRef")) {}
15 
16  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
17 
18 private:
20 };
21 
23 
24 double SCAndECALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
25  double dist = -1.0;
26  const reco::PFBlockElementCluster* ecalelem(nullptr);
27  const reco::PFBlockElementSuperCluster* scelem(nullptr);
28  if (elem1->type() < elem2->type()) {
29  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
30  scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem2);
31  } else {
32  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
33  scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem1);
34  }
35  const reco::PFClusterRef& clus = ecalelem->clusterRef();
36  const reco::SuperClusterRef& sclus = scelem->superClusterRef();
37  if (sclus.isNull()) {
38  throw cms::Exception("BadRef") << "SuperClusterRef is invalid!";
39  }
40 
42  if (sclus == ecalelem->superClusterRef())
43  dist = 0.001;
44  } else {
45  if (ClusterClusterMapping::overlap(*sclus, *clus)) {
47  sclus->position().eta(), sclus->position().phi(), clus->positionREP().Eta(), clus->positionREP().Phi());
48  }
49  }
50  return dist;
51 }
const SuperClusterRef & superClusterRef() const
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 bool overlap(const reco::CaloCluster &sc1, const reco::CaloCluster &sc, float minfrac=0.01, bool debug=false)
Type type() const
const SuperClusterRef & superClusterRef() const
bool isNull() const
Checks for null.
Definition: Ref.h:235
const PFClusterRef & clusterRef() const override
SCAndECALLinker(const edm::ParameterSet &conf)
#define DEFINE_EDM_PLUGIN(factory, type, name)
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override