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 PreshowerAndECALLinker : public BlockElementLinkerBase {
0007 public:
0008 PreshowerAndECALLinker(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, PreshowerAndECALLinker, "PreshowerAndECALLinker");
0022
0023 bool PreshowerAndECALLinker::linkPrefilter(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0024 bool result = false;
0025
0026 switch (elem1->type()) {
0027 case reco::PFBlockElement::PS1:
0028 case reco::PFBlockElement::PS2:
0029 result = (elem1->isMultilinksValide(elem2->type()) && !elem1->getMultilinks(elem2->type()).empty() &&
0030 elem2->isMultilinksValide(elem1->type()));
0031 break;
0032 case reco::PFBlockElement::ECAL:
0033 result = (elem2->isMultilinksValide(elem1->type()) && !elem2->getMultilinks(elem1->type()).empty() &&
0034 elem1->isMultilinksValide(elem2->type()));
0035 break;
0036 default:
0037 break;
0038 }
0039 return (useKDTree_ ? result : true);
0040 }
0041
0042 double PreshowerAndECALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0043 const reco::PFBlockElementCluster *pselem(nullptr), *ecalelem(nullptr);
0044 double dist(-1.0);
0045 if (elem1->type() < elem2->type()) {
0046 pselem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0047 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0048 } else {
0049 pselem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0050 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0051 }
0052 const reco::PFClusterRef& psref = pselem->clusterRef();
0053 const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
0054 if (psref.isNull() || ecalref.isNull()) {
0055 throw cms::Exception("BadClusterRefs") << "PFBlockElementCluster's refs are null!";
0056 }
0057
0058
0059 if (useKDTree_ && pselem->isMultilinksValide(ecalelem->type())) {
0060 const reco::PFMultilinksType& multilinks = pselem->getMultilinks(ecalelem->type());
0061 const math::XYZPoint& ecalxyzpos = ecalref->position();
0062 const math::XYZPoint& psxyzpos = psref->position();
0063
0064
0065 reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
0066 for (; mlit != multilinks.end(); ++mlit)
0067 if (mlit->clusterRef == ecalref)
0068 break;
0069
0070
0071 if (mlit != multilinks.end()) {
0072 dist = LinkByRecHit::computeDist(
0073 ecalxyzpos.X() / 1000., ecalxyzpos.Y() / 1000., psxyzpos.X() / 1000., psxyzpos.Y() / 1000., false);
0074 }
0075 } else {
0076 dist = LinkByRecHit::testECALAndPSByRecHit(*ecalref, *psref, debug_);
0077 }
0078 return dist;
0079 }