CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ECALAndHCALLinker.cc
Go to the documentation of this file.
5 
7 public:
10  minAbsEtaEcal_(conf.getParameter<double>("minAbsEtaEcal")),
11  useKDTree_(conf.getParameter<bool>("useKDTree")),
12  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
13 
14  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
15 
16 private:
19 };
20 
22 
23 double ECALAndHCALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
24  const reco::PFBlockElementCluster *hcalelem(nullptr), *ecalelem(nullptr);
25  double dist(-1.0);
26  if (elem1->type() < elem2->type()) {
27  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
28  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
29  } else {
30  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
31  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
32  }
33  const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
34  const reco::PFClusterRef& hcalref = hcalelem->clusterRef();
35  const reco::PFCluster::REPPoint& ecalreppos = ecalref->positionREP();
36  if (hcalref.isNull() || ecalref.isNull()) {
37  throw cms::Exception("BadClusterRefs") << "PFBlockElementCluster's refs are null!";
38  }
39  dist = (std::abs(ecalreppos.Eta()) > minAbsEtaEcal_
41  ecalreppos.Eta(), ecalreppos.Phi(), hcalref->positionREP().Eta(), hcalref->positionREP().Phi())
42  : -1.0);
43  return (dist < 0.2 ? dist : -1.0);
44 }
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
Type type() const
ECALAndHCALLinker(const edm::ParameterSet &conf)
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isNull() const
Checks for null.
Definition: Ref.h:235
const PFClusterRef & clusterRef() const override
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
#define DEFINE_EDM_PLUGIN(factory, type, name)