CMS 3D CMS Logo

SpecialClusterImporter.cc
Go to the documentation of this file.
10 
11 // NOTE! This should come *after* and importers that bring in super clusters
12 // of their own (like electron seeds or photons)
13 // otherwise ECAL <-> ECAL linking will not work correctly
14 template<reco::PFBlockElement::Type T>
16 public:
18  edm::ConsumesCollector& sumes) :
19  BlockElementImporterBase(conf,sumes),
20  _src(sumes.consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("source"))),
21  _assoc(sumes.consumes<edm::ValueMap<reco::CaloClusterPtr> >(conf.getParameter<edm::InputTag>("BCtoPFCMap"))) {}
22 
23  void importToBlock( const edm::Event& ,
24  ElementList& ) const override;
25 
26 private:
29 };
30 
31 
32 template<reco::PFBlockElement::Type T>
39  e.getByToken(_src,clusters);
40  e.getByToken(_assoc,assoc);
41  auto bclus = clusters->cbegin();
42  auto eclus = clusters->cend();
43  // get all the SCs in the element list
44  auto sc_end = std::partition(elems.begin(),elems.end(),
45  [](const ElementList::value_type& o){
46  return o->type() == reco::PFBlockElement::SC;
47  });
48  ecals.reserve(clusters->size());
49  for( auto clus = bclus; clus != eclus; ++clus ) {
50  reco::PFClusterRef tempref(clusters, std::distance(bclus,clus));
51  reco::PFBlockElementCluster* newelem =
52  new reco::PFBlockElementCluster(tempref,T);
53  for( auto scelem = elems.begin(); scelem != sc_end; ++scelem ) {
54  const reco::PFBlockElementSuperCluster* elem_as_sc =
55  static_cast<const reco::PFBlockElementSuperCluster*>(scelem->get());
56  const reco::SuperClusterRef& this_sc = elem_as_sc->superClusterRef();
57 
58  const bool in_sc = ( elem_as_sc->fromPFSuperCluster() ?
59  // use association map if from PFSC
61  *this_sc,
62  *assoc) :
63  // match by overlapping rechit otherwise
65  *this_sc) );
66  if( in_sc ) {
67  newelem->setSuperClusterRef(this_sc);
68  break;
69  }
70  }
71  ecals.emplace_back(newelem);
72  }
73  elems.reserve(elems.size()+ecals.size());
74  for( auto& ecal : ecals ) {
75  elems.emplace_back(ecal.release());
76  }
77 }
78 
79 
83  "ECALClusterImporter");
84 
88  "HGCalClusterImporter");
const SuperClusterRef & superClusterRef() const
static bool overlap(const reco::CaloCluster &sc1, const reco::CaloCluster &sc, float minfrac=0.01, bool debug=false)
edm::EDGetTokenT< edm::ValueMap< reco::CaloClusterPtr > > _assoc
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
edm::Ptr< CaloCluster > CaloClusterPtr
edm::EDGetTokenT< reco::PFClusterCollection > _src
SpecialClusterImporter< reco::PFBlockElement::HGCAL > HGCalClusterImporter
SpecialClusterImporter< reco::PFBlockElement::ECAL > ECALClusterImporter
fixed size matrix
HLT enums.
SpecialClusterImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
void importToBlock(const edm::Event &, ElementList &) const override
void setSuperClusterRef(const SuperClusterRef &ref)
long double T