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 {
23  usesResource();
24 
26 
27  produces<reco::GenParticleCollection>("neutrinos");
28  produces<reco::GenParticleCollection>("photons");
29  produces<reco::GenJetCollection>("leptons");
30  produces<reco::GenJetCollection>("jets");
31  produces<reco::GenJetCollection>("fatjets");
32  produces<reco::GenParticleCollection>("consts");
33  produces<reco::GenParticleCollection>("tags");
34  produces<reco::METCollection>("mets");
35 
36  analysisHandler_.addAnalysis(rivetAnalysis_);
37 }
38 
39 void ParticleLevelProducer::addGenJet(Rivet::Jet jet, std::unique_ptr<reco::GenJetCollection> &jets,
40  std::unique_ptr<reco::GenParticleCollection> &consts, edm::RefProd<reco::GenParticleCollection>& constsRefHandle, int &iConstituent,
41  std::unique_ptr<reco::GenParticleCollection> &tags, edm::RefProd<reco::GenParticleCollection>& tagsRefHandle, int &iTag)
42 {
43  const auto pjet = jet.pseudojet();
44 
45  reco::GenJet genJet;
46  genJet.setP4(p4(jet));
47  genJet.setVertex(genVertex_);
48  if ( jet.bTagged() ) genJet.setPdgId(5);
49  else if ( jet.cTagged() ) genJet.setPdgId(4);
50  genJet.setJetArea(pjet.has_area() ? pjet.area() : 0);
51 
52  for ( auto const & p : jet.particles()) {
53  auto pp4 = p4(p);
54  bool match = false; int iMatch = -1;
55  for ( auto const & q : *consts ) {
56  ++iMatch;
57  if (q.p4() == pp4) { match = true; break; }
58  }
59  if (match){
60  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, iMatch)));
61  }
62  else {
63  consts->push_back(reco::GenParticle(p.charge(), pp4, genVertex_, p.pdgId(), 1, true));
64  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
65  }
66  }
67  for ( auto const & p : jet.tags()) {
68  // The tag particles are accessible as jet daughters, so scale down p4 for safety.
69  // p4 needs to be multiplied by 1e20 for fragmentation analysis.
70  auto pp4 = p4(p)*1e-20;
71  bool match = false; int iMatch = -1;
72  for ( auto const & q : *tags ) {
73  ++iMatch;
74  if (q.p4() == pp4) { match = true; break; }
75  }
76  if (match){
77  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, iMatch)));
78  }
79  else {
80  tags->push_back(reco::GenParticle(p.charge(), p4(p)*1e-20, genVertex_, p.pdgId(), 2, true));
81  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, ++iTag)));
82  }
83  }
84 
85  jets->push_back(genJet);
86 }
87 
89 {
90  using namespace Rivet;
92 
93  std::unique_ptr<reco::GenParticleCollection> neutrinos(new reco::GenParticleCollection);
94  std::unique_ptr<reco::GenParticleCollection> photons(new reco::GenParticleCollection);
95  std::unique_ptr<reco::GenJetCollection> leptons(new reco::GenJetCollection);
96  std::unique_ptr<reco::GenJetCollection> jets(new reco::GenJetCollection);
97  std::unique_ptr<reco::GenJetCollection> fatjets(new reco::GenJetCollection);
98  std::unique_ptr<reco::GenParticleCollection> consts(new reco::GenParticleCollection);
99  std::unique_ptr<reco::GenParticleCollection> tags(new reco::GenParticleCollection);
100  std::unique_ptr<reco::METCollection> mets(new reco::METCollection);
101  auto constsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("consts");
102  auto tagsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("tags");
103 
104  edm::Handle<HepMCProduct> srcHandle;
105  event.getByToken(srcToken_, srcHandle);
106 
107  const HepMC::GenEvent* genEvent = srcHandle->GetEvent();
108  analysisHandler_.analyze(*genEvent);
109 
110  // Convert into edm objects
111  // Prompt neutrinos
112  for ( auto const & p : rivetAnalysis_->neutrinos() ) {
113  neutrinos->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
114  }
115  std::sort(neutrinos->begin(), neutrinos->end(), GreaterByPt<reco::Candidate>());
116 
117  // Photons
118  for ( auto const & p : rivetAnalysis_->photons() ) {
119  photons->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
120  }
121  std::sort(photons->begin(), photons->end(), GreaterByPt<reco::Candidate>());
122 
123  // Prompt leptons
124  int iConstituent = -1;
125  for ( auto const & lepton : rivetAnalysis_->leptons() ) {
126  reco::GenJet lepJet;
127  lepJet.setP4(p4(lepton));
128  lepJet.setVertex(genVertex_);
129  lepJet.setPdgId(lepton.pdgId());
130  lepJet.setCharge(lepton.charge());
131 
132  const auto& cl = lepton.constituentLepton();
133  consts->push_back(reco::GenParticle(cl.charge(), p4(cl), genVertex_, cl.pdgId(), 1, true));
134  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
135 
136  for ( auto const & p : lepton.constituentPhotons()) {
137  consts->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pdgId(), 1, true));
138  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
139  }
140 
141  leptons->push_back(lepJet);
142  }
143  std::sort(leptons->begin(), leptons->end(), GreaterByPt<reco::GenJet>());
144 
145  // Jets with constituents and tag particles
146  int iTag = -1;
147  for ( auto jet : rivetAnalysis_->jets() ) {
148  addGenJet(jet, jets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
149  }
150  for ( auto jet : rivetAnalysis_->fatjets() ) {
151  addGenJet(jet, fatjets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
152  }
153 
154  // MET
156  mets->push_back(reco::MET(metP4, genVertex_));
157 
158  event.put(std::move(neutrinos), "neutrinos");
159  event.put(std::move(photons), "photons");
160  event.put(std::move(leptons), "leptons");
161  event.put(std::move(jets), "jets");
162  event.put(std::move(fatjets), "fatjets");
163  event.put(std::move(consts), "consts");
164  event.put(std::move(tags), "tags");
165  event.put(std::move(mets), "mets");
166 }
167 
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
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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
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
std::vector< DressedLepton > leptons() const
Definition: RivetAnalysis.h:26
Jets made from MC generator particles.
Definition: GenJet.h:24
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
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:510
Definition: event.py:1
math::PtEtaPhiELorentzVectorF LorentzVector
reco::Particle::Point genVertex_