CMS 3D CMS Logo

ParticleLevelProducer.cc
Go to the documentation of this file.
2 
11 
12 #include "Rivet/Analysis.hh"
13 
14 using namespace std;
15 using namespace edm;
16 using namespace reco;
17 using namespace Rivet;
18 
20  : srcToken_(consumes<edm::HepMCProduct>(pset.getParameter<edm::InputTag>("src"))),
21  rivetAnalysis_(new Rivet::RivetAnalysis(pset)) {
22  usesResource("Rivet");
23 
25 
26  produces<reco::GenParticleCollection>("neutrinos");
27  produces<reco::GenParticleCollection>("photons");
28  produces<reco::GenJetCollection>("leptons");
29  produces<reco::GenJetCollection>("jets");
30  produces<reco::GenJetCollection>("fatjets");
31  produces<reco::GenParticleCollection>("consts");
32  produces<reco::GenParticleCollection>("tags");
33  produces<reco::METCollection>("mets");
34 
35  analysisHandler_.setIgnoreBeams(true);
36  analysisHandler_.addAnalysis(rivetAnalysis_);
37 }
38 
40  std::unique_ptr<reco::GenJetCollection>& jets,
41  std::unique_ptr<reco::GenParticleCollection>& consts,
43  int& iConstituent,
44  std::unique_ptr<reco::GenParticleCollection>& tags,
46  int& iTag) {
47  const auto pjet = jet.pseudojet();
48 
49  reco::GenJet genJet;
50  genJet.setP4(p4(jet));
51  genJet.setVertex(genVertex_);
52  if (jet.bTagged())
53  genJet.setPdgId(5);
54  else if (jet.cTagged())
55  genJet.setPdgId(4);
56  genJet.setJetArea(pjet.has_area() ? pjet.area() : 0);
57 
58  for (auto const& p : jet.particles()) {
59  auto pp4 = p4(p);
60  bool match = false;
61  int iMatch = -1;
62  for (auto const& q : *consts) {
63  ++iMatch;
64  if (q.p4() == pp4) {
65  match = true;
66  break;
67  }
68  }
69  if (match) {
70  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, iMatch)));
71  } else {
72  consts->push_back(reco::GenParticle(p.charge(), pp4, genVertex_, p.pdgId(), 1, true));
73  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
74  }
75  }
76  for (auto const& p : jet.tags()) {
77  // The tag particles are accessible as jet daughters, so scale down p4 for safety.
78  // p4 needs to be multiplied by 1e20 for fragmentation analysis.
79  auto pp4 = p4(p) * 1e-20;
80  bool match = false;
81  int iMatch = -1;
82  for (auto const& q : *tags) {
83  ++iMatch;
84  if (q.p4() == pp4) {
85  match = true;
86  break;
87  }
88  }
89  if (match) {
90  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, iMatch)));
91  } else {
92  tags->push_back(reco::GenParticle(p.charge(), p4(p) * 1e-20, genVertex_, p.pdgId(), 2, true));
93  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, ++iTag)));
94  // Also save lepton+neutrino daughters of tag particles
95  int iTagMother = iTag;
96  for (auto const& d : p.constituents()) {
97  ++iTag;
98  int d_status = (d.isStable()) ? 1 : 2;
99  tags->push_back(reco::GenParticle(d.charge(), p4(d) * 1e-20, genVertex_, d.pid(), d_status, true));
100  tags->at(iTag).addMother(reco::GenParticleRef(tagsRefHandle, iTagMother));
101  tags->at(iTagMother).addDaughter(reco::GenParticleRef(tagsRefHandle, iTag));
102  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, iTag)));
103  }
104  }
105  }
106 
107  jets->push_back(genJet);
108 }
109 
111  using namespace Rivet;
113 
114  std::unique_ptr<reco::GenParticleCollection> neutrinos(new reco::GenParticleCollection);
115  std::unique_ptr<reco::GenParticleCollection> photons(new reco::GenParticleCollection);
116  std::unique_ptr<reco::GenJetCollection> leptons(new reco::GenJetCollection);
117  std::unique_ptr<reco::GenJetCollection> jets(new reco::GenJetCollection);
118  std::unique_ptr<reco::GenJetCollection> fatjets(new reco::GenJetCollection);
119  std::unique_ptr<reco::GenParticleCollection> consts(new reco::GenParticleCollection);
120  std::unique_ptr<reco::GenParticleCollection> tags(new reco::GenParticleCollection);
121  std::unique_ptr<reco::METCollection> mets(new reco::METCollection);
122  auto constsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("consts");
123  auto tagsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("tags");
124 
125  edm::Handle<HepMCProduct> srcHandle;
126  event.getByToken(srcToken_, srcHandle);
127 
128  const HepMC::GenEvent* genEvent = srcHandle->GetEvent();
129  analysisHandler_.analyze(*genEvent);
130 
131  // Convert into edm objects
132  // Prompt neutrinos
133  for (auto const& p : rivetAnalysis_->neutrinos()) {
134  neutrinos->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
135  }
136  std::sort(neutrinos->begin(), neutrinos->end(), GreaterByPt<reco::Candidate>());
137 
138  // Photons
139  for (auto const& p : rivetAnalysis_->photons()) {
140  photons->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
141  }
142  std::sort(photons->begin(), photons->end(), GreaterByPt<reco::Candidate>());
143 
144  // Prompt leptons
145  int iConstituent = -1;
146  int iTag = -1;
147  for (auto const& lepton : rivetAnalysis_->leptons()) {
148  reco::GenJet lepJet;
149  lepJet.setP4(p4(lepton));
150  lepJet.setVertex(genVertex_);
151  lepJet.setPdgId(lepton.pdgId());
152  lepJet.setCharge(lepton.charge());
153 
154  for (auto const & p : lepton.constituents()) {
155  // ghost taus (momentum scaled with 10e-20 in RivetAnalysis.h already)
156  if (p.abspid() == 15) {
157  tags->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 2, true));
158  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, ++iTag)));
159  }
160  // electrons, muons, photons
161  else {
162  consts->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
163  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
164  }
165  }
166 
167  leptons->push_back(lepJet);
168  }
169  std::sort(leptons->begin(), leptons->end(), GreaterByPt<reco::GenJet>());
170 
171  // Jets with constituents and tag particles
172  for (auto jet : rivetAnalysis_->jets()) {
173  addGenJet(jet, jets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
174  }
175  for (auto jet : rivetAnalysis_->fatjets()) {
176  addGenJet(jet, fatjets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
177  }
178 
179  // MET
181  rivetAnalysis_->met().y(),
182  0.,
183  sqrt(pow(rivetAnalysis_->met().x(), 2) + pow(rivetAnalysis_->met().y(), 2)));
184  mets->push_back(reco::MET(metP4, genVertex_));
185 
186  event.put(std::move(neutrinos), "neutrinos");
187  event.put(std::move(photons), "photons");
188  event.put(std::move(leptons), "leptons");
189  event.put(std::move(jets), "jets");
190  event.put(std::move(fatjets), "fatjets");
191  event.put(std::move(consts), "consts");
192  event.put(std::move(tags), "tags");
193  event.put(std::move(mets), "mets");
194 }
195 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
Vector3 met() const
Definition: RivetAnalysis.h:31
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
std::vector< GenJet > GenJetCollection
collection of GenJet objects
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
Jets fatjets() const
Definition: RivetAnalysis.h:30
void setVertex(const Point &vertex) override
set vertex
ParticleVector neutrinos() const
Definition: RivetAnalysis.h:28
void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
Jets jets() const
Definition: RivetAnalysis.h:29
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:103
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParticleVector photons() const
Definition: RivetAnalysis.h:27
reco::Candidate::LorentzVector p4(const T &p) const
Definition: MET.h:42
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Jets made from MC generator particles.
Definition: GenJet.h:25
void addGenJet(Rivet::Jet jet, std::unique_ptr< reco::GenJetCollection > &jets, std::unique_ptr< reco::GenParticleCollection > &consts, edm::RefProd< reco::GenParticleCollection > &constsRefHandle, int &iConstituent, std::unique_ptr< reco::GenParticleCollection > &tags, edm::RefProd< reco::GenParticleCollection > &tagsRefHandle, int &iTag)
Rivet::RivetAnalysis * rivetAnalysis_
ParticleLevelProducer(const edm::ParameterSet &pset)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
const edm::EDGetTokenT< edm::HepMCProduct > srcToken_
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
ParticleVector leptons() const
Definition: RivetAnalysis.h:26
fixed size matrix
HLT enums.
Rivet::AnalysisHandler analysisHandler_
void addDaughter(const CandidatePtr &)
add a daughter via a reference
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
void setPdgId(int pdgId) final
void setP4(const LorentzVector &p4) final
set 4-momentum
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
math::PtEtaPhiELorentzVectorF LorentzVector
reco::Particle::Point genVertex_