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/PFBlockElementGsfTrack.h"
0005 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
0006
0007 class GSFAndECALLinker : public BlockElementLinkerBase {
0008 public:
0009 GSFAndECALLinker(const edm::ParameterSet& conf)
0010 : BlockElementLinkerBase(conf),
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 bool useKDTree_, debug_;
0018 };
0019
0020 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, GSFAndECALLinker, "GSFAndECALLinker");
0021
0022 double GSFAndECALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0023 constexpr reco::PFTrajectoryPoint::LayerType ECALShowerMax = reco::PFTrajectoryPoint::ECALShowerMax;
0024 const reco::PFBlockElementCluster* ecalelem(nullptr);
0025 const reco::PFBlockElementGsfTrack* gsfelem(nullptr);
0026 double dist(-1.0);
0027 if (elem1->type() < elem2->type()) {
0028 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
0029 gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem2);
0030 } else {
0031 ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
0032 gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem1);
0033 }
0034 const reco::PFRecTrack& track = gsfelem->GsftrackPF();
0035 const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
0036 const reco::PFTrajectoryPoint& tkAtECAL = track.extrapolatedPoint(ECALShowerMax);
0037 if (tkAtECAL.isValid()) {
0038 dist = LinkByRecHit::testTrackAndClusterByRecHit(track, *clusterref, false, debug_);
0039 }
0040 if (debug_) {
0041 if (dist > 0.) {
0042 std::cout << " Here a link has been established"
0043 << " between a GSF track an Ecal with dist " << dist << std::endl;
0044 } else {
0045 if (debug_)
0046 std::cout << " No link found " << std::endl;
0047 }
0048 }
0049
0050 return dist;
0051 }