CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackFromParentImporter.h
Go to the documentation of this file.
1 #ifndef __TrackFromParentImporter_H__
2 #define __TrackFromParentImporter_H__
3 
9 
10 namespace pflow {
11  namespace noop {
12  // this adaptor class gets redefined later to match the
13  // needs of the collection and importing cuts we are using
14  template <class Collection>
16  public:
17  static bool check_importable(const typename Collection::value_type&) { return true; }
18  static const std::vector<reco::PFRecTrackRef>& get_track_refs(const typename Collection::value_type&) {
19  return empty_;
20  }
22  static const std::vector<reco::PFRecTrackRef> empty_;
23  };
24  } // namespace noop
25  namespace importers {
26  template <class Collection, class Adaptor = noop::ParentCollectionAdaptor<Collection>>
28  public:
30  : BlockElementImporterBase(conf, cc),
31  src_(cc.consumes<Collection>(conf.getParameter<edm::InputTag>("source"))),
32  vetoEndcap_(conf.getParameter<bool>("vetoEndcap")) {
33  if (vetoEndcap_) {
34  vetoMode_ = conf.getParameter<unsigned>("vetoMode");
35  switch (vetoMode_) {
38  break;
39  case ticlSeedingRegion:
41  cc.consumes<std::vector<TICLSeedingRegion>>(conf.getParameter<edm::InputTag>("vetoSrc"));
43  break;
47  break;
48  } // switch
49  } // vetoEndcap_
50  }
51 
52  void importToBlock(const edm::Event&, ElementList&) const override;
53 
54  private:
57  const bool vetoEndcap_;
58  unsigned int vetoMode_;
63  };
64 
65  template <class Collection, class Adaptor>
67  const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
69  auto pfparents = e.getHandle(src_);
70  //
71  // Store tracks to be vetoed
72  typedef std::pair<edm::ProductID, unsigned> TrackProdIDKey;
73  std::vector<TrackProdIDKey> vetoed;
74  edm::ProductID prodIdForVeto;
75  if (vetoEndcap_) {
76  switch (vetoMode_) {
77  case pfRecTrackCollection: {
78  const auto& vetoes = e.get(vetoPFTracksSrc_);
79  for (const auto& veto : vetoes) {
80  vetoed.emplace_back(veto.trackRef().id(), veto.trackRef().key());
81  }
82  break;
83  }
84  case ticlSeedingRegion: {
85  const auto& vetoes = e.get(vetoTICLSeedingSrc_);
86  auto tracksH = e.getHandle(tracksSrc_);
87  for (const auto& veto : vetoes) {
88  assert(veto.collectionID == tracksH.id());
89  vetoed.emplace_back(tracksH.id(), veto.index); // track prod id and key
90  }
91  break;
92  }
93  case pfCandidateCollection: {
94  const auto& vetoes = e.get(vetoPFCandidatesSrc_);
95  for (const auto& veto : vetoes) {
96  if (veto.trackRef().isNull())
97  continue;
98  vetoed.emplace_back(veto.trackRef().id(), veto.trackRef().key());
99  }
100  break;
101  }
102  } // switch
103  std::sort(vetoed.begin(), vetoed.end());
104  }
105  //
106  elems.reserve(elems.size() + 2 * pfparents->size());
107  //
108  auto TKs_end = std::partition(
109  elems.begin(), elems.end(), [](const ElementType& a) { return a->type() == reco::PFBlockElement::TRACK; });
110  // insert tracks into the element list, updating tracks that exist already
111  auto bpar = pfparents->cbegin();
112  auto epar = pfparents->cend();
113  edm::Ref<Collection> parentRef;
114  reco::PFBlockElement* trkElem = nullptr;
115  for (auto pfparent = bpar; pfparent != epar; ++pfparent) {
116  if (Adaptor::check_importable(*pfparent)) {
117  parentRef = edm::Ref<Collection>(pfparents, std::distance(bpar, pfparent));
118  const auto& pftracks = Adaptor::get_track_refs(*pfparent);
119  for (const auto& pftrack : pftracks) {
120  if (vetoEndcap_) { // vetoEndcap flag
121  TrackProdIDKey trk = std::make_pair(pftrack->trackRef().id(), pftrack->trackRef().key());
122  auto lower = std::lower_bound(vetoed.begin(), vetoed.end(), trk);
123  bool inVetoList = (lower != vetoed.end() && *lower == trk);
124  if (inVetoList)
125  continue; // found a track in a veto list
126  }
127  //
128  // Now try to update an entry in pfblock or import
129  auto tk_elem = std::find_if(
130  elems.begin(), TKs_end, [&](const ElementType& a) { return (a->trackRef() == pftrack->trackRef()); });
131  if (tk_elem != TKs_end) { // if found flag the track, otherwise import
132  Adaptor::set_element_info(tk_elem->get(), parentRef);
133  } else {
134  trkElem = new reco::PFBlockElementTrack(pftrack);
135  Adaptor::set_element_info(trkElem, parentRef);
136  TKs_end = elems.insert(TKs_end, ElementType(trkElem));
137  ++TKs_end;
138  }
139  } // daughter track loop ends
140  } // end of importable check
141  } // loop on tracking coming from common parent
142  elems.shrink_to_fit();
143  } // end of importToBlock
144  } // namespace importers
145 } // namespace pflow
146 #endif
static bool check_importable(const typename Collection::value_type &)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Abstract base class for a PFBlock element (track, cluster...)
static const std::vector< reco::PFRecTrackRef > empty_
edm::EDGetTokenT< reco::PFCandidateCollection > vetoPFCandidatesSrc_
edm::EDGetTokenT< reco::PFRecTrackCollection > vetoPFTracksSrc_
static void set_element_info(reco::PFBlockElement *, const typename edm::Ref< Collection > &)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::TrackCollection > tracksSrc_
assert(be >=bs)
tuple pfCandidateCollection
TrackFromParentImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
static const std::vector< reco::PFRecTrackRef > & get_track_refs(const typename Collection::value_type &)
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void importToBlock(const edm::Event &, ElementList &) const override
__host__ __device__ constexpr RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
double a
Definition: hdecay.h:119
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
edm::EDGetTokenT< std::vector< TICLSeedingRegion > > vetoTICLSeedingSrc_