CMS 3D CMS Logo

SCAndECALLinker.cc
Go to the documentation of this file.
7 
9 public:
12  _useKDTree(conf.getParameter<bool>("useKDTree")),
13  _debug(conf.getUntrackedParameter<bool>("debug",false)),
14  _superClusterMatchByRef(conf.getParameter<bool>("SuperClusterMatchByRef")){}
15 
16  double testLink
17  ( const reco::PFBlockElement*,
18  const reco::PFBlockElement* ) const override;
19 
20 private:
22 };
23 
26  "SCAndECALLinker");
27 
29  ( const reco::PFBlockElement* elem1,
30  const reco::PFBlockElement* elem2) const {
31  double dist = -1.0;
32  const reco::PFBlockElementCluster* ecalelem(NULL);
33  const reco::PFBlockElementSuperCluster* scelem(NULL);
34  if( elem1->type() < elem2->type() ) {
35  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
36  scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem2);
37  } else {
38  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
39  scelem = static_cast<const reco::PFBlockElementSuperCluster*>(elem1);
40  }
41  const reco::PFClusterRef& clus = ecalelem->clusterRef();
42  const reco::SuperClusterRef& sclus = scelem->superClusterRef();
43  if( sclus.isNull() ) {
44  throw cms::Exception("BadRef")
45  << "SuperClusterRef is invalid!";
46  }
47 
49  if( sclus == ecalelem->superClusterRef() ) dist = 0.001;
50  } else {
51  if( ClusterClusterMapping::overlap(*sclus,*clus) ) {
52  dist = LinkByRecHit::computeDist( sclus->position().eta(),
53  sclus->position().phi(),
54  clus->positionREP().Eta(),
55  clus->positionREP().Phi() );
56  }
57  }
58  return dist;
59 }
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
#define NULL
Definition: scimark2.h:8
const PFClusterRef & clusterRef() const
const SuperClusterRef & superClusterRef() const
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
bool isNull() const
Checks for null.
Definition: Ref.h:250
SCAndECALLinker(const edm::ParameterSet &conf)
#define DEFINE_EDM_PLUGIN(factory, type, name)