CMS 3D CMS Logo

PFTauSecondaryVertexProducer.cc
Go to the documentation of this file.
1 /* class PFTauSecondaryVertexProducer
2  * EDProducer of the
3  * authors: Ian M. Nugent
4  * This work is based on the impact parameter work by Rosamaria Venditti and reconstructing the 3 prong taus.
5  * The idea of the fully reconstructing the tau using a kinematic fit comes from
6  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
7  * work was continued by Ian M. Nugent and Vladimir Cherepanov.
8  * Thanks goes to Christian Veelken and Evan Klose Friis for their help and suggestions.
9  */
10 
14 
24 
27 
32 
43 
47 
48 #include <memory>
49 
50 using namespace reco;
51 using namespace edm;
52 using namespace std;
53 
55 public:
56  enum Alg { useInputPV = 0, usePVwithMaxSumPt, useTauPV };
57 
58  explicit PFTauSecondaryVertexProducer(const edm::ParameterSet& iConfig);
59  ~PFTauSecondaryVertexProducer() override;
60  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
61  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
62 
63 private:
66 };
67 
69  : PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
70  PFTauToken_(consumes<std::vector<reco::PFTau>>(iConfig.getParameter<edm::InputTag>("PFTauTag"))) {
71  produces<edm::AssociationVector<PFTauRefProd, std::vector<std::vector<reco::VertexRef>>>>();
72  produces<VertexCollection>("PFTauSecondaryVertices");
73 }
74 
76 
77 namespace {
78  const reco::Track* getTrack(const reco::Candidate& cand) {
79  const reco::PFCandidate* pfCand = dynamic_cast<const reco::PFCandidate*>(&cand);
80  if (pfCand != nullptr) {
81  if (pfCand->trackRef().isNonnull())
82  return &*pfCand->trackRef();
83  else if (pfCand->gsfTrackRef().isNonnull())
84  return &*pfCand->gsfTrackRef();
85  }
86  const pat::PackedCandidate* pCand = dynamic_cast<const pat::PackedCandidate*>(&cand);
87  if (pCand != nullptr && pCand->hasTrackDetails())
88  return &pCand->pseudoTrack();
89  return nullptr;
90  }
91 } // namespace
93  // Obtain
94  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
95  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", transTrackBuilder);
96 
98  iEvent.getByToken(PFTauToken_, Tau);
99 
100  // Set Association Map
101  auto AVPFTauSV = std::make_unique<edm::AssociationVector<PFTauRefProd, std::vector<std::vector<reco::VertexRef>>>>(
102  PFTauRefProd(Tau));
103  auto VertexCollection_out = std::make_unique<VertexCollection>();
104  reco::VertexRefProd VertexRefProd_out = iEvent.getRefBeforePut<reco::VertexCollection>("PFTauSecondaryVertices");
105 
106  // For each Tau Run Algorithim
107  if (Tau.isValid()) {
108  for (reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
109  reco::PFTauRef RefPFTau(Tau, iPFTau);
110  std::vector<reco::VertexRef> SV;
111  if (RefPFTau->decayMode() >= 5) {
113  // Get tracks form PFTau daugthers
114  std::vector<reco::TransientTrack> transTrk;
115  TransientVertex transVtx;
116  const std::vector<edm::Ptr<reco::Candidate>> cands = RefPFTau->signalChargedHadrCands();
117  for (const auto& cand : cands) {
118  if (cand.isNull())
119  continue;
120  const reco::Track* track = getTrack(*cand);
121  if (track != nullptr)
122  transTrk.push_back(transTrackBuilder->build(*track));
123  }
125  // Fit the secondary vertex
126  bool FitOk(true);
127  KalmanVertexFitter kvf(true);
128  if (transTrk.size() > 1) {
129  transVtx = kvf.vertex(transTrk); //KalmanVertexFitter
130  } else {
131  FitOk = false;
132  }
133  if (!transVtx.hasRefittedTracks())
134  FitOk = false;
135  if (transVtx.refittedTracks().size() != transTrk.size())
136  FitOk = false;
137  if (FitOk) {
138  SV.push_back(reco::VertexRef(VertexRefProd_out, VertexCollection_out->size()));
139  VertexCollection_out->push_back(transVtx);
140  }
141  }
142  AVPFTauSV->setValue(iPFTau, SV);
143  }
144  }
145  iEvent.put(std::move(VertexCollection_out), "PFTauSecondaryVertices");
146  iEvent.put(std::move(AVPFTauSV));
147 }
148 
150  // PFTauSecondaryVertexProducer
152  desc.add<edm::InputTag>("PFTauTag", edm::InputTag("hpsPFTauProducer"));
153  descriptions.add("PFTauSecondaryVertexProducer", desc);
154 }
155 
ConfigurationDescriptions.h
reco::PFCandidate::trackRef
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:408
edm::RefProd< VertexCollection >
RefProd.h
edm::StreamID
Definition: StreamID.h:30
pat::PackedCandidate::hasTrackDetails
bool hasTrackDetails() const
Return true if a bestTrack can be extracted from this Candidate.
Definition: PackedCandidate.h:787
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
TransientVertex::hasRefittedTracks
bool hasRefittedTracks() const
Definition: TransientVertex.h:206
PFTauFwd.h
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
MessageLogger.h
TransientVertex::refittedTracks
std::vector< reco::TransientTrack > const & refittedTracks() const
Definition: TransientVertex.h:211
KalmanVertexFitter.h
ESHandle.h
PFTauSecondaryVertexProducer_cfi.PFTauSecondaryVertexProducer
PFTauSecondaryVertexProducer
Definition: PFTauSecondaryVertexProducer_cfi.py:3
edm::EDGetTokenT
Definition: EDGetToken.h:33
L1TRate_Offline_cfi.Tau
Tau
Definition: L1TRate_Offline_cfi.py:43
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
reco::PFTauRefProd
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
PFTauSecondaryVertexProducer::~PFTauSecondaryVertexProducer
~PFTauSecondaryVertexProducer() override
Definition: PFTauSecondaryVertexProducer.cc:75
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
reco::PFTau
Definition: PFTau.h:36
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
Tau
Definition: Tau.py:1
Association.h
edm::Handle
Definition: AssociativeIterator.h:50
edm::Ref< PFTauCollection >
PFTauSecondaryVertexProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PFTauSecondaryVertexProducer.cc:149
MakerMacros.h
trigger::size_type
uint16_t size_type
Definition: TriggerTypeDefs.h:18
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TransientTrackRecord
Definition: TransientTrackRecord.h:11
reco::Track
Definition: Track.h:27
edm::ESHandle< TransientTrackBuilder >
ParameterSetDescription.h
RefToBase.h
PFTauSecondaryVertexProducer::Alg
Alg
Definition: PFTauSecondaryVertexProducer.cc:56
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
Vertex.h
pat::PackedCandidate::pseudoTrack
virtual const reco::Track & pseudoTrack() const
Definition: PackedCandidate.h:772
HLT_FULL_cff.cands
cands
Definition: HLT_FULL_cff.py:15142
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:47
reco::PFCandidate::gsfTrackRef
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:440
Event.h
PFTauSecondaryVertexProducer::PFTauSecondaryVertexProducer
PFTauSecondaryVertexProducer(const edm::ParameterSet &iConfig)
Definition: PFTauSecondaryVertexProducer.cc:68
pat::PackedCandidate
Definition: PackedCandidate.h:22
PackedCandidate.h
edm::Ref::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
cand
Definition: decayParser.h:32
nanoDQM_cfi.SV
SV
Definition: nanoDQM_cfi.py:616
iEvent
int iEvent
Definition: GenABIO.cc:224
GsfTrack.h
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
TransientTrackRecord.h
get
#define get
reco::Candidate
Definition: Candidate.h:27
VertexFwd.h
getTrack
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
Definition: GhostTrackState.cc:49
TransientVertex.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
PFTauSecondaryVertexProducer::PFTauToken_
const edm::EDGetTokenT< std::vector< reco::PFTau > > PFTauToken_
Definition: PFTauSecondaryVertexProducer.cc:65
GsfTrackFwd.h
Frameworkfwd.h
PFTauSecondaryVertexProducer
Definition: PFTauSecondaryVertexProducer.cc:54
PFTau.h
EventSetup.h
Exception.h
PFTauTagInfo.h
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
AssociationVector.h
ParameterSet.h
EDProducer.h
PFTauSecondaryVertexProducer::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: PFTauSecondaryVertexProducer.cc:92
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
PFTauSecondaryVertexProducer::PFTauTag_
const edm::InputTag PFTauTag_
Definition: PFTauSecondaryVertexProducer.cc:64
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22