CMS 3D CMS Logo

SpecialClusterImporter.cc
Go to the documentation of this file.
8 
9 // NOTE! This should come *after* and importers that bring in super clusters
10 // of their own (like electron seeds or photons)
11 // otherwise ECAL <-> ECAL linking will not work correctly
12 template <reco::PFBlockElement::Type T>
14 public:
17  _src(cc.consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("source"))),
18  _assoc(cc.consumes<edm::ValueMap<reco::CaloClusterPtr> >(conf.getParameter<edm::InputTag>("BCtoPFCMap"))) {}
19 
20  void importToBlock(const edm::Event&, ElementList&) const override;
21 
22 private:
25 };
26 
27 template <reco::PFBlockElement::Type T>
30  auto clusters = e.getHandle(_src);
31  auto assoc = e.getHandle(_assoc);
32  auto bclus = clusters->cbegin();
33  auto eclus = clusters->cend();
34  // get all the SCs in the element list
35  auto sc_end = std::partition(elems.begin(), elems.end(), [](const ElementList::value_type& o) {
36  return o->type() == reco::PFBlockElement::SC;
37  });
38  ecals.reserve(clusters->size());
39  for (auto clus = bclus; clus != eclus; ++clus) {
40  reco::PFClusterRef tempref(clusters, std::distance(bclus, clus));
42  for (auto scelem = elems.begin(); scelem != sc_end; ++scelem) {
43  const reco::PFBlockElementSuperCluster* elem_as_sc =
44  static_cast<const reco::PFBlockElementSuperCluster*>(scelem->get());
45  const reco::SuperClusterRef& this_sc = elem_as_sc->superClusterRef();
46 
47  const bool in_sc = (elem_as_sc->fromPFSuperCluster() ?
48  // use association map if from PFSC
49  ClusterClusterMapping::overlap(tempref, *this_sc, *assoc)
50  :
51  // match by overlapping rechit otherwise
52  ClusterClusterMapping::overlap(*tempref, *this_sc));
53  if (in_sc) {
54  newelem->setSuperClusterRef(this_sc);
55  break;
56  }
57  }
58  ecals.emplace_back(newelem);
59  }
60  elems.reserve(elems.size() + ecals.size());
61  for (auto& ecal : ecals) {
62  elems.emplace_back(ecal.release());
63  }
64 }
65 
68 
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
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
const SuperClusterRef & superClusterRef() const
edm::Ptr< CaloCluster > CaloClusterPtr
edm::EDGetTokenT< reco::PFClusterCollection > _src
SpecialClusterImporter< reco::PFBlockElement::HGCAL > HGCalClusterImporter
void importToBlock(const edm::Event &, ElementList &) const override
SpecialClusterImporter< reco::PFBlockElement::ECAL > ECALClusterImporter
std::vector< l1t::PFCluster > PFClusterCollection
Definition: PFCluster.h:87
fixed size matrix
HLT enums.
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
SpecialClusterImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
void setSuperClusterRef(const SuperClusterRef &ref)
long double T