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>
37  auto clusters = e.getHandle(_src);
38  auto assoc = e.getHandle(_assoc);
39  auto bclus = clusters->cbegin();
40  auto eclus = clusters->cend();
41  // get all the SCs in the element list
42  auto sc_end = std::partition(elems.begin(),elems.end(),
43  [](const ElementList::value_type& o){
44  return o->type() == reco::PFBlockElement::SC;
45  });
46  ecals.reserve(clusters->size());
47  for( auto clus = bclus; clus != eclus; ++clus ) {
48  reco::PFClusterRef tempref(clusters, std::distance(bclus,clus));
49  reco::PFBlockElementCluster* newelem =
50  new reco::PFBlockElementCluster(tempref,T);
51  for( auto scelem = elems.begin(); scelem != sc_end; ++scelem ) {
52  const reco::PFBlockElementSuperCluster* elem_as_sc =
53  static_cast<const reco::PFBlockElementSuperCluster*>(scelem->get());
54  const reco::SuperClusterRef& this_sc = elem_as_sc->superClusterRef();
55 
56  const bool in_sc = ( elem_as_sc->fromPFSuperCluster() ?
57  // use association map if from PFSC
59  *this_sc,
60  *assoc) :
61  // match by overlapping rechit otherwise
63  *this_sc) );
64  if( in_sc ) {
65  newelem->setSuperClusterRef(this_sc);
66  break;
67  }
68  }
69  ecals.emplace_back(newelem);
70  }
71  elems.reserve(elems.size()+ecals.size());
72  for( auto& ecal : ecals ) {
73  elems.emplace_back(ecal.release());
74  }
75 }
76 
77 
81  "ECALClusterImporter");
82 
86  "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
edm::Ptr< CaloCluster > CaloClusterPtr
edm::EDGetTokenT< reco::PFClusterCollection > _src
SpecialClusterImporter< reco::PFBlockElement::HGCAL > HGCalClusterImporter
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:539
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