Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /* class PFTauSecondaryVertexProducer
0002  * EDProducer of the 
0003  * authors: Ian M. Nugent
0004  * This work is based on the impact parameter work by Rosamaria Venditti and reconstructing the 3 prong taus.
0005  * The idea of the fully reconstructing the tau using a kinematic fit comes from
0006  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
0007  * work was continued by Ian M. Nugent and Vladimir Cherepanov.
0008  * Thanks goes to Christian Veelken and Evan Klose Friis for their help and suggestions.
0009  */
0010 
0011 #include "DataFormats/TauReco/interface/PFTauTagInfo.h"
0012 #include "DataFormats/VertexReco/interface/Vertex.h"
0013 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0014 
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/Framework/interface/global/EDProducer.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/Utilities/interface/Exception.h"
0023 
0024 #include <FWCore/ParameterSet/interface/ConfigurationDescriptions.h>
0025 #include <FWCore/ParameterSet/interface/ParameterSetDescription.h>
0026 
0027 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0028 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0029 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0030 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0031 
0032 #include "DataFormats/TauReco/interface/PFTau.h"
0033 #include "DataFormats/TauReco/interface/PFTauFwd.h"
0034 #include "DataFormats/VertexReco/interface/Vertex.h"
0035 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0036 #include "DataFormats/TrackReco/interface/Track.h"
0037 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0038 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0039 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0040 #include "DataFormats/Common/interface/RefToBase.h"
0041 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0042 
0043 #include "DataFormats/Common/interface/Association.h"
0044 #include "DataFormats/Common/interface/AssociationVector.h"
0045 #include "DataFormats/Common/interface/RefProd.h"
0046 
0047 #include <memory>
0048 
0049 using namespace reco;
0050 using namespace edm;
0051 using namespace std;
0052 
0053 class PFTauSecondaryVertexProducer : public edm::global::EDProducer<> {
0054 public:
0055   enum Alg { useInputPV = 0, usePVwithMaxSumPt, useTauPV };
0056 
0057   explicit PFTauSecondaryVertexProducer(const edm::ParameterSet& iConfig);
0058   ~PFTauSecondaryVertexProducer() override;
0059   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0060   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0061 
0062 private:
0063   const edm::InputTag PFTauTag_;
0064   const edm::EDGetTokenT<std::vector<reco::PFTau>> PFTauToken_;
0065   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> transTrackBuilderToken_;
0066 };
0067 
0068 PFTauSecondaryVertexProducer::PFTauSecondaryVertexProducer(const edm::ParameterSet& iConfig)
0069     : PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
0070       PFTauToken_(consumes<std::vector<reco::PFTau>>(iConfig.getParameter<edm::InputTag>("PFTauTag"))),
0071       transTrackBuilderToken_(esConsumes(edm::ESInputTag{"", "TransientTrackBuilder"})) {
0072   produces<edm::AssociationVector<PFTauRefProd, std::vector<std::vector<reco::VertexRef>>>>();
0073   produces<VertexCollection>("PFTauSecondaryVertices");
0074 }
0075 
0076 PFTauSecondaryVertexProducer::~PFTauSecondaryVertexProducer() {}
0077 
0078 namespace {
0079   const reco::Track* getTrack(const reco::Candidate& cand) {
0080     const reco::PFCandidate* pfCand = dynamic_cast<const reco::PFCandidate*>(&cand);
0081     if (pfCand != nullptr) {
0082       if (pfCand->trackRef().isNonnull())
0083         return &*pfCand->trackRef();
0084       else if (pfCand->gsfTrackRef().isNonnull())
0085         return &*pfCand->gsfTrackRef();
0086     }
0087     const pat::PackedCandidate* pCand = dynamic_cast<const pat::PackedCandidate*>(&cand);
0088     if (pCand != nullptr && pCand->hasTrackDetails())
0089       return &pCand->pseudoTrack();
0090     return nullptr;
0091   }
0092 }  // namespace
0093 void PFTauSecondaryVertexProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0094   // Obtain
0095   auto const& transTrackBuilder = iSetup.getData(transTrackBuilderToken_);
0096 
0097   edm::Handle<std::vector<reco::PFTau>> Tau;
0098   iEvent.getByToken(PFTauToken_, Tau);
0099 
0100   // Set Association Map
0101   auto AVPFTauSV = std::make_unique<edm::AssociationVector<PFTauRefProd, std::vector<std::vector<reco::VertexRef>>>>(
0102       PFTauRefProd(Tau));
0103   auto VertexCollection_out = std::make_unique<VertexCollection>();
0104   reco::VertexRefProd VertexRefProd_out = iEvent.getRefBeforePut<reco::VertexCollection>("PFTauSecondaryVertices");
0105 
0106   // For each Tau Run Algorithim
0107   if (Tau.isValid()) {
0108     for (reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
0109       reco::PFTauRef RefPFTau(Tau, iPFTau);
0110       std::vector<reco::VertexRef> SV;
0111       if (RefPFTau->decayMode() >= 5) {
0112         ///////////////////////////////////////////////////////////////////////////////////////////////
0113         // Get tracks form PFTau daugthers
0114         std::vector<reco::TransientTrack> transTrk;
0115         TransientVertex transVtx;
0116         const std::vector<edm::Ptr<reco::Candidate>> cands = RefPFTau->signalChargedHadrCands();
0117         for (const auto& cand : cands) {
0118           if (cand.isNull())
0119             continue;
0120           const reco::Track* track = getTrack(*cand);
0121           if (track != nullptr)
0122             transTrk.push_back(transTrackBuilder.build(*track));
0123         }
0124         ///////////////////////////////////////////////////////////////////////////////////////////////
0125         // Fit the secondary vertex
0126         bool FitOk(true);
0127         KalmanVertexFitter kvf(true);
0128         if (transTrk.size() > 1) {
0129           transVtx = kvf.vertex(transTrk);  //KalmanVertexFitter
0130         } else {
0131           FitOk = false;
0132         }
0133         if (!transVtx.hasRefittedTracks())
0134           FitOk = false;
0135         if (transVtx.refittedTracks().size() != transTrk.size())
0136           FitOk = false;
0137         if (FitOk) {
0138           SV.push_back(reco::VertexRef(VertexRefProd_out, VertexCollection_out->size()));
0139           VertexCollection_out->push_back(transVtx);
0140         }
0141       }
0142       AVPFTauSV->setValue(iPFTau, SV);
0143     }
0144   }
0145   iEvent.put(std::move(VertexCollection_out), "PFTauSecondaryVertices");
0146   iEvent.put(std::move(AVPFTauSV));
0147 }
0148 
0149 void PFTauSecondaryVertexProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0150   // PFTauSecondaryVertexProducer
0151   edm::ParameterSetDescription desc;
0152   desc.add<edm::InputTag>("PFTauTag", edm::InputTag("hpsPFTauProducer"));
0153   descriptions.add("PFTauSecondaryVertexProducer", desc);
0154 }
0155 
0156 DEFINE_FWK_MODULE(PFTauSecondaryVertexProducer);