CMS 3D CMS Logo

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 
21 
22 //
23 // class declaration
24 //
25 
27 public:
28  explicit GenJetGenPartMerger(const edm::ParameterSet&);
29  ~GenJetGenPartMerger() override;
30 
31 private:
32  void beginStream(edm::StreamID) override;
33  void produce(edm::Event&, const edm::EventSetup&) override;
34  void endStream() override;
35 
40 };
41 
42 //
43 // constants, enums and typedefs
44 //
45 
46 //
47 // static data member definitions
48 //
49 
50 //
51 // constructors and destructor
52 //
54  : jetToken_(consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("srcJet"))),
55  partToken_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("srcPart"))),
56  cut_(iConfig.getParameter<std::string>("cut")),
57  tauAncToken_(consumes<edm::ValueMap<bool>>(iConfig.getParameter<edm::InputTag>("hasTauAnc"))) {
58  produces<reco::GenJetCollection>("merged");
59  produces<edm::ValueMap<bool>>("hasTauAnc");
60 }
61 
63 
64 //
65 // member functions
66 //
67 
68 // ------------ method called to produce the data ------------
70  using namespace edm;
71  std::unique_ptr<reco::GenJetCollection> merged(new reco::GenJetCollection);
72 
73  std::vector<bool> hasTauAncValues;
74 
76  iEvent.getByToken(jetToken_, jetHandle);
77 
79  iEvent.getByToken(partToken_, partHandle);
80 
81  edm::Handle<edm::ValueMap<bool>> tauAncHandle;
82  iEvent.getByToken(tauAncToken_, tauAncHandle);
83 
84  for (unsigned int ijet = 0; ijet < jetHandle->size(); ++ijet) {
85  auto jet = jetHandle->at(ijet);
86  if (cut_(jet)) {
87  merged->push_back(reco::GenJet(jet));
88  reco::GenJetRef jetRef(jetHandle, ijet);
89  hasTauAncValues.push_back((*tauAncHandle)[jetRef]);
90  }
91  }
92 
93  for (auto& part : *partHandle) {
95  jet.setP4(part.p4());
96  jet.setPdgId(part.pdgId());
97  jet.setCharge(part.charge());
98  merged->push_back(jet);
99  hasTauAncValues.push_back(false);
100  }
101 
102  auto newmerged = iEvent.put(std::move(merged), "merged");
103 
104  std::unique_ptr<edm::ValueMap<bool>> out(new edm::ValueMap<bool>());
106  filler.insert(newmerged, hasTauAncValues.begin(), hasTauAncValues.end());
107  filler.fill();
108 
109  iEvent.put(std::move(out), "hasTauAnc");
110 }
111 
112 // ------------ method called once each stream before processing any runs, lumis or events ------------
114 
115 // ------------ method called once each stream after processing all runs, lumis and events ------------
117 
118 //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_
const StringCutObjectSelector< reco::Candidate > cut_
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