File indexing completed on 2023-10-25 10:01:43
0001 #include "RecoParticleFlow/PFProducer/interface/BlockElementImporterBase.h"
0002 #include "RecoParticleFlow/PFProducer/interface/PhotonSelectorAlgo.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementSuperCluster.h"
0005 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0006 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0007 #include "RecoParticleFlow/PFProducer/interface/PFBlockElementSCEqual.h"
0008
0009 #include <memory>
0010
0011 #include <unordered_map>
0012
0013 class EGPhotonImporter : public BlockElementImporterBase {
0014 public:
0015 enum SelectionChoices { SeparateDetectorIso, CombinedDetectorIso };
0016
0017 EGPhotonImporter(const edm::ParameterSet&, edm::ConsumesCollector&);
0018
0019 void importToBlock(const edm::Event&, ElementList&) const override;
0020
0021 private:
0022 edm::EDGetTokenT<reco::PhotonCollection> _src;
0023 const std::unordered_map<std::string, SelectionChoices> _selectionTypes;
0024 SelectionChoices _selectionChoice;
0025 std::unique_ptr<const PhotonSelectorAlgo> _selector;
0026 bool _superClustersArePF;
0027 };
0028
0029 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, EGPhotonImporter, "EGPhotonImporter");
0030
0031 EGPhotonImporter::EGPhotonImporter(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0032 : BlockElementImporterBase(conf, cc),
0033 _src(cc.consumes<reco::PhotonCollection>(conf.getParameter<edm::InputTag>("source"))),
0034 _selectionTypes({{"SeparateDetectorIso", EGPhotonImporter::SeparateDetectorIso},
0035 {"CombinedDetectorIso", EGPhotonImporter::CombinedDetectorIso}}),
0036 _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {
0037 const std::string& selChoice = conf.getParameter<std::string>("SelectionChoice");
0038 _selectionChoice = _selectionTypes.at(selChoice);
0039 const edm::ParameterSet& selDef = conf.getParameterSet("SelectionDefinition");
0040 const float minEt = selDef.getParameter<double>("minEt");
0041 const float trackIso_const = selDef.getParameter<double>("trackIsoConstTerm");
0042 const float trackIso_slope = selDef.getParameter<double>("trackIsoSlopeTerm");
0043 const float ecalIso_const = selDef.getParameter<double>("ecalIsoConstTerm");
0044 const float ecalIso_slope = selDef.getParameter<double>("ecalIsoSlopeTerm");
0045 const float hcalIso_const = selDef.getParameter<double>("hcalIsoConstTerm");
0046 const float hcalIso_slope = selDef.getParameter<double>("hcalIsoSlopeTerm");
0047 const float hoe = selDef.getParameter<double>("HoverE");
0048 const float loose_hoe = selDef.getParameter<double>("LooseHoverE");
0049 const float combIso = selDef.getParameter<double>("combIsoConstTerm");
0050 _selector = std::make_unique<PhotonSelectorAlgo>((float)_selectionChoice,
0051 minEt,
0052 trackIso_const,
0053 trackIso_slope,
0054 ecalIso_const,
0055 ecalIso_slope,
0056 hcalIso_const,
0057 hcalIso_slope,
0058 hoe,
0059 combIso,
0060 loose_hoe);
0061 }
0062
0063 void EGPhotonImporter::importToBlock(const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
0064 typedef BlockElementImporterBase::ElementList::value_type ElementType;
0065 auto photons = e.getHandle(_src);
0066 elems.reserve(elems.size() + photons->size());
0067
0068 auto SCs_end = std::partition(
0069 elems.begin(), elems.end(), [](const ElementType& a) { return a->type() == reco::PFBlockElement::SC; });
0070
0071 auto bphoton = photons->cbegin();
0072 auto ephoton = photons->cend();
0073 reco::PFBlockElementSuperCluster* scbe = nullptr;
0074 reco::PhotonRef phoref;
0075 for (auto photon = bphoton; photon != ephoton; ++photon) {
0076 if (_selector->passPhotonSelection(*photon)) {
0077 phoref = reco::PhotonRef(photons, std::distance(bphoton, photon));
0078 const reco::SuperClusterRef& scref = photon->superCluster();
0079 PFBlockElementSCEqual myEqual(scref);
0080 auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
0081 if (sc_elem != SCs_end) {
0082 scbe = static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
0083 scbe->setFromPhoton(true);
0084 scbe->setPhotonRef(phoref);
0085 scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
0086 scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
0087 scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
0088 scbe->setHoE(photon->hadronicOverEm());
0089 } else {
0090 scbe = new reco::PFBlockElementSuperCluster(scref);
0091 scbe->setFromPhoton(true);
0092 scbe->setFromPFSuperCluster(_superClustersArePF);
0093 scbe->setPhotonRef(phoref);
0094 scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
0095 scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
0096 scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
0097 scbe->setHoE(photon->hadronicOverEm());
0098 SCs_end = elems.insert(SCs_end, ElementType(scbe));
0099 ++SCs_end;
0100 }
0101 }
0102 }
0103 elems.shrink_to_fit();
0104 }