CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  edm::Handle<Collection> pfparents;
52  e.getByToken(_src,pfparents);
53  elems.reserve(elems.size() + 2*pfparents->size());
54  // setup our elements so that all the SCs are grouped together
55  auto TKs_end = std::partition(elems.begin(),elems.end(),
56  [](const ElementType& a){
57  return a->type() == reco::PFBlockElement::TRACK;
58  });
59  // insert tracks into the element list, updating tracks that exist already
60  auto bpar = pfparents->cbegin();
61  auto epar = pfparents->cend();
62  edm::Ref<Collection> parentRef;
63  reco::PFBlockElement* trkElem = NULL;
64  for( auto pfparent = bpar; pfparent != epar; ++pfparent ) {
65  if( Adaptor::check_importable(*pfparent) ) {
66  parentRef = edm::Ref<Collection>(pfparents,std::distance(bpar,pfparent));
67  const auto& pftracks = Adaptor::get_track_refs(*pfparent);
68  for( const auto& pftrack : pftracks ) {
69  auto tk_elem = std::find_if(elems.begin(),TKs_end,
70  [&](const ElementType& a){
71  return ( a->trackRef() ==
72  pftrack->trackRef() );
73  });
74  if( tk_elem != TKs_end ) { // if found flag the track, otherwise import
75  Adaptor::set_element_info(tk_elem->get(),parentRef);
76  } else {
77  trkElem = new reco::PFBlockElementTrack(pftrack);
78  Adaptor::set_element_info(trkElem,parentRef);
79  TKs_end = elems.insert(TKs_end,ElementType(trkElem));
80  ++TKs_end;
81  }
82  }
83  }
84  }// loop on tracking coming from common parent
85  elems.shrink_to_fit();
86  }
87  } // importers
88 } // pflow
89 #endif
static bool check_importable(const typename Collection::value_type &)
Abstract base class for a PFBlock element (track, cluster...)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
static void set_element_info(reco::PFBlockElement *, const typename edm::Ref< Collection > &)
#define NULL
Definition: scimark2.h:8
static const std::vector< reco::PFRecTrackRef > & get_track_refs(const typename Collection::value_type &)
static const std::vector< reco::PFRecTrackRef > _empty
void importToBlock(const edm::Event &, ElementList &) const override
TrackFromParentImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
double a
Definition: hdecay.h:121
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList