CMS 3D CMS Logo

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