Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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 TrackAndTrackLinker : public BlockElementLinkerBase {
0008 public:
0009   TrackAndTrackLinker(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   bool useKDTree_, debug_;
0020 };
0021 
0022 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, TrackAndTrackLinker, "TrackAndTrackLinker");
0023 
0024 bool TrackAndTrackLinker::linkPrefilter(const reco::PFBlockElement* e1, const reco::PFBlockElement* e2) const {
0025   return (e1->isLinkedToDisplacedVertex() || e2->isLinkedToDisplacedVertex());
0026 }
0027 
0028 double TrackAndTrackLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0029   constexpr reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
0030   constexpr reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
0031   double dist = -1.0;
0032 
0033   const reco::PFDisplacedTrackerVertexRef& ni1_TO_DISP = elem1->displacedVertexRef(T_TO_DISP);
0034   const reco::PFDisplacedTrackerVertexRef& ni2_TO_DISP = elem2->displacedVertexRef(T_TO_DISP);
0035   const reco::PFDisplacedTrackerVertexRef& ni1_FROM_DISP = elem1->displacedVertexRef(T_FROM_DISP);
0036   const reco::PFDisplacedTrackerVertexRef& ni2_FROM_DISP = elem2->displacedVertexRef(T_FROM_DISP);
0037 
0038   if (ni1_TO_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
0039     if (ni1_TO_DISP == ni2_FROM_DISP) {
0040       dist = 1.0;
0041     }
0042 
0043   if (ni1_FROM_DISP.isNonnull() && ni2_TO_DISP.isNonnull())
0044     if (ni1_FROM_DISP == ni2_TO_DISP) {
0045       dist = 1.0;
0046     }
0047 
0048   if (ni1_FROM_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
0049     if (ni1_FROM_DISP == ni2_FROM_DISP) {
0050       dist = 1.0;
0051     }
0052 
0053   if (elem1->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) &&
0054       elem2->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) {
0055     for (const auto& conv1 : elem1->convRefs()) {
0056       for (const auto& conv2 : elem2->convRefs()) {
0057         if (conv1.isNonnull() && conv2.isNonnull() && conv1 == conv2) {
0058           dist = 1.0;
0059           break;
0060         }
0061       }
0062     }
0063   }
0064 
0065   if (elem1->trackType(reco::PFBlockElement::T_FROM_V0) && elem2->trackType(reco::PFBlockElement::T_FROM_V0)) {
0066     if (elem1->V0Ref().isNonnull() && elem2->V0Ref().isNonnull()) {
0067       if (elem1->V0Ref() == elem2->V0Ref()) {
0068         dist = 1.0;
0069       }
0070     }
0071   }
0072 
0073   return dist;
0074 }