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 
13 
16 
19 
23 
24 //
25 // class declaration
26 //
27 
29 public:
30  explicit GenJetGenPartMerger(const edm::ParameterSet&);
31  ~GenJetGenPartMerger() override;
32  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
33 
34 private:
35  void beginStream(edm::StreamID) override;
36  void produce(edm::Event&, const edm::EventSetup&) override;
37  void endStream() override;
38 
43 };
44 
45 //
46 // constants, enums and typedefs
47 //
48 
49 //
50 // static data member definitions
51 //
52 
53 //
54 // constructors and destructor
55 //
57  : jetToken_(consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("srcJet"))),
58  partToken_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("srcPart"))),
59  cut_(iConfig.getParameter<std::string>("cut")),
60  tauAncToken_(consumes<edm::ValueMap<bool>>(iConfig.getParameter<edm::InputTag>("hasTauAnc"))) {
61  produces<reco::GenJetCollection>("merged");
62  produces<edm::ValueMap<bool>>("hasTauAnc");
63 }
64 
66 
69  desc.add<edm::InputTag>("srcJet")->setComment("reco::GenJetCollection input collection");
70  desc.add<edm::InputTag>("srcPart")->setComment("reco::GenParticleCollection input collection");
71  desc.add<std::string>("cut")->setComment("a selection to apply to GenJet");
72  desc.add<edm::InputTag>("hasTauAnc")->setComment("value map defining GenJet tau origin");
73  descriptions.addWithDefaultLabel(desc);
74 }
75 //
76 // member functions
77 //
78 
79 // ------------ method called to produce the data ------------
81  using namespace edm;
82  auto merged = std::make_unique<reco::GenJetCollection>();
83 
84  std::vector<bool> hasTauAncValues;
85 
86  auto jetHandle = iEvent.getHandle(jetToken_);
87  const auto& partProd = iEvent.get(partToken_);
88  const auto& tauAncProd = iEvent.get(tauAncToken_);
89 
90  for (unsigned int ijet = 0; ijet < jetHandle->size(); ++ijet) {
91  auto jet = jetHandle->at(ijet);
92  if (cut_(jet)) {
93  merged->push_back(reco::GenJet(jet));
94  reco::GenJetRef jetRef(jetHandle, ijet);
95  hasTauAncValues.push_back(tauAncProd[jetRef]);
96  }
97  }
98 
99  for (const auto& part : partProd) {
101  jet.setP4(part.p4());
102  jet.setPdgId(part.pdgId());
103  jet.setCharge(part.charge());
104  merged->push_back(jet);
105  hasTauAncValues.push_back(false);
106  }
107 
108  auto newmerged = iEvent.put(std::move(merged), "merged");
109 
110  auto out = std::make_unique<edm::ValueMap<bool>>();
112  filler.insert(newmerged, hasTauAncValues.begin(), hasTauAncValues.end());
113  filler.fill();
114 
115  iEvent.put(std::move(out), "hasTauAnc");
116 }
117 
118 // ------------ method called once each stream before processing any runs, lumis or events ------------
120 
121 // ------------ method called once each stream after processing all runs, lumis and events ------------
123 
124 //define this as a plug-in
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::vector< GenJet > GenJetCollection
collection of GenJet objects
void beginStream(edm::StreamID) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
GenJetGenPartMerger(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
Jets made from MC generator particles.
Definition: GenJet.h:23
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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_
def move(src, dest)
Definition: eostools.py:511