CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ParticleLevelProducer Class Reference

#include <ParticleLevelProducer.h>

Inheritance diagram for ParticleLevelProducer:
edm::one::EDProducer< edm::one::SharedResources > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 ParticleLevelProducer (const edm::ParameterSet &pset)
 
void produce (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 ~ParticleLevelProducer () override
 
- Public Member Functions inherited from edm::one::EDProducer< edm::one::SharedResources >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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)
 
template<typename T >
reco::Candidate::LorentzVector p4 (const T &p) const
 

Private Attributes

Rivet::AnalysisHandler analysisHandler_
 
reco::Particle::Point genVertex_
 
Rivet::RivetAnalysisrivetAnalysis_
 
const edm::EDGetTokenT< edm::HepMCProductsrcToken_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 18 of file ParticleLevelProducer.h.

Constructor & Destructor Documentation

ParticleLevelProducer::ParticleLevelProducer ( const edm::ParameterSet pset)

Definition at line 19 of file ParticleLevelProducer.cc.

References analysisHandler_, genVertex_, and rivetAnalysis_.

20  : srcToken_(consumes<edm::HepMCProduct>(pset.getParameter<edm::InputTag>("src"))),
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 }
T getParameter(std::string const &) const
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Rivet::RivetAnalysis * rivetAnalysis_
const edm::EDGetTokenT< edm::HepMCProduct > srcToken_
Rivet::AnalysisHandler analysisHandler_
reco::Particle::Point genVertex_
ParticleLevelProducer::~ParticleLevelProducer ( )
inlineoverride

Member Function Documentation

void ParticleLevelProducer::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 
)
private

Definition at line 39 of file ParticleLevelProducer.cc.

References reco::CompositePtrCandidate::addDaughter(), MillePedeFileConverter_cfg::e, genVertex_, match(), AlCaHLTBitMon_ParallelJobs::p, p4(), data-class-funcs::q, edm::refToPtr(), reco::Jet::setJetArea(), reco::LeafCandidate::setP4(), reco::LeafCandidate::setPdgId(), and reco::LeafCandidate::setVertex().

Referenced by produce(), and ~ParticleLevelProducer().

46  {
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.pid(), 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.pid(), 2, true));
93  genJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, ++iTag)));
94  }
95  }
96 
97  jets->push_back(genJet);
98 }
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
void setVertex(const Point &vertex) override
set vertex
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:101
reco::Candidate::LorentzVector p4(const T &p) const
Jets made from MC generator particles.
Definition: GenJet.h:23
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
reco::Particle::Point genVertex_
template<typename T >
reco::Candidate::LorentzVector ParticleLevelProducer::p4 ( const T p) const
inlineprivate

Definition at line 35 of file ParticleLevelProducer.h.

Referenced by addGenJet(), Tau.Tau::dxy_approx(), Tau.Tau::dz(), Lepton.Lepton::p4WithFSR(), and produce().

35  {
36  return reco::Candidate::LorentzVector(p.px(), p.py(), p.pz(), p.energy());
37  }
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void ParticleLevelProducer::produce ( edm::Event event,
const edm::EventSetup eventSetup 
)
override

Definition at line 100 of file ParticleLevelProducer.cc.

References reco::CompositePtrCandidate::addDaughter(), addGenJet(), analysisHandler_, DEFINE_FWK_MODULE, Rivet::RivetAnalysis::fatjets(), nano_cff::genEvent, genVertex_, metsig::jet, Rivet::RivetAnalysis::jets(), singleTopDQM_cfi::jets, Rivet::RivetAnalysis::leptons(), HLT_2018_cff::leptons, Rivet::RivetAnalysis::met(), singleTopDQM_cfi::mets, eostools::move(), TtSemiLepHypHitFit_cfi::neutrinos, Rivet::RivetAnalysis::neutrinos(), AlCaHLTBitMon_ParallelJobs::p, p4(), Rivet::RivetAnalysis::photons(), BPHMonitor_cfi::photons, funct::pow(), edm::refToPtr(), rivetAnalysis_, reco::LeafCandidate::setCharge(), reco::LeafCandidate::setP4(), reco::LeafCandidate::setPdgId(), reco::LeafCandidate::setVertex(), mathSSE::sqrt(), srcToken_, and triggerMatcherToHLTDebug_cfi::tags.

Referenced by ~ParticleLevelProducer().

100  {
101  using namespace Rivet;
103 
104  std::unique_ptr<reco::GenParticleCollection> neutrinos(new reco::GenParticleCollection);
105  std::unique_ptr<reco::GenParticleCollection> photons(new reco::GenParticleCollection);
106  std::unique_ptr<reco::GenJetCollection> leptons(new reco::GenJetCollection);
107  std::unique_ptr<reco::GenJetCollection> jets(new reco::GenJetCollection);
108  std::unique_ptr<reco::GenJetCollection> fatjets(new reco::GenJetCollection);
109  std::unique_ptr<reco::GenParticleCollection> consts(new reco::GenParticleCollection);
110  std::unique_ptr<reco::GenParticleCollection> tags(new reco::GenParticleCollection);
111  std::unique_ptr<reco::METCollection> mets(new reco::METCollection);
112  auto constsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("consts");
113  auto tagsRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("tags");
114 
115  edm::Handle<HepMCProduct> srcHandle;
116  event.getByToken(srcToken_, srcHandle);
117 
118  const HepMC::GenEvent* genEvent = srcHandle->GetEvent();
119  analysisHandler_.analyze(*genEvent);
120 
121  // Convert into edm objects
122  // Prompt neutrinos
123  for (auto const& p : rivetAnalysis_->neutrinos()) {
124  neutrinos->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pid(), 1, true));
125  }
126  std::sort(neutrinos->begin(), neutrinos->end(), GreaterByPt<reco::Candidate>());
127 
128  // Photons
129  for (auto const& p : rivetAnalysis_->photons()) {
130  photons->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pid(), 1, true));
131  }
132  std::sort(photons->begin(), photons->end(), GreaterByPt<reco::Candidate>());
133 
134  // Prompt leptons
135  int iConstituent = -1;
136  int iTag = -1;
137  for (auto const& lepton : rivetAnalysis_->leptons()) {
138  reco::GenJet lepJet;
139  lepJet.setP4(p4(lepton));
140  lepJet.setVertex(genVertex_);
141  lepJet.setPdgId(lepton.pid());
142  lepJet.setCharge(lepton.charge());
143 
144  for (auto const& p : lepton.constituents()) {
145  // ghost taus (momentum scaled with 10e-20 in RivetAnalysis.h already)
146  if (p.abspid() == 15) {
147  tags->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pid(), 2, true));
148  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(tagsRefHandle, ++iTag)));
149  }
150  // electrons, muons, photons
151  else {
152  consts->push_back(reco::GenParticle(p.charge(), p4(p), genVertex_, p.pid(), 1, true));
153  lepJet.addDaughter(edm::refToPtr(reco::GenParticleRef(constsRefHandle, ++iConstituent)));
154  }
155  }
156 
157  leptons->push_back(lepJet);
158  }
159  std::sort(leptons->begin(), leptons->end(), GreaterByPt<reco::GenJet>());
160 
161  // Jets with constituents and tag particles
162  for (auto jet : rivetAnalysis_->jets()) {
163  addGenJet(jet, jets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
164  }
165  for (auto jet : rivetAnalysis_->fatjets()) {
166  addGenJet(jet, fatjets, consts, constsRefHandle, iConstituent, tags, tagsRefHandle, iTag);
167  }
168 
169  // MET
171  rivetAnalysis_->met().y(),
172  0.,
173  sqrt(pow(rivetAnalysis_->met().x(), 2) + pow(rivetAnalysis_->met().y(), 2)));
174  mets->push_back(reco::MET(metP4, genVertex_));
175 
176  event.put(std::move(neutrinos), "neutrinos");
177  event.put(std::move(photons), "photons");
178  event.put(std::move(leptons), "leptons");
179  event.put(std::move(jets), "jets");
180  event.put(std::move(fatjets), "fatjets");
181  event.put(std::move(consts), "consts");
182  event.put(std::move(tags), "tags");
183  event.put(std::move(mets), "mets");
184 }
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
std::vector< GenJet > GenJetCollection
collection of GenJet objects
Particles leptons() const
Definition: RivetAnalysis.h:25
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
Jets jets() const
Definition: RivetAnalysis.h:28
reco::Candidate::LorentzVector p4(const T &p) const
Definition: MET.h:41
T sqrt(T t)
Definition: SSEVec.h:19
Jets made from MC generator particles.
Definition: GenJet.h:23
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_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
const edm::EDGetTokenT< edm::HepMCProduct > srcToken_
Rivet::AnalysisHandler analysisHandler_
void addDaughter(const CandidatePtr &)
add a daughter via a reference
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:30
def move(src, dest)
Definition: eostools.py:511
math::PtEtaPhiELorentzVectorF LorentzVector
reco::Particle::Point genVertex_

Member Data Documentation

Rivet::AnalysisHandler ParticleLevelProducer::analysisHandler_
private

Definition at line 44 of file ParticleLevelProducer.h.

Referenced by ParticleLevelProducer(), and produce().

reco::Particle::Point ParticleLevelProducer::genVertex_
private

Definition at line 41 of file ParticleLevelProducer.h.

Referenced by addGenJet(), ParticleLevelProducer(), and produce().

Rivet::RivetAnalysis* ParticleLevelProducer::rivetAnalysis_
private

Definition at line 43 of file ParticleLevelProducer.h.

Referenced by ParticleLevelProducer(), and produce().

const edm::EDGetTokenT<edm::HepMCProduct> ParticleLevelProducer::srcToken_
private

Definition at line 39 of file ParticleLevelProducer.h.

Referenced by produce().