CMS 3D CMS Logo

PFEGammaToCandidate.cc
Go to the documentation of this file.
1 
22 
24 public:
25  explicit PFEGammaToCandidate(const edm::ParameterSet &iConfig);
26  ~PFEGammaToCandidate() override = default;
27 
28  void produce(edm::StreamID iID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override;
29  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
30 
31 private:
35 
36  template <typename T>
38  const edm::EDGetTokenT<edm::View<T>> &colltoken,
39  const edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> &oldmaptoken,
40  const std::string &name) const {
41  auto const &coll = iEvent.get(colltoken);
42 
44  iEvent.getByToken(oldmaptoken, oldmap);
45 
46  auto result = std::make_unique<std::vector<reco::VertexCompositePtrCandidate>>();
47  for (auto const &obj : coll) {
49  for (const reco::PFCandidateRef &pfRef : (*oldmap)[obj.originalObjectRef()])
50  result->back().addDaughter(refToPtr(pfRef));
51  }
52  iEvent.put(std::move(result), name);
53  }
54 };
55 
57  : photons_(consumes<edm::View<pat::Photon>>(iConfig.getParameter<edm::InputTag>("photons"))),
58  electrons_(consumes<edm::View<pat::Electron>>(iConfig.getParameter<edm::InputTag>("electrons"))),
59  photon2pf_(
60  consumes<edm::ValueMap<std::vector<reco::PFCandidateRef>>>(iConfig.getParameter<edm::InputTag>("photon2pf"))),
61  electron2pf_(consumes<edm::ValueMap<std::vector<reco::PFCandidateRef>>>(
62  iConfig.getParameter<edm::InputTag>("electron2pf"))) {
63  produces<std::vector<reco::VertexCompositePtrCandidate>>("photons");
64  produces<std::vector<reco::VertexCompositePtrCandidate>>("electrons");
65 }
66 
68  run<pat::Photon>(iEvent, photons_, photon2pf_, "photons");
69  run<pat::Electron>(iEvent, electrons_, electron2pf_, "electrons");
70 }
71 
74  desc.add<edm::InputTag>("photons", edm::InputTag("selectedPatPhotons"));
75  desc.add<edm::InputTag>("electrons", edm::InputTag("selectedPatElectrons"));
76  desc.add<edm::InputTag>("photon2pf", edm::InputTag("particleBasedIsolation", "gedPhotons"));
77  desc.add<edm::InputTag>("electron2pf", edm::InputTag("particleBasedIsolation", "gedGsfElectrons"));
78  descriptions.addWithDefaultLabel(desc);
79 }
80 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
Definition: Photon.py:1
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
const edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > photon2pf_
const edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > electron2pf_
void produce(edm::StreamID iID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
void run(edm::Event &iEvent, const edm::EDGetTokenT< edm::View< T >> &colltoken, const edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef >>> &oldmaptoken, const std::string &name) const
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PFEGammaToCandidate(const edm::ParameterSet &iConfig)
fixed size matrix
HLT enums.
edm::Ref< l1t::PFCandidateCollection > PFCandidateRef
Definition: PFCandidate.h:83
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~PFEGammaToCandidate() override=default
const edm::EDGetTokenT< edm::View< pat::Photon > > photons_
const edm::EDGetTokenT< edm::View< pat::Electron > > electrons_
def move(src, dest)
Definition: eostools.py:511