CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenJetGenPartMerger.cc
Go to the documentation of this file.
1 
2 // system include files
3 #include <memory>
4 
5 // user include files
8 
11 
14 
17 
19 
20 //
21 // class declaration
22 //
23 
25 public:
26  explicit GenJetGenPartMerger(const edm::ParameterSet&);
27  ~GenJetGenPartMerger() override;
28 
29 private:
30  void beginStream(edm::StreamID) override;
31  void produce(edm::Event&, const edm::EventSetup&) override;
32  void endStream() override;
33 
37 };
38 
39 //
40 // constants, enums and typedefs
41 //
42 
43 //
44 // static data member definitions
45 //
46 
47 //
48 // constructors and destructor
49 //
51  : jetToken_(consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("srcJet"))),
52  partToken_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("srcPart"))),
53  tauAncToken_(consumes<edm::ValueMap<bool>>(iConfig.getParameter<edm::InputTag>("hasTauAnc"))) {
54  produces<reco::GenJetCollection>("merged");
55  produces<edm::ValueMap<bool>>("hasTauAnc");
56 }
57 
59 
60 //
61 // member functions
62 //
63 
64 // ------------ method called to produce the data ------------
66  using namespace edm;
67  std::unique_ptr<reco::GenJetCollection> merged(new reco::GenJetCollection);
68 
69  std::vector<bool> hasTauAncValues;
70 
72  iEvent.getByToken(jetToken_, jetHandle);
73 
75  iEvent.getByToken(partToken_, partHandle);
76 
77  edm::Handle<edm::ValueMap<bool>> tauAncHandle;
78  iEvent.getByToken(tauAncToken_, tauAncHandle);
79 
80  for (unsigned int ijet = 0; ijet < jetHandle->size(); ++ijet) {
81  auto jet = jetHandle->at(ijet);
82  merged->push_back(reco::GenJet(jet));
83  reco::GenJetRef jetRef(jetHandle, ijet);
84  hasTauAncValues.push_back((*tauAncHandle)[jetRef]);
85  }
86 
87  for (auto& part : *partHandle) {
89  jet.setP4(part.p4());
90  jet.setPdgId(part.pdgId());
91  jet.setCharge(part.charge());
92  merged->push_back(jet);
93  hasTauAncValues.push_back(false);
94  }
95 
96  auto newmerged = iEvent.put(std::move(merged), "merged");
97 
98  std::unique_ptr<edm::ValueMap<bool>> out(new edm::ValueMap<bool>());
100  filler.insert(newmerged, hasTauAncValues.begin(), hasTauAncValues.end());
101  filler.fill();
102 
103  iEvent.put(std::move(out), "hasTauAnc");
104 }
105 
106 // ------------ method called once each stream before processing any runs, lumis or events ------------
108 
109 // ------------ method called once each stream after processing all runs, lumis and events ------------
111 
112 //define this as a plug-in
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< GenJet > GenJetCollection
collection of GenJet objects
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
void beginStream(edm::StreamID) override
void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
GenJetGenPartMerger(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jets made from MC generator particles.
Definition: GenJet.h:25
void produce(edm::Event &, const edm::EventSetup &) override
part
Definition: HCALResponse.h:20
const edm::EDGetTokenT< reco::GenJetCollection > jetToken_
fixed size matrix
HLT enums.
const edm::EDGetTokenT< reco::GenParticleCollection > partToken_
const edm::EDGetTokenT< edm::ValueMap< bool > > tauAncToken_
void setPdgId(int pdgId) final
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511