CMS 3D CMS Logo

EGPhotonImporter.cc
Go to the documentation of this file.
10 
11 #include <unordered_map>
12 
14 public:
16 
18 
19  void importToBlock( const edm::Event& ,
20  ElementList& ) const override;
21 
22 private:
24  const std::unordered_map<std::string,SelectionChoices> _selectionTypes;
26  std::unique_ptr<const PhotonSelectorAlgo> _selector;
28 
29 };
30 
33  "EGPhotonImporter");
34 
36  edm::ConsumesCollector& sumes) :
37  BlockElementImporterBase(conf,sumes),
38  _src(sumes.consumes<reco::PhotonCollection>(conf.getParameter<edm::InputTag>("source"))),
40  {"CombinedDetectorIso",EGPhotonImporter::CombinedDetectorIso} }),
41  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {
42  const std::string& selChoice =
43  conf.getParameter<std::string>("SelectionChoice");
44  _selectionChoice = _selectionTypes.at(selChoice);
45  const edm::ParameterSet& selDef =
46  conf.getParameterSet("SelectionDefinition");
47  const float minEt = selDef.getParameter<double>("minEt");
48  const float trackIso_const = selDef.getParameter<double>("trackIsoConstTerm");
49  const float trackIso_slope = selDef.getParameter<double>("trackIsoSlopeTerm");
50  const float ecalIso_const = selDef.getParameter<double>("ecalIsoConstTerm");
51  const float ecalIso_slope = selDef.getParameter<double>("ecalIsoSlopeTerm");
52  const float hcalIso_const = selDef.getParameter<double>("hcalIsoConstTerm");
53  const float hcalIso_slope = selDef.getParameter<double>("hcalIsoSlopeTerm");
54  const float hoe = selDef.getParameter<double>("HoverE");
55  const float loose_hoe = selDef.getParameter<double>("LooseHoverE");
56  const float combIso = selDef.getParameter<double>("combIsoConstTerm");
57  _selector.reset(new PhotonSelectorAlgo((float)_selectionChoice,
58  minEt,
59  trackIso_const, trackIso_slope,
60  ecalIso_const, ecalIso_slope,
61  hcalIso_const, hcalIso_slope,
62  hoe,
63  combIso,
64  loose_hoe));
65 }
66 
71  auto photons = e.getHandle(_src);
72  elems.reserve(elems.size()+photons->size());
73  // setup our elements so that all the SCs are grouped together
74  auto SCs_end = std::partition(elems.begin(),elems.end(),
75  [](const ElementType& a){
76  return a->type() == reco::PFBlockElement::SC;
77  });
78  //now add the photons
79  auto bphoton = photons->cbegin();
80  auto ephoton = photons->cend();
81  reco::PFBlockElementSuperCluster* scbe = nullptr;
82  reco::PhotonRef phoref;
83  for( auto photon = bphoton; photon != ephoton; ++photon ) {
84  if( _selector->passPhotonSelection(*photon) ) {
85  phoref = reco::PhotonRef(photons,std::distance(bphoton,photon));
86  const reco::SuperClusterRef& scref = photon->superCluster();
87  PFBlockElementSCEqual myEqual(scref);
88  auto sc_elem = std::find_if(elems.begin(),SCs_end,myEqual);
89  if( sc_elem != SCs_end ) {
90  scbe = static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
91  scbe->setFromPhoton(true);
92  scbe->setPhotonRef(phoref);
93  scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
94  scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
95  scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
96  scbe->setHoE(photon->hadronicOverEm());
97  } else {
98  scbe = new reco::PFBlockElementSuperCluster(scref);
99  scbe->setFromPhoton(true);
101  scbe->setPhotonRef(phoref);
102  scbe->setTrackIso(photon->trkSumPtHollowConeDR04());
103  scbe->setEcalIso(photon->ecalRecHitSumEtConeDR04());
104  scbe->setHcalIso(photon->hcalTowerSumEtConeDR04());
105  scbe->setHoE(photon->hadronicOverEm());
106  SCs_end = elems.insert(SCs_end,ElementType(scbe));
107  ++SCs_end; // point to element *after* the new one
108  }
109  }
110  }// loop on photons
111  elems.shrink_to_fit();
112 }
T getParameter(std::string const &) const
SelectionChoices _selectionChoice
void setFromPhoton(bool val)
set provenance
void setHcalIso(float val)
set the had Iso
EGPhotonImporter(const edm::ParameterSet &, edm::ConsumesCollector &)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:539
void setTrackIso(float val)
set the track Iso
void setPhotonRef(const PhotonRef &ref)
set photonRef
void importToBlock(const edm::Event &, ElementList &) const override
ParameterSet const & getParameterSet(std::string const &) const
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
const std::unordered_map< std::string, SelectionChoices > _selectionTypes