CMS 3D CMS Logo

PFTauMiniAODPrimaryVertexProducer.cc
Go to the documentation of this file.
5 
8 
9  public:
12 
13  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
14  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
15 
16  protected:
18  const std::vector<edm::Ptr<reco::TrackBase> >&,
19  std::vector<const reco::Track*>&) override;
20 
21  private:
22  void nonTauTracksInPVFromPackedCands(const size_t&,
24  const std::vector<edm::Ptr<reco::TrackBase> >&,
25  std::vector<const reco::Track*> &);
26 
29 
30 };
31 
34  packedCandsToken_(consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("packedCandidatesTag"))),
35  lostCandsToken_(consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("lostCandidatesTag"))) {}
36 
38 
40 
41  //Get candidate collections
44 }
45 
47  const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
48  std::vector<const reco::Track*> &nonTauTracks){
49 
50  //Find non-tau tracks associated to thePV
51  //PackedCandidates first...
52  if(packedCands_.isValid()) {
53  nonTauTracksInPVFromPackedCands(thePVRef.key(),*packedCands_,tauTracks,nonTauTracks);
54  }
55  //then lostCandidates
56  if(lostCands_.isValid()) {
57  nonTauTracksInPVFromPackedCands(thePVRef.key(),*lostCands_,tauTracks,nonTauTracks);
58  }
59 }
60 
63  const std::vector<edm::Ptr<reco::TrackBase> > &tauTracks,
64  std::vector<const reco::Track*> &nonTauTracks){
65 
66  //Find candidates/tracks associated to thePV
67  for(const auto& cand: cands){
68  if(cand.vertexRef().isNull()) continue;
69  int quality = cand.pvAssociationQuality();
70  if(cand.vertexRef().key()!=thePVkey ||
72  quality!=pat::PackedCandidate::UsedInFitLoose)) continue;
73  const reco::Track *track = cand.bestTrack();
74  if(track == nullptr) 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  ){
84  matched = true;
85  break;
86  }
87  }
88  if( !matched ) nonTauTracks.push_back(track);
89  }
90 }
91 
92 void
95  desc.add<edm::InputTag>("lostCandidatesTag", edm::InputTag("lostTracks"));
96  desc.add<edm::InputTag>("packedCandidatesTag", edm::InputTag("packedPFCandidates"));
97 
98  descriptions.add("pfTauMiniAODPrimaryVertexProducer", desc);
99 }
100 
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:517
edm::EDGetTokenT< pat::PackedCandidateCollection > packedCandsToken_
std::vector< pat::PackedCandidate > PackedCandidateCollection
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:678
key_type key() const
Accessor for product key.
Definition: Ref.h:263
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:684
double pt() const
track transverse momentum
Definition: TrackBase.h:654
MiniAOD implementation of the PFTauPrimaryVertexProducer plugin.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:74
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_