CMS 3D CMS Logo

GeneralTracksImporter.cc
Go to the documentation of this file.
11 
13 public:
15  : BlockElementImporterBase(conf, cc),
16  src_(cc.consumes<reco::PFRecTrackCollection>(conf.getParameter<edm::InputTag>("source"))),
17  vetoEndcap_(conf.getParameter<bool>("vetoEndcap")),
18  muons_(cc.consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
19  trackQuality_(reco::TrackBase::qualityByName(conf.getParameter<std::string>("trackQuality"))),
20  DPtovPtCut_(conf.getParameter<std::vector<double>>("DPtOverPtCuts_byTrackAlgo")),
21  NHitCut_(conf.getParameter<std::vector<unsigned>>("NHitCuts_byTrackAlgo")),
22  useIterTracking_(conf.getParameter<bool>("useIterativeTracking")),
23  cleanBadConvBrems_(conf.getParameter<bool>("cleanBadConvertedBrems")),
24  muonMaxDPtOPt_(conf.getParameter<double>("muonMaxDPtOPt")) {
25  if (vetoEndcap_) {
26  vetoMode_ = conf.getParameter<unsigned>("vetoMode");
27  switch (vetoMode_) {
30  break;
31  case ticlSeedingRegion:
33  cc.consumes<std::vector<TICLSeedingRegion>>(conf.getParameter<edm::InputTag>("vetoSrc"));
35  break;
38  break;
39  } // switch
40  }
41  }
42 
43  void importToBlock(const edm::Event&, ElementList&) const override;
44 
45 private:
47  const bool vetoEndcap_;
50  const std::vector<double> DPtovPtCut_;
51  const std::vector<unsigned> NHitCut_;
53  const double muonMaxDPtOPt_;
54  unsigned int vetoMode_;
59 };
60 
62 
65  auto tracks = e.getHandle(src_);
66 
67  typedef std::pair<edm::ProductID, unsigned> TrackProdIDKey;
68  std::vector<TrackProdIDKey> vetoed;
69  if (vetoEndcap_) {
70  switch (vetoMode_) {
71  case pfRecTrackCollection: {
72  const auto& vetoes = e.get(vetoPFTracksSrc_);
73  for (const auto& veto : vetoes)
74  vetoed.emplace_back(veto.trackRef().id(), veto.trackRef().key());
75  break;
76  }
77  case ticlSeedingRegion: {
78  const auto& vetoes = e.get(vetoTICLSeedingSrc_);
79  auto tracksH = e.getHandle(tracksSrc_);
80  for (const auto& veto : vetoes) {
81  assert(veto.collectionID == tracksH.id());
82  reco::TrackRef trkref = reco::TrackRef(tracksH, veto.index);
83  vetoed.emplace_back(tracksH.id(), veto.index); // track prod id and key
84  }
85  break;
86  }
87  case pfCandidateCollection: {
88  const auto& vetoes = e.get(vetoPFCandidatesSrc_);
89  for (const auto& veto : vetoes) {
90  if (veto.trackRef().isNull())
91  continue;
92  vetoed.emplace_back(veto.trackRef().id(), veto.trackRef().key());
93  }
94  break;
95  }
96  } // switch
97  std::sort(vetoed.begin(), vetoed.end());
98  }
99  const auto muonH = e.getHandle(muons_);
100  const auto muons = *muonH;
101  elems.reserve(elems.size() + tracks->size());
102  std::vector<bool> mask(tracks->size(), true);
103  reco::MuonRef muonref;
104 
105  // remove converted brems with bad pT resolution if requested
106  // this reproduces the old behavior of PFBlockAlgo
107  if (cleanBadConvBrems_) {
108  auto itr = elems.begin();
109  while (itr != elems.end()) {
110  if ((*itr)->type() == reco::PFBlockElement::TRACK) {
111  const reco::PFBlockElementTrack* trkel = static_cast<reco::PFBlockElementTrack*>(itr->get());
112  const reco::ConversionRefVector& cRef = trkel->convRefs();
114  const reco::VertexCompositeCandidateRef& v0Ref = trkel->V0Ref();
115  // if there is no displaced vertex reference and it is marked
116  // as a conversion it's gotta be a converted brem
117  if (trkel->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) && cRef.empty() && dvRef.isNull() &&
118  v0Ref.isNull()) {
119  // if the Pt resolution is bad we kill this element
122  itr = elems.erase(itr);
123  continue;
124  }
125  }
126  }
127  ++itr;
128  } // loop on existing elements
129  }
130  // preprocess existing tracks in the element list and create a mask
131  // so that we do not import tracks twice, tag muons we find
132  // in this collection
133  auto TKs_end = std::partition(
134  elems.begin(), elems.end(), [](const ElementType& a) { return a->type() == reco::PFBlockElement::TRACK; });
135  auto btk_elems = elems.begin();
136  auto btrack = tracks->cbegin();
137  auto etrack = tracks->cend();
138  for (auto track = btrack; track != etrack; ++track) {
139  auto tk_elem =
140  std::find_if(btk_elems, TKs_end, [&](const ElementType& a) { return (a->trackRef() == track->trackRef()); });
141  if (tk_elem != TKs_end) {
142  mask[std::distance(tracks->cbegin(), track)] = false;
143  // check and update if this track is a muon
144  const int muId = PFMuonAlgo::muAssocToTrack((*tk_elem)->trackRef(), muons);
145  if (muId != -1) {
146  muonref = reco::MuonRef(muonH, muId);
147  if (PFMuonAlgo::isLooseMuon(muonref) || PFMuonAlgo::isMuon(muonref)) {
148  static_cast<reco::PFBlockElementTrack*>(tk_elem->get())->setMuonRef(muonref);
149  }
150  }
151  }
152  }
153  // now we actually insert tracks, again tagging muons along the way
154  reco::PFRecTrackRef pftrackref;
155  reco::PFBlockElementTrack* trkElem = nullptr;
156  for (auto track = btrack; track != etrack; ++track) {
157  const unsigned idx = std::distance(btrack, track);
158  // since we already set muon refs in the previously imported tracks,
159  // here we can skip everything that is already imported
160  if (!mask[idx])
161  continue;
162  muonref = reco::MuonRef();
163  pftrackref = reco::PFRecTrackRef(tracks, idx);
164  // Get the eventual muon associated to this track
165  const int muId = PFMuonAlgo::muAssocToTrack(pftrackref->trackRef(), muons);
166  bool thisIsAPotentialMuon = false;
167  if (muId != -1) {
168  muonref = reco::MuonRef(muonH, muId);
169  thisIsAPotentialMuon =
170  ((PFMuonAlgo::hasValidTrack(muonref, true, muonMaxDPtOPt_) && PFMuonAlgo::isLooseMuon(muonref)) ||
171  (PFMuonAlgo::hasValidTrack(muonref, false, muonMaxDPtOPt_) && PFMuonAlgo::isMuon(muonref)));
172  }
173  if (thisIsAPotentialMuon || PFTrackAlgoTools::goodPtResolution(
174  pftrackref->trackRef(), DPtovPtCut_, NHitCut_, useIterTracking_, trackQuality_)) {
175  trkElem = new reco::PFBlockElementTrack(pftrackref);
176  if (thisIsAPotentialMuon) {
177  LogDebug("GeneralTracksImporter")
178  << "Potential Muon P " << pftrackref->trackRef()->p() << " pt " << pftrackref->trackRef()->p() << std::endl;
179  }
180  if (muId != -1)
181  trkElem->setMuonRef(muonref);
182  //
183  // if _vetoEndcap is false, add this trk automatically.
184  // if _vetoEndcap is true, veto against the hgcal region tracks or charged PF candidates.
185  // when simPF is used, we don't veto tracks with muonref even if they are in the hgcal region.
186  //
187  if (!vetoEndcap_)
188  elems.emplace_back(trkElem);
189  else {
190  TrackProdIDKey trk = std::make_pair(pftrackref->trackRef().id(), pftrackref->trackRef().key());
191  auto lower = std::lower_bound(vetoed.begin(), vetoed.end(), trk);
192  bool inVetoList = (lower != vetoed.end() && *lower == trk);
193  if (!inVetoList || (vetoMode_ == pfRecTrackCollection && muonref.isNonnull())) {
194  elems.emplace_back(trkElem);
195  } else
196  delete trkElem;
197  }
198  }
199  }
200  elems.shrink_to_fit();
201 }
GeneralTracksImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
static int muAssocToTrack(const reco::TrackRef &trackref, const reco::MuonCollection &muons)
Definition: PFMuonAlgo.cc:479
const ConversionRefVector & convRefs() const override
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const PFDisplacedTrackerVertexRef & displacedVertexRef(TrackType trType) const override
bool goodPtResolution(const reco::TrackRef &, const std::vector< double > &DPtovPtCut, const std::vector< unsigned > &NHitCut, bool useIterTracking, const reco::TrackBase::TrackQuality trackQuality)
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
const std::vector< double > DPtovPtCut_
Quality qualityByName(std::string const &name)
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:48
TrackQuality
track quality
Definition: TrackBase.h:150
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::TrackCollection > tracksSrc_
edm::EDGetTokenT< std::vector< TICLSeedingRegion > > vetoTICLSeedingSrc_
assert(be >=bs)
key_type key() const
Accessor for product key.
Definition: Ref.h:250
constexpr uint32_t mask
Definition: gpuClustering.h:24
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
const VertexCompositeCandidateRef & V0Ref() const override
void importToBlock(const edm::Event &, ElementList &) const override
const edm::EDGetTokenT< reco::MuonCollection > muons_
bool trackType(TrackType trType) const override
edm::Ref< PFRecTrackCollection > PFRecTrackRef
persistent reference to PFRecTrack objects
Definition: PFRecTrackFwd.h:15
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
const std::vector< unsigned > NHitCut_
bool isNull() const
Checks for null.
Definition: Ref.h:235
const reco::TrackBase::TrackQuality trackQuality_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
auto const & tracks
cannot be loose
edm::EDGetTokenT< reco::PFRecTrackCollection > vetoPFTracksSrc_
edm::EDGetTokenT< reco::PFCandidateCollection > vetoPFCandidatesSrc_
static bool isLooseMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:57
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
const reco::TrackRef & trackRef() const override
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:119
void setMuonRef(const MuonRef &muref) override
reference to the Muon
static bool hasValidTrack(const reco::MuonRef &muonRef, bool loose, double maxDPtOPt)
Definition: PFMuonAlgo.cc:348
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
const edm::EDGetTokenT< reco::PFRecTrackCollection > src_
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
#define LogDebug(id)