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 "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
0004 
0005 class GenericClusterImporter : public BlockElementImporterBase {
0006 public:
0007   GenericClusterImporter(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0008       : BlockElementImporterBase(conf, cc),
0009         _src(cc.consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("source"))) {}
0010 
0011   void importToBlock(const edm::Event&, ElementList&) const override;
0012 
0013 private:
0014   edm::EDGetTokenT<reco::PFClusterCollection> _src;
0015 };
0016 
0017 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, GenericClusterImporter, "GenericClusterImporter");
0018 
0019 void GenericClusterImporter::importToBlock(const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
0020   auto clusters = e.getHandle(_src);
0021   auto cbegin = clusters->cbegin();
0022   auto cend = clusters->cend();
0023   for (auto clus = cbegin; clus != cend; ++clus) {
0024     reco::PFBlockElement::Type type = reco::PFBlockElement::NONE;
0025     reco::PFClusterRef cref(clusters, std::distance(cbegin, clus));
0026     switch (clus->layer()) {
0027       case PFLayer::PS1:
0028         type = reco::PFBlockElement::PS1;
0029         break;
0030       case PFLayer::PS2:
0031         type = reco::PFBlockElement::PS2;
0032         break;
0033       case PFLayer::ECAL_BARREL:
0034       case PFLayer::ECAL_ENDCAP:
0035         type = reco::PFBlockElement::ECAL;
0036         break;
0037       case PFLayer::HCAL_BARREL1:
0038       case PFLayer::HCAL_ENDCAP:
0039         type = reco::PFBlockElement::HCAL;
0040         break;
0041       case PFLayer::HCAL_BARREL2:
0042         type = reco::PFBlockElement::HO;
0043         break;
0044       case PFLayer::HF_EM:
0045         type = reco::PFBlockElement::HFEM;
0046         break;
0047       case PFLayer::HF_HAD:
0048         type = reco::PFBlockElement::HFHAD;
0049         break;
0050       case PFLayer::HGCAL:
0051         type = reco::PFBlockElement::HGCAL;
0052         break;
0053       default:
0054         throw cms::Exception("InvalidPFLayer") << "Layer given, " << clus->layer() << " is not a valid PFLayer!";
0055     }
0056     reco::PFBlockElement* cptr = new reco::PFBlockElementCluster(cref, type);
0057     elems.emplace_back(cptr);
0058   }
0059 }