CMS 3D CMS Logo

PFTauMiniAODPrimaryVertexProducer.cc
Go to the documentation of this file.
5 
8 public:
11 
12  void beginEvent(const edm::Event &, const edm::EventSetup &) override;
13  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
14 
15 protected:
16  void nonTauTracksInPV(const reco::VertexRef &,
18  std::vector<const reco::Track *> &) override;
19 
20 private:
21  void nonTauTracksInPVFromPackedCands(const size_t &,
24  std::vector<const reco::Track *> &);
25 
28 };
29 
32  packedCandsToken_(
33  consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("packedCandidatesTag"))),
34  lostCandsToken_(
35  consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("lostCandidatesTag"))) {}
36 
38 
40  //Get candidate collections
42  iEvent.getByToken(lostCandsToken_, lostCands_);
43 }
44 
46  const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
47  std::vector<const reco::Track *> &nonTauTracks) {
48  //Find non-tau tracks associated to thePV
49  //PackedCandidates first...
50  if (packedCands_.isValid()) {
51  nonTauTracksInPVFromPackedCands(thePVRef.key(), *packedCands_, tauTracks, nonTauTracks);
52  }
53  //then lostCandidates
54  if (lostCands_.isValid()) {
55  nonTauTracksInPVFromPackedCands(thePVRef.key(), *lostCands_, tauTracks, nonTauTracks);
56  }
57 }
58 
60  const size_t &thePVkey,
62  const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
63  std::vector<const reco::Track *> &nonTauTracks) {
64  //Find candidates/tracks associated to thePV
65  for (const auto &cand : cands) {
66  //MB: Skip candidates with ill-defined momentum as they return ill-defined tracks (why it happens?)
67  if (!std::isfinite(cand.pt())) //MB: it is enough to check just pt (?)
68  continue;
69  if (cand.vertexRef().isNull())
70  continue;
71  int quality = cand.pvAssociationQuality();
72  if (cand.vertexRef().key() != thePVkey ||
74  continue;
75  const reco::Track *track = cand.bestTrack();
76  if (track == nullptr)
77  continue;
78  //Remove signal (tau) tracks
79  //MB: Only deltaR deltaPt overlap removal possible (?)
80  //MB: It should be fine as pat objects stores same track info with same presision
81  bool matched = false;
82  for (const auto &tauTrack : tauTracks) {
83  if (std::abs(tauTrack->eta() - track->eta()) < 0.005 &&
84  std::abs(deltaPhi(tauTrack->phi(), track->phi())) < 0.005 &&
85  std::abs(tauTrack->pt() / track->pt() - 1.) < 0.005) {
86  matched = true;
87  break;
88  }
89  }
90  if (!matched)
91  nonTauTracks.push_back(track);
92  }
93 }
94 
97  desc.add<edm::InputTag>("lostCandidatesTag", edm::InputTag("lostTracks"));
98  desc.add<edm::InputTag>("packedCandidatesTag", edm::InputTag("packedPFCandidates"));
99 
100  descriptions.add("pfTauMiniAODPrimaryVertexProducer", desc);
101 }
102 
void beginEvent(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< pat::PackedCandidateCollection > lostCandsToken_
edm::EDGetTokenT< pat::PackedCandidateCollection > packedCandsToken_
std::vector< pat::PackedCandidate > PackedCandidateCollection
void nonTauTracksInPVFromPackedCands(const size_t &, const pat::PackedCandidateCollection &, const std::vector< edm::Ptr< reco::TrackBase > > &, std::vector< const reco::Track *> &)
key_type key() const
Accessor for product key.
Definition: Ref.h:250
PFTauMiniAODPrimaryVertexProducer(const edm::ParameterSet &iConfig)
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
MiniAOD implementation of the PFTauPrimaryVertexProducer plugin.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Handle< pat::PackedCandidateCollection > packedCands_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void nonTauTracksInPV(const reco::VertexRef &, const std::vector< edm::Ptr< reco::TrackBase > > &, std::vector< const reco::Track *> &) override
bool isValid() const
Definition: HandleBase.h:70
static edm::ParameterSetDescription getDescriptionsBase()
HLT enums.
edm::Handle< pat::PackedCandidateCollection > lostCands_
string quality