Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:28

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 ECALAndECALLinker : public BlockElementLinkerBase {
0007 public:
0008   ECALAndECALLinker(const edm::ParameterSet& conf)
0009       : BlockElementLinkerBase(conf),
0010         useKDTree_(conf.getParameter<bool>("useKDTree")),
0011         debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
0012 
0013   bool linkPrefilter(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0014 
0015   double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0016 
0017 private:
0018   bool useKDTree_, debug_;
0019 };
0020 
0021 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, ECALAndECALLinker, "ECALAndECALLinker");
0022 
0023 bool ECALAndECALLinker::linkPrefilter(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0024   const reco::PFBlockElementCluster* ecal1 = static_cast<const reco::PFBlockElementCluster*>(elem1);
0025   const reco::PFBlockElementCluster* ecal2 = static_cast<const reco::PFBlockElementCluster*>(elem2);
0026   return (ecal1->superClusterRef().isNonnull() && ecal2->superClusterRef().isNonnull());
0027 }
0028 
0029 double ECALAndECALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0030   double dist = -1.0;
0031 
0032   const reco::PFBlockElementCluster* ecal1 = static_cast<const reco::PFBlockElementCluster*>(elem1);
0033   const reco::PFBlockElementCluster* ecal2 = static_cast<const reco::PFBlockElementCluster*>(elem2);
0034 
0035   const reco::SuperClusterRef& sc1 = ecal1->superClusterRef();
0036   const reco::SuperClusterRef& sc2 = ecal2->superClusterRef();
0037 
0038   const reco::PFClusterRef& clus1 = ecal1->clusterRef();
0039   const reco::PFClusterRef& clus2 = ecal2->clusterRef();
0040 
0041   if (sc1.isNonnull() && sc2.isNonnull() && sc1 == sc2) {
0042     dist = LinkByRecHit::computeDist(
0043         clus1->positionREP().Eta(), clus1->positionREP().Phi(), clus2->positionREP().Eta(), clus2->positionREP().Phi());
0044   }
0045 
0046   return dist;
0047 }