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/PFRecTrack.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
0005 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
0006 
0007 class TrackAndGSFLinker : public BlockElementLinkerBase {
0008 public:
0009   TrackAndGSFLinker(const edm::ParameterSet& conf)
0010       : BlockElementLinkerBase(conf),
0011         useKDTree_(conf.getParameter<bool>("useKDTree")),
0012         useConvertedBrems_(conf.getParameter<bool>("useConvertedBrems")),
0013         debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
0014 
0015   double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
0016 
0017 private:
0018   bool useKDTree_, useConvertedBrems_, debug_;
0019 };
0020 
0021 DEFINE_EDM_PLUGIN(BlockElementLinkerFactory, TrackAndGSFLinker, "TrackAndGSFLinker");
0022 
0023 double TrackAndGSFLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
0024   constexpr reco::PFBlockElement::TrackType T_FROM_GAMMACONV = reco::PFBlockElement::T_FROM_GAMMACONV;
0025   double dist = -1.0;
0026   const reco::PFBlockElementGsfTrack* gsfelem(nullptr);
0027   const reco::PFBlockElementTrack* tkelem(nullptr);
0028   if (elem1->type() < elem2->type()) {
0029     tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
0030     gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem2);
0031   } else {
0032     tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
0033     gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem1);
0034   }
0035 
0036   const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
0037   const reco::GsfPFRecTrackRef& gsfref = gsfelem->GsftrackRefPF();
0038   const reco::TrackRef& kftrackref = trackref->trackRef();
0039   const reco::TrackBaseRef kftrackrefbase(kftrackref);
0040   const reco::PFRecTrackRef& refkf = gsfref->kfPFRecTrackRef();
0041   if (refkf.isNonnull()) {
0042     const reco::TrackRef& gsftrackref = refkf->trackRef();
0043     if (gsftrackref.isNonnull() && kftrackref.isNonnull() && kftrackref == gsftrackref) {
0044       dist = 0.001;
0045     }
0046   }
0047 
0048   //override for converted brems
0049   if (useConvertedBrems_) {
0050     if (tkelem->isLinkedToDisplacedVertex()) {
0051       const std::vector<reco::PFRecTrackRef>& convbrems = gsfref->convBremPFRecTrackRef();
0052       for (const auto& convbrem : convbrems) {
0053         if (tkelem->trackType(T_FROM_GAMMACONV) && kftrackref == convbrem->trackRef()) {
0054           dist = 0.001;
0055         } else {  // check the base ref as well (for dedicated conversions?)
0056           const reco::TrackBaseRef convbrembase(convbrem->trackRef());
0057           if (convbrembase == kftrackrefbase) {
0058             dist = 0.001;
0059           }
0060         }
0061       }
0062     }
0063   }
0064   return dist;
0065 }