CMS 3D CMS Logo

GeneralTracksImporterWithVeto.cc
Go to the documentation of this file.
12 
14 public:
16  : BlockElementImporterBase(conf, sumes),
17  src_(sumes.consumes<reco::PFRecTrackCollection>(conf.getParameter<edm::InputTag>("source"))),
18  veto_(sumes.consumes<reco::PFRecTrackCollection>(conf.getParameter<edm::InputTag>("veto"))),
19  muons_(sumes.consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
20  trackQuality_((conf.existsAs<std::string>("trackQuality"))
21  ? reco::TrackBase::qualityByName(conf.getParameter<std::string>("trackQuality"))
22  : reco::TrackBase::highPurity),
23  DPtovPtCut_(conf.getParameter<std::vector<double> >("DPtOverPtCuts_byTrackAlgo")),
24  NHitCut_(conf.getParameter<std::vector<unsigned> >("NHitCuts_byTrackAlgo")),
25  useIterTracking_(conf.getParameter<bool>("useIterativeTracking")),
27  conf.existsAs<bool>("cleanBadConvertedBrems") ? conf.getParameter<bool>("cleanBadConvertedBrems") : false) {
28  bool postMuonCleaning =
29  conf.existsAs<bool>("postMuonCleaning") ? conf.getParameter<bool>("postMuonCleaning") : false;
30  pfmu_ = std::unique_ptr<PFMuonAlgo>(new PFMuonAlgo(conf, postMuonCleaning));
31  }
32 
33  void importToBlock(const edm::Event&, ElementList&) const override;
34 
35 private:
36  int muAssocToTrack(const reco::TrackRef& trackref, const edm::Handle<reco::MuonCollection>& muonh) const;
37 
41  const std::vector<double> DPtovPtCut_;
42  const std::vector<unsigned> NHitCut_;
44 
45  std::unique_ptr<PFMuonAlgo> pfmu_;
46 };
47 
49 
53  auto tracks = e.getHandle(src_);
54  auto vetosH = e.getHandle(veto_);
55  const auto& vetos = *vetosH;
56  std::unordered_set<unsigned> vetoed;
57  for (unsigned i = 0; i < vetos.size(); ++i) {
58  vetoed.insert(vetos[i].trackRef().key());
59  }
60  auto muons = e.getHandle(muons_);
61  elems.reserve(elems.size() + tracks->size());
62  std::vector<bool> mask(tracks->size(), true);
63  reco::MuonRef muonref;
64 
65  // remove converted brems with bad pT resolution if requested
66  // this reproduces the old behavior of PFBlockAlgo
67  if (cleanBadConvBrems_) {
68  auto itr = elems.begin();
69  while (itr != elems.end()) {
70  if ((*itr)->type() == reco::PFBlockElement::TRACK) {
71  const reco::PFBlockElementTrack* trkel = static_cast<reco::PFBlockElementTrack*>(itr->get());
72  const reco::ConversionRefVector& cRef = trkel->convRefs();
74  const reco::VertexCompositeCandidateRef& v0Ref = trkel->V0Ref();
75  // if there is no displaced vertex reference and it is marked
76  // as a conversion it's gotta be a converted brem
77  if (trkel->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) && cRef.empty() && dvRef.isNull() &&
78  v0Ref.isNull()) {
79  // if the Pt resolution is bad we kill this element
82  itr = elems.erase(itr);
83  continue;
84  }
85  }
86  }
87  ++itr;
88  } // loop on existing elements
89  }
90  // preprocess existing tracks in the element list and create a mask
91  // so that we do not import tracks twice, tag muons we find
92  // in this collection
93  auto TKs_end = std::partition(
94  elems.begin(), elems.end(), [](const ElementType& a) { return a->type() == reco::PFBlockElement::TRACK; });
95  auto btk_elems = elems.begin();
96  auto btrack = tracks->cbegin();
97  auto etrack = tracks->cend();
98  for (auto track = btrack; track != etrack; ++track) {
99  auto tk_elem =
100  std::find_if(btk_elems, TKs_end, [&](const ElementType& a) { return (a->trackRef() == track->trackRef()); });
101  if (tk_elem != TKs_end) {
102  mask[std::distance(tracks->cbegin(), track)] = false;
103  // check and update if this track is a muon
104  const int muId = muAssocToTrack((*tk_elem)->trackRef(), muons);
105  if (muId != -1) {
106  muonref = reco::MuonRef(muons, muId);
107  if (PFMuonAlgo::isLooseMuon(muonref) || PFMuonAlgo::isMuon(muonref)) {
108  static_cast<reco::PFBlockElementTrack*>(tk_elem->get())->setMuonRef(muonref);
109  }
110  }
111  }
112  }
113  // now we actually insert tracks, again tagging muons along the way
114  reco::PFRecTrackRef pftrackref;
115  reco::PFBlockElementTrack* trkElem = nullptr;
116  for (auto track = btrack; track != etrack; ++track) {
117  const unsigned idx = std::distance(btrack, track);
118  // since we already set muon refs in the previously imported tracks,
119  // here we can skip everything that is already imported
120  if (!mask[idx])
121  continue;
122  muonref = reco::MuonRef();
123  pftrackref = reco::PFRecTrackRef(tracks, idx);
124  // Get the eventual muon associated to this track
125  const int muId = muAssocToTrack(pftrackref->trackRef(), muons);
126  bool thisIsAPotentialMuon = false;
127  if (muId != -1) {
128  muonref = reco::MuonRef(muons, muId);
129  thisIsAPotentialMuon = ((pfmu_->hasValidTrack(muonref, true) && PFMuonAlgo::isLooseMuon(muonref)) ||
130  (pfmu_->hasValidTrack(muonref, false) && PFMuonAlgo::isMuon(muonref)));
131  }
132  if (thisIsAPotentialMuon || PFTrackAlgoTools::goodPtResolution(
133  pftrackref->trackRef(), DPtovPtCut_, NHitCut_, useIterTracking_, trackQuality_)) {
134  trkElem = new reco::PFBlockElementTrack(pftrackref);
135  if (thisIsAPotentialMuon) {
136  LogDebug("GeneralTracksImporterWithVeto")
137  << "Potential Muon P " << pftrackref->trackRef()->p() << " pt " << pftrackref->trackRef()->p() << std::endl;
138  }
139  if (muId != -1)
140  trkElem->setMuonRef(muonref);
141  if (vetoed.count(pftrackref->trackRef().key()) == 0 || muonref.isNonnull()) {
142  elems.emplace_back(trkElem);
143  } else
144  delete trkElem;
145  }
146  }
147  elems.shrink_to_fit();
148 }
149 
151  const edm::Handle<reco::MuonCollection>& muonh) const {
152  auto muon = std::find_if(muonh->cbegin(), muonh->cend(), [&](const reco::Muon& m) {
153  return (m.track().isNonnull() && m.track() == trackref);
154  });
155  return (muon != muonh->cend() ? std::distance(muonh->cbegin(), muon) : -1);
156 }
#define LogDebug(id)
T getParameter(std::string const &) const
int muAssocToTrack(const reco::TrackRef &trackref, const edm::Handle< reco::MuonCollection > &muonh) const
bool goodPtResolution(const reco::TrackRef &, const std::vector< double > &DPtovPtCut, const std::vector< unsigned > &NHitCut, bool useIterTracking, const reco::TrackBase::TrackQuality trackQuality)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:160
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:56
TrackQuality
track quality
Definition: TrackBase.h:150
edm::EDGetTokenT< reco::MuonCollection > muons_
const reco::TrackRef & trackRef() const override
const ConversionRefVector & convRefs() const override
bool trackType(TrackType trType) const override
edm::EDGetTokenT< reco::PFRecTrackCollection > veto_
key_type key() const
Accessor for product key.
Definition: Ref.h:250
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:547
edm::EDGetTokenT< reco::PFRecTrackCollection > src_
edm::Ref< PFRecTrackCollection > PFRecTrackRef
persistent reference to PFRecTrack objects
Definition: PFRecTrackFwd.h:15
const VertexCompositeCandidateRef & V0Ref() const override
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
bool isNull() const
Checks for null.
Definition: Ref.h:235
const PFDisplacedTrackerVertexRef & displacedVertexRef(TrackType trType) const override
GeneralTracksImporterWithVeto(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
void importToBlock(const edm::Event &, ElementList &) const override
static bool isLooseMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:65
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:119
void setMuonRef(const MuonRef &muref) override
reference to the Muon
const reco::TrackBase::TrackQuality trackQuality_
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
const std::vector< unsigned > NHitCut_