Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:27

0001 #include "RecoParticleFlow/PFProducer/interface/BlockElementImporterBase.h"
0002 #include "RecoParticleFlow/PFClusterTools/interface/ClusterClusterMapping.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementSuperCluster.h"
0006 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0007 #include "DataFormats/Common/interface/ValueMap.h"
0008 
0009 // NOTE! This should come *after* and importers that bring in super clusters
0010 // of their own (like electron seeds or photons)
0011 // otherwise ECAL <-> ECAL linking will not work correctly
0012 template <reco::PFBlockElement::Type T>
0013 class SpecialClusterImporter : public BlockElementImporterBase {
0014 public:
0015   SpecialClusterImporter(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0016       : BlockElementImporterBase(conf, cc),
0017         _src(cc.consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("source"))),
0018         _assoc(cc.consumes<edm::ValueMap<reco::CaloClusterPtr> >(conf.getParameter<edm::InputTag>("BCtoPFCMap"))) {}
0019 
0020   void importToBlock(const edm::Event&, ElementList&) const override;
0021 
0022 private:
0023   edm::EDGetTokenT<reco::PFClusterCollection> _src;
0024   edm::EDGetTokenT<edm::ValueMap<reco::CaloClusterPtr> > _assoc;
0025 };
0026 
0027 template <reco::PFBlockElement::Type T>
0028 void SpecialClusterImporter<T>::importToBlock(const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
0029   BlockElementImporterBase::ElementList ecals;
0030   auto clusters = e.getHandle(_src);
0031   auto assoc = e.getHandle(_assoc);
0032   auto bclus = clusters->cbegin();
0033   auto eclus = clusters->cend();
0034   // get all the SCs in the element list
0035   auto sc_end = std::partition(elems.begin(), elems.end(), [](const ElementList::value_type& o) {
0036     return o->type() == reco::PFBlockElement::SC;
0037   });
0038   ecals.reserve(clusters->size());
0039   for (auto clus = bclus; clus != eclus; ++clus) {
0040     reco::PFClusterRef tempref(clusters, std::distance(bclus, clus));
0041     reco::PFBlockElementCluster* newelem = new reco::PFBlockElementCluster(tempref, T);
0042     for (auto scelem = elems.begin(); scelem != sc_end; ++scelem) {
0043       const reco::PFBlockElementSuperCluster* elem_as_sc =
0044           static_cast<const reco::PFBlockElementSuperCluster*>(scelem->get());
0045       const reco::SuperClusterRef& this_sc = elem_as_sc->superClusterRef();
0046 
0047       const bool in_sc = (elem_as_sc->fromPFSuperCluster() ?
0048                                                            // use association map if from PFSC
0049                               ClusterClusterMapping::overlap(tempref, *this_sc, *assoc)
0050                                                            :
0051                                                            // match by overlapping rechit otherwise
0052                               ClusterClusterMapping::overlap(*tempref, *this_sc));
0053       if (in_sc) {
0054         newelem->setSuperClusterRef(this_sc);
0055         break;
0056       }
0057     }
0058     ecals.emplace_back(newelem);
0059   }
0060   elems.reserve(elems.size() + ecals.size());
0061   for (auto& ecal : ecals) {
0062     elems.emplace_back(ecal.release());
0063   }
0064 }
0065 
0066 typedef SpecialClusterImporter<reco::PFBlockElement::ECAL> ECALClusterImporter;
0067 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, ECALClusterImporter, "ECALClusterImporter");
0068 
0069 typedef SpecialClusterImporter<reco::PFBlockElement::HGCAL> HGCalClusterImporter;
0070 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, HGCalClusterImporter, "HGCalClusterImporter");