File indexing completed on 2024-04-06 12:27:50
0001 #include "RecoTauTag/RecoTau/interface/PFTauPrimaryVertexProducerBase.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0004 #include "DataFormats/Math/interface/deltaR.h"
0005
0006
0007 class PFTauMiniAODPrimaryVertexProducer final : public PFTauPrimaryVertexProducerBase {
0008 public:
0009 explicit PFTauMiniAODPrimaryVertexProducer(const edm::ParameterSet &iConfig);
0010 ~PFTauMiniAODPrimaryVertexProducer() override;
0011
0012 void beginEvent(const edm::Event &, const edm::EventSetup &) override;
0013 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0014
0015 protected:
0016 void nonTauTracksInPV(const reco::VertexRef &,
0017 const std::vector<edm::Ptr<reco::TrackBase> > &,
0018 std::vector<const reco::Track *> &) override;
0019
0020 private:
0021 void nonTauTracksInPVFromPackedCands(const size_t &,
0022 const pat::PackedCandidateCollection &,
0023 const std::vector<edm::Ptr<reco::TrackBase> > &,
0024 std::vector<const reco::Track *> &);
0025
0026 edm::EDGetTokenT<pat::PackedCandidateCollection> packedCandsToken_, lostCandsToken_;
0027 edm::Handle<pat::PackedCandidateCollection> packedCands_, lostCands_;
0028 };
0029
0030 PFTauMiniAODPrimaryVertexProducer::PFTauMiniAODPrimaryVertexProducer(const edm::ParameterSet &iConfig)
0031 : PFTauPrimaryVertexProducerBase::PFTauPrimaryVertexProducerBase(iConfig),
0032 packedCandsToken_(
0033 consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("packedCandidatesTag"))),
0034 lostCandsToken_(
0035 consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("lostCandidatesTag"))) {}
0036
0037 PFTauMiniAODPrimaryVertexProducer::~PFTauMiniAODPrimaryVertexProducer() {}
0038
0039 void PFTauMiniAODPrimaryVertexProducer::beginEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0040
0041 iEvent.getByToken(packedCandsToken_, packedCands_);
0042 iEvent.getByToken(lostCandsToken_, lostCands_);
0043 }
0044
0045 void PFTauMiniAODPrimaryVertexProducer::nonTauTracksInPV(const reco::VertexRef &thePVRef,
0046 const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
0047 std::vector<const reco::Track *> &nonTauTracks) {
0048
0049
0050 if (packedCands_.isValid()) {
0051 nonTauTracksInPVFromPackedCands(thePVRef.key(), *packedCands_, tauTracks, nonTauTracks);
0052 }
0053
0054 if (lostCands_.isValid()) {
0055 nonTauTracksInPVFromPackedCands(thePVRef.key(), *lostCands_, tauTracks, nonTauTracks);
0056 }
0057 }
0058
0059 void PFTauMiniAODPrimaryVertexProducer::nonTauTracksInPVFromPackedCands(
0060 const size_t &thePVkey,
0061 const pat::PackedCandidateCollection &cands,
0062 const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
0063 std::vector<const reco::Track *> &nonTauTracks) {
0064
0065 for (const auto &cand : cands) {
0066
0067 if (!std::isfinite(cand.pt()))
0068 continue;
0069 if (cand.vertexRef().isNull())
0070 continue;
0071 int quality = cand.pvAssociationQuality();
0072 if (cand.vertexRef().key() != thePVkey ||
0073 (quality != pat::PackedCandidate::UsedInFitTight && quality != pat::PackedCandidate::UsedInFitLoose))
0074 continue;
0075 const reco::Track *track = cand.bestTrack();
0076 if (track == nullptr)
0077 continue;
0078
0079
0080
0081 bool matched = false;
0082 for (const auto &tauTrack : tauTracks) {
0083 if (std::abs(tauTrack->eta() - track->eta()) < 0.005 &&
0084 std::abs(deltaPhi(tauTrack->phi(), track->phi())) < 0.005 &&
0085 std::abs(tauTrack->pt() / track->pt() - 1.) < 0.005) {
0086 matched = true;
0087 break;
0088 }
0089 }
0090 if (!matched)
0091 nonTauTracks.push_back(track);
0092 }
0093 }
0094
0095 void PFTauMiniAODPrimaryVertexProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0096 auto desc = PFTauPrimaryVertexProducerBase::getDescriptionsBase();
0097 desc.add<edm::InputTag>("lostCandidatesTag", edm::InputTag("lostTracks"));
0098 desc.add<edm::InputTag>("packedCandidatesTag", edm::InputTag("packedPFCandidates"));
0099
0100 descriptions.add("pfTauMiniAODPrimaryVertexProducer", desc);
0101 }
0102
0103 DEFINE_FWK_MODULE(PFTauMiniAODPrimaryVertexProducer);