Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:21:09

0001 #include "RecoParticleFlow/PFProducer/interface/BlockElementLinkerBase.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
0004 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
0005 
0006 class ECALAndHCALLinker : public BlockElementLinkerBase {
0007 public:
0008   ECALAndHCALLinker(const edm::ParameterSet& conf)
0009       : BlockElementLinkerBase(conf),
0010         minAbsEtaEcal_(conf.getParameter<double>("minAbsEtaEcal")),
0011         useKDTree_(conf.getParameter<bool>("useKDTree")),
0012         debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
0013 
0014   double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0015 
0016 private:
0017   double minAbsEtaEcal_;
0018   bool useKDTree_, debug_;
0019 };
0020 
0021 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, ECALAndHCALLinker, "ECALAndHCALLinker");
0022 
0023 double ECALAndHCALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0024   const reco::PFBlockElementCluster *hcalelem(nullptr), *ecalelem(nullptr);
0025   double dist(-1.0);
0026   if (elem1->type() < elem2->type()) {
0027     ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0028     hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0029   } else {
0030     ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0031     hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0032   }
0033   const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
0034   const reco::PFClusterRef& hcalref = hcalelem->clusterRef();
0035   const reco::PFCluster::REPPoint& ecalreppos = ecalref->positionREP();
0036   if (hcalref.isNull() || ecalref.isNull()) {
0037     throw cms::Exception("BadClusterRefs") << "PFBlockElementCluster's refs are null!";
0038   }
0039   dist = (std::abs(ecalreppos.Eta()) > minAbsEtaEcal_
0040               ? LinkByRecHit::computeDist(
0041                     ecalreppos.Eta(), ecalreppos.Phi(), hcalref->positionREP().Eta(), hcalref->positionREP().Phi())
0042               : -1.0);
0043   return (dist < 0.2 ? dist : -1.0);
0044 }