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 &,
17  const std::vector<edm::Ptr<reco::TrackBase> > &,
18  std::vector<const reco::Track *> &) override;
19 
20 private:
21  void nonTauTracksInPVFromPackedCands(const size_t &,
23  const std::vector<edm::Ptr<reco::TrackBase> > &,
24  std::vector<const reco::Track *> &);
25 
28 };
29 
33  consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("packedCandidatesTag"))),
35  consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("lostCandidatesTag"))) {}
36 
38 
40  //Get candidate collections
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  if (cand.vertexRef().isNull())
67  continue;
68  int quality = cand.pvAssociationQuality();
69  if (cand.vertexRef().key() != thePVkey ||
71  continue;
72  const reco::Track *track = cand.bestTrack();
73  if (track == nullptr)
74  continue;
75  //Remove signal (tau) tracks
76  //MB: Only deltaR deltaPt overlap removal possible (?)
77  //MB: It should be fine as pat objects stores same track info with same presision
78  bool matched = false;
79  for (const auto &tauTrack : tauTracks) {
80  if (std::abs(tauTrack->eta() - track->eta()) < 0.005 &&
81  std::abs(deltaPhi(tauTrack->phi(), track->phi())) < 0.005 &&
82  std::abs(tauTrack->pt() / track->pt() - 1.) < 0.005) {
83  matched = true;
84  break;
85  }
86  }
87  if (!matched)
88  nonTauTracks.push_back(track);
89  }
90 }
91 
94  desc.add<edm::InputTag>("lostCandidatesTag", edm::InputTag("lostTracks"));
95  desc.add<edm::InputTag>("packedCandidatesTag", edm::InputTag("packedPFCandidates"));
96 
97  descriptions.add("pfTauMiniAODPrimaryVertexProducer", desc);
98 }
99 
void beginEvent(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< pat::PackedCandidateCollection > lostCandsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< pat::PackedCandidateCollection > packedCandsToken_
std::vector< pat::PackedCandidate > PackedCandidateCollection
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
key_type key() const
Accessor for product key.
Definition: Ref.h:250
PFTauMiniAODPrimaryVertexProducer(const edm::ParameterSet &iConfig)
Definition: HeavyIon.h:7
void nonTauTracksInPVFromPackedCands(const size_t &, const pat::PackedCandidateCollection &, const std::vector< edm::Ptr< reco::TrackBase > > &, std::vector< const reco::Track * > &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:617
double pt() const
track transverse momentum
Definition: TrackBase.h:602
MiniAOD implementation of the PFTauPrimaryVertexProducer plugin.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
void nonTauTracksInPV(const reco::VertexRef &, const std::vector< edm::Ptr< reco::TrackBase > > &, std::vector< const reco::Track * > &) override
edm::Handle< pat::PackedCandidateCollection > packedCands_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static edm::ParameterSetDescription getDescriptionsBase()
HLT enums.
edm::Handle< pat::PackedCandidateCollection > lostCands_