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 "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
0005 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
0006
0007 class TrackAndECALLinker : public BlockElementLinkerBase {
0008 public:
0009 TrackAndECALLinker(const edm::ParameterSet& conf)
0010 : BlockElementLinkerBase(conf),
0011 useKDTree_(conf.getParameter<bool>("useKDTree")),
0012 debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
0013
0014 bool linkPrefilter(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0015
0016 double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0017
0018 private:
0019 const bool useKDTree_, debug_;
0020 };
0021
0022 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, TrackAndECALLinker, "TrackAndECALLinker");
0023
0024 bool TrackAndECALLinker::linkPrefilter(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0025 bool result = false;
0026
0027 switch (elem1->type()) {
0028 case reco::PFBlockElement::TRACK:
0029 result = (elem2->isMultilinksValide(elem1->type()) && !elem2->getMultilinks(elem1->type()).empty() &&
0030 elem1->isMultilinksValide(elem2->type()));
0031 break;
0032 case reco::PFBlockElement::ECAL:
0033 result = (elem1->isMultilinksValide(elem2->type()) && !elem1->getMultilinks(elem2->type()).empty() &&
0034 elem2->isMultilinksValide(elem1->type()));
0035 default:
0036 break;
0037 }
0038 return (useKDTree_ ? result : true);
0039 }
0040
0041 double TrackAndECALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0042 constexpr reco::PFTrajectoryPoint::LayerType ECALShowerMax = reco::PFTrajectoryPoint::ECALShowerMax;
0043 const reco::PFBlockElementCluster* ecalelem(nullptr);
0044 const reco::PFBlockElementTrack* tkelem(nullptr);
0045 double dist(-1.0);
0046 if (elem1->type() < elem2->type()) {
0047 tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
0048 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0049 } else {
0050 tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
0051 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0052 }
0053 const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
0054 const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
0055 const reco::PFCluster::REPPoint& ecalreppos = clusterref->positionREP();
0056 const reco::PFTrajectoryPoint& tkAtECAL = trackref->extrapolatedPoint(ECALShowerMax);
0057
0058
0059
0060 if (useKDTree_ && ecalelem->isMultilinksValide(tkelem->type())) {
0061 const reco::PFMultilinksType& multilinks = ecalelem->getMultilinks(tkelem->type());
0062 const double tracketa = tkAtECAL.positionREP().Eta();
0063 const double trackphi = tkAtECAL.positionREP().Phi();
0064
0065 reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
0066 for (; mlit != multilinks.end(); ++mlit)
0067 if (mlit->trackRef == trackref)
0068 break;
0069
0070
0071 if (mlit != multilinks.end()) {
0072 dist = LinkByRecHit::computeDist(ecalreppos.Eta(), ecalreppos.Phi(), tracketa, trackphi);
0073 }
0074
0075 } else {
0076 if (tkAtECAL.isValid())
0077 dist = LinkByRecHit::testTrackAndClusterByRecHit(*trackref, *clusterref, false, debug_);
0078 }
0079
0080 if (debug_) {
0081 if (dist > 0.) {
0082 std::cout << " Here a link has been established"
0083 << " between a track an Ecal with dist " << dist << std::endl;
0084 } else
0085 std::cout << " No link found " << std::endl;
0086 }
0087 return dist;
0088 }