CMS 3D CMS Logo

TrackFromParentImporter.h
Go to the documentation of this file.
1 #ifndef __TrackFromParentImporter_H__
2 #define __TrackFromParentImporter_H__
3 
8 
9 namespace pflow {
10  namespace noop {
11  // this adaptor class gets redefined later to match the
12  // needs of the collection and importing cuts we are using
13  template<class Collection>
15  public:
16  static bool check_importable(const typename Collection::value_type&) {
17  return true;
18  }
19  static const std::vector<reco::PFRecTrackRef>&
21  return _empty;
22  }
24  const typename edm::Ref<Collection>&) {
25  }
26  static const std::vector<reco::PFRecTrackRef> _empty;
27  };
28  }
29  namespace importers {
30  template<class Collection,class Adaptor=noop::ParentCollectionAdaptor<Collection> >
32  public:
34  edm::ConsumesCollector& sumes) :
35  BlockElementImporterBase(conf,sumes),
36  _src(sumes.consumes<Collection>(conf.getParameter<edm::InputTag>("source"))) {}
37 
38  void importToBlock( const edm::Event& ,
39  ElementList& ) const override;
40 
41  private:
43  };
44 
45 
46  template<class Collection, class Adaptor>
51  auto pfparents = e.getHandle(_src);
52  elems.reserve(elems.size() + 2*pfparents->size());
53  // setup our elements so that all the SCs are grouped together
54  auto TKs_end = std::partition(elems.begin(),elems.end(),
55  [](const ElementType& a){
56  return a->type() == reco::PFBlockElement::TRACK;
57  });
58  // insert tracks into the element list, updating tracks that exist already
59  auto bpar = pfparents->cbegin();
60  auto epar = pfparents->cend();
61  edm::Ref<Collection> parentRef;
62  reco::PFBlockElement* trkElem = nullptr;
63  for( auto pfparent = bpar; pfparent != epar; ++pfparent ) {
64  if( Adaptor::check_importable(*pfparent) ) {
65  parentRef = edm::Ref<Collection>(pfparents,std::distance(bpar,pfparent));
66  const auto& pftracks = Adaptor::get_track_refs(*pfparent);
67  for( const auto& pftrack : pftracks ) {
68  auto tk_elem = std::find_if(elems.begin(),TKs_end,
69  [&](const ElementType& a){
70  return ( a->trackRef() ==
71  pftrack->trackRef() );
72  });
73  if( tk_elem != TKs_end ) { // if found flag the track, otherwise import
74  Adaptor::set_element_info(tk_elem->get(),parentRef);
75  } else {
76  trkElem = new reco::PFBlockElementTrack(pftrack);
77  Adaptor::set_element_info(trkElem,parentRef);
78  TKs_end = elems.insert(TKs_end,ElementType(trkElem));
79  ++TKs_end;
80  }
81  }
82  }
83  }// loop on tracking coming from common parent
84  elems.shrink_to_fit();
85  }
86  } // importers
87 } // pflow
88 #endif
static bool check_importable(const typename Collection::value_type &)
Abstract base class for a PFBlock element (track, cluster...)
static void set_element_info(reco::PFBlockElement *, const typename edm::Ref< Collection > &)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:539
static const std::vector< reco::PFRecTrackRef > & get_track_refs(const typename Collection::value_type &)
static const std::vector< reco::PFRecTrackRef > _empty
HLT enums.
TrackFromParentImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
double a
Definition: hdecay.h:121
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList