CMS 3D CMS Logo

EGPhotonImporter.cc
Go to the documentation of this file.
8 
9 #include <memory>
10 
11 #include <unordered_map>
12 
14 public:
16 
18 
19  void importToBlock(const edm::Event&, ElementList&) const override;
20 
21 private:
23  const std::unordered_map<std::string, SelectionChoices> _selectionTypes;
25  std::unique_ptr<const PhotonSelectorAlgo> _selector;
27 };
28 
30 
33  _src(cc.consumes<reco::PhotonCollection>(conf.getParameter<edm::InputTag>("source"))),
34  _selectionTypes({{"SeparateDetectorIso", EGPhotonImporter::SeparateDetectorIso},
35  {"CombinedDetectorIso", EGPhotonImporter::CombinedDetectorIso}}),
36  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {
37  const std::string& selChoice = conf.getParameter<std::string>("SelectionChoice");
38  _selectionChoice = _selectionTypes.at(selChoice);
39  const edm::ParameterSet& selDef = conf.getParameterSet("SelectionDefinition");
40  const float minEt = selDef.getParameter<double>("minEt");
41  const float trackIso_const = selDef.getParameter<double>("trackIsoConstTerm");
42  const float trackIso_slope = selDef.getParameter<double>("trackIsoSlopeTerm");
43  const float ecalIso_const = selDef.getParameter<double>("ecalIsoConstTerm");
44  const float ecalIso_slope = selDef.getParameter<double>("ecalIsoSlopeTerm");
45  const float hcalIso_const = selDef.getParameter<double>("hcalIsoConstTerm");
46  const float hcalIso_slope = selDef.getParameter<double>("hcalIsoSlopeTerm");
47  const float hoe = selDef.getParameter<double>("HoverE");
48  const float loose_hoe = selDef.getParameter<double>("LooseHoverE");
49  const float combIso = selDef.getParameter<double>("combIsoConstTerm");
50  _selector = std::make_unique<PhotonSelectorAlgo>((float)_selectionChoice,
51  minEt,
52  trackIso_const,
53  trackIso_slope,
54  ecalIso_const,
55  ecalIso_slope,
56  hcalIso_const,
57  hcalIso_slope,
58  hoe,
59  combIso,
60  loose_hoe);
61 }
62 
65  auto photons = e.getHandle(_src);
66  elems.reserve(elems.size() + photons->size());
67  // setup our elements so that all the SCs are grouped together
68  auto SCs_end = std::partition(
69  elems.begin(), elems.end(), [](const ElementType& a) { return a->type() == reco::PFBlockElement::SC; });
70  //now add the photons
71  auto bphoton = photons->cbegin();
72  auto ephoton = photons->cend();
73  reco::PFBlockElementSuperCluster* scbe = nullptr;
74  reco::PhotonRef phoref;
75  for (auto photon = bphoton; photon != ephoton; ++photon) {
76  if (_selector->passPhotonSelection(*photon)) {
77  phoref = reco::PhotonRef(photons, std::distance(bphoton, photon));
78  const reco::SuperClusterRef& scref = photon->superCluster();
79  PFBlockElementSCEqual myEqual(scref);
80  auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
81  if (sc_elem != SCs_end) {
82  scbe = static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
83  scbe->setFromPhoton(true);
84  scbe->setPhotonRef(phoref);
85  scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
86  scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
87  scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
88  scbe->setHoE(photon->hadronicOverEm());
89  } else {
90  scbe = new reco::PFBlockElementSuperCluster(scref);
91  scbe->setFromPhoton(true);
93  scbe->setPhotonRef(phoref);
94  scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
95  scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
96  scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
97  scbe->setHoE(photon->hadronicOverEm());
98  SCs_end = elems.insert(SCs_end, ElementType(scbe));
99  ++SCs_end; // point to element *after* the new one
100  }
101  }
102  } // loop on photons
103  elems.shrink_to_fit();
104 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const std::unordered_map< std::string, SelectionChoices > _selectionTypes
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
ParameterSet const & getParameterSet(std::string const &) const
SelectionChoices _selectionChoice
void setFromPhoton(bool val)
set provenance
void setHcalIso(float val)
set the had Iso
void importToBlock(const edm::Event &, ElementList &) const override
EGPhotonImporter(const edm::ParameterSet &, edm::ConsumesCollector &)
void setTrackIso(float val)
set the track Iso
void setPhotonRef(const PhotonRef &ref)
set photonRef
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
edm::Ref< PhotonCollection > PhotonRef
reference to an object in a collection of Photon objects
Definition: PhotonFwd.h:15
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
edm::EDGetTokenT< reco::PhotonCollection > _src
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
void setEcalIso(float val)
set the ecal Iso
std::unique_ptr< const PhotonSelectorAlgo > _selector