CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes
pat::GenPlusSimParticleProducer Class Reference

Produces reco::GenParticle from SimTracks. More...

Inheritance diagram for pat::GenPlusSimParticleProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

struct  LessById
 

Public Member Functions

 GenPlusSimParticleProducer (const edm::ParameterSet &)
 
 ~GenPlusSimParticleProducer () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () 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 Types

typedef StringCutObjectSelector< reco::GenParticleStrFilter
 

Private Member Functions

void addGenParticle (const SimTrack &stMom, const SimTrack &stDau, unsigned int momidx, const edm::SimTrackContainer &simtks, const edm::SimVertexContainer &simvtxs, reco::GenParticleCollection &mergedGens, const reco::GenParticleRefProd ref, std::vector< int > &genBarcodes, bool &barcodesAreSorted) const
 Try to link the GEANT particle to the generator particle it came from. More...
 
void endJob () override
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::unique_ptr< StrFilterfilter_
 
bool firstEvent_
 
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
 
edm::EDGetTokenT< reco::GenParticleCollectiongensToken_
 Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of barcodes, aligned with them. More...
 
std::set< int > pdgIds_
 
std::vector< PdtEntrypdts_
 
int setStatus_
 
edm::EDGetTokenT< edm::SimTrackContainersimtracksToken_
 
edm::EDGetTokenT< edm::SimVertexContainersimverticesToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer 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::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- 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

Produces reco::GenParticle from SimTracks.

The GenPlusSimParticleProducer produces GenParticles from SimTracks, so they can be used for MC matching. A copy of the original genParticle list is made and the genParticles created from the GEANT/FAMOS particles are added to the list including all ancestors and correct mother/daughter references

Sample useage in cfg.py file: process.genParticlePlusGEANT = cms.EDProducer("GenPlusSimParticleProducer", src = cms.InputTag("g4SimHits"), # use "fastSimProducer" for FastSim setStatus = cms.int32(8), # set status = 8 for GEANT GPs particleTypes = cms.vstring("pi+"), # also picks pi- (optional) filter = cms.vstring("pt > 0.0"), # just for testing genParticles = cms.InputTag("genParticles") # original genParticle list )

Author
Jordan Tucker (original module), Keith Ulmer (generalization)

Definition at line 39 of file GenPlusSimParticleProducer.cc.

Member Typedef Documentation

Definition at line 55 of file GenPlusSimParticleProducer.cc.

Constructor & Destructor Documentation

GenPlusSimParticleProducer::GenPlusSimParticleProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 98 of file GenPlusSimParticleProducer.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), ALCARECOTkAlBeamHalo_cff::filter, filter_, edm::ParameterSet::getParameter(), and pdts_.

99  : firstEvent_(true),
100  simtracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
102  consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
103  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
104  gensToken_(consumes<GenParticleCollection>(
105  cfg.getParameter<InputTag>("genParticles"))), // get the genParticles to add GEANT particles to
106  genBarcodesToken_(consumes<std::vector<int>>(
107  cfg.getParameter<InputTag>("genParticles"))) // get the genParticles to add GEANT particles to
108 {
109  // Possibly allow a list of particle types
110  if (cfg.exists("particleTypes")) {
111  pdts_ = cfg.getParameter<vector<PdtEntry>>("particleTypes");
112  }
113 
114  // Possibly allow a string cut
115  if (cfg.existsAs<string>("filter")) {
116  string filter = cfg.getParameter<string>("filter");
117  if (!filter.empty()) {
118  filter_ = std::make_unique<StrFilter>(filter);
119  }
120  }
121  produces<GenParticleCollection>();
122  produces<vector<int>>();
123 }
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:160
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::unique_ptr< StrFilter > filter_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::GenParticleCollection > gensToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_
pat::GenPlusSimParticleProducer::~GenPlusSimParticleProducer ( )
inlineoverride

Definition at line 42 of file GenPlusSimParticleProducer.cc.

References produce().

42 {}

Member Function Documentation

void GenPlusSimParticleProducer::addGenParticle ( const SimTrack stMom,
const SimTrack stDau,
unsigned int  momidx,
const edm::SimTrackContainer simtks,
const edm::SimVertexContainer simvtxs,
reco::GenParticleCollection mergedGens,
const reco::GenParticleRefProd  ref,
std::vector< int > &  genBarcodes,
bool &  barcodesAreSorted 
) const
private

Try to link the GEANT particle to the generator particle it came from.

Arguments: – Specific – gp: GenParticle made from the GEANT particle st: The GEANT simTrack for which we create a genParticle

– GEANT related – simtks: A list of GEANT tracks, sorted by track id simvtxs: The list of GEANT vertices, in their natural order (skimtks have pointers into this vector!)

GenParticle related – gens : Handle to the GenParticles, to make the ref to genBarcodes : Barcodes for each GenParticle, in a vector aligned to the GenParticleCollection. barcodesAreSorted: true if the barcodes are sorted (which means I can use binary_search on them)

Definition at line 125 of file GenPlusSimParticleProducer.cc.

References reco::CompositeRefCandidateT< D >::addDaughter(), CoreSimTrack::charge(), ALCARECOTkAlJpsiMuMu_cff::charge, filter_, EgammaValidation_cff::genp, training_settings::idx, CoreSimTrack::momentum(), SimVertex::noParent(), SimTrack::noVertex(), p4, SimVertex::parentIndex(), setStatus_, CoreSimTrack::trackId(), CoreSimTrack::type(), SimTrack::vertIndex(), and badGlobalMuonTaggersAOD_cff::vtx.

Referenced by produce().

133  {
134  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
135  // Make up a GenParticleCandidate from the GEANT track info.
136  int charge = static_cast<int>(stDau.charge());
137  const Particle::LorentzVector &p4 = stDau.momentum();
138  Particle::Point vtx; // = (0,0,0) by default
139  if (!stDau.noVertex())
140  vtx = simvertices[stDau.vertIndex()].position();
141  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
142 
143  // Maybe apply filter on the particle
144  if (filter_.get() != nullptr) {
145  if (!(*filter_)(genp))
146  return;
147  }
148 
149  reco::GenParticleRef parentRef(ref, momidx);
150  genp.addMother(parentRef);
151  mergedGens.push_back(genp);
152  // get the index for the daughter just added
153  unsigned int dauidx = mergedGens.size() - 1;
154 
155  // update add daughter relationship
156  reco::GenParticle &cand = mergedGens[momidx];
157  cand.addDaughter(GenParticleRef(ref, dauidx));
158 
159  //look for simtrack daughters of stDau to see if we need to recur further down the chain
160 
161  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin(); isimtrk != simtracksSorted.end();
162  ++isimtrk) {
163  if (!isimtrk->noVertex()) {
164  // Pick the vertex (isimtrk.vertIndex() is really an index)
165  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
166 
167  // Check if the vertex has a parent track (otherwise, we're lost)
168  if (!vtx.noParent()) {
169  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
170  unsigned int idx = vtx.parentIndex();
171  SimTrackContainer::const_iterator it =
172  std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
173  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
174  if (it->trackId() == stDau.trackId()) {
175  //need the genparticle index of stDau which is dauidx
177  stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
178  }
179  }
180  }
181  }
182  }
183 }
genp
produce generated paricles in acceptance #
std::unique_ptr< StrFilter > filter_
float charge() const
charge
Definition: CoreSimTrack.cc:17
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
int parentIndex() const
Definition: SimVertex.h:29
double p4[4]
Definition: TauolaWrapper.h:92
bool noVertex() const
Definition: SimTrack.h:31
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:30
unsigned int trackId() const
Definition: CoreSimTrack.h:31
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
void addGenParticle(const SimTrack &stMom, const SimTrack &stDau, unsigned int momidx, const edm::SimTrackContainer &simtks, const edm::SimVertexContainer &simvtxs, reco::GenParticleCollection &mergedGens, const reco::GenParticleRefProd ref, std::vector< int > &genBarcodes, bool &barcodesAreSorted) const
Try to link the GEANT particle to the generator particle it came from.
bool noParent() const
Definition: SimVertex.h:30
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
void pat::GenPlusSimParticleProducer::endJob ( void  )
inlineoverrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 46 of file GenPlusSimParticleProducer.cc.

46 {}
void GenPlusSimParticleProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 185 of file GenPlusSimParticleProducer.cc.

References funct::abs(), reco::CompositeRefCandidateT< D >::addDaughter(), addGenParticle(), reco::CompositeRefCandidateT< D >::addMother(), HLT_2018_cff::cands, ztail::d, reco::CompositeRefCandidateT< D >::daughterRefVector(), DEFINE_FWK_MODULE, spr::find(), firstEvent_, genBarcodesToken_, gensToken_, mps_fire::i, edm::RefProd< C >::id(), training_settings::idx, crabWrapper::key, visualization-live-secondInstance_cfg::m, reco::CompositeRefCandidateT< D >::motherRefVector(), eostools::move(), SimVertex::noParent(), reco::CompositeRefCandidateT< D >::numberOfDaughters(), reco::CompositeRefCandidateT< D >::numberOfMothers(), SimVertex::parentIndex(), pdgIds_, pdts_, reco::CompositeRefCandidateT< D >::resetDaughters(), reco::CompositeRefCandidateT< D >::resetMothers(), simtracksToken_, simverticesToken_, and badGlobalMuonTaggersAOD_cff::vtx.

Referenced by ~GenPlusSimParticleProducer().

185  {
186  if (firstEvent_) {
187  if (!pdts_.empty()) {
188  pdgIds_.clear();
189  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
190  itp->setup(iSetup); // decode string->pdgId and vice-versa
191  pdgIds_.insert(std::abs(itp->pdgId()));
192  }
193  pdts_.clear();
194  }
195  firstEvent_ = false;
196  }
197 
198  // Simulated tracks (i.e. GEANT particles).
199  Handle<SimTrackContainer> simtracks;
200  event.getByToken(simtracksToken_, simtracks);
201 
202  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
203  std::unique_ptr<SimTrackContainer> simtracksTmp;
204  const SimTrackContainer *simtracksSorted = &*simtracks;
205  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
206  simtracksTmp.reset(new SimTrackContainer(*simtracks));
207  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
208  simtracksSorted = &*simtracksTmp;
209  }
210 
211  // Get the associated vertices
212  Handle<SimVertexContainer> simvertices;
213  event.getByToken(simverticesToken_, simvertices);
214 
215  // Get the GenParticles and barcodes, if needed to set mother and daughter links
217  Handle<std::vector<int>> genBarcodes;
218  bool barcodesAreSorted = true;
219  event.getByToken(gensToken_, gens);
220  event.getByToken(genBarcodesToken_, genBarcodes);
221  if (gens->size() != genBarcodes->size())
222  throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
223 
224  // make the output collection
225  auto candsPtr = std::make_unique<GenParticleCollection>();
226  GenParticleCollection &cands = *candsPtr;
227 
228  const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
229 
230  // add the original genParticles to the merged output list
231  for (size_t i = 0; i < gens->size(); ++i) {
232  reco::GenParticle cand((*gens)[i]);
233  cands.push_back(cand);
234  }
235 
236  // make new barcodes vector and fill it with the original list
237  auto newGenBarcodes = std::make_unique<vector<int>>();
238  for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
239  newGenBarcodes->push_back((*genBarcodes)[i]);
240  }
241  barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
242 
243  for (size_t i = 0; i < cands.size(); ++i) {
244  reco::GenParticle &cand = cands[i];
245  size_t nDaus = cand.numberOfDaughters();
247  cand.resetDaughters(ref.id());
248  for (size_t d = 0; d < nDaus; ++d) {
249  cand.addDaughter(GenParticleRef(ref, daus[d].key()));
250  }
251 
252  size_t nMoms = cand.numberOfMothers();
254  cand.resetMothers(ref.id());
255  for (size_t m = 0; m < nMoms; ++m) {
256  cand.addMother(GenParticleRef(ref, moms[m].key()));
257  }
258  }
259 
260  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin(); isimtrk != simtracks->end(); ++isimtrk) {
261  // Skip PYTHIA tracks.
262  if (isimtrk->genpartIndex() != -1)
263  continue;
264 
265  // Maybe apply the PdgId filter
266  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
267  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end())
268  continue;
269  }
270 
271  // find simtrack that has a genParticle match to its parent
272  // Look at the production vertex. If there is no vertex, I can do nothing...
273  if (!isimtrk->noVertex()) {
274  // Pick the vertex (isimtrk.vertIndex() is really an index)
275  const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
276 
277  // Check if the vertex has a parent track (otherwise, we're lost)
278  if (!vtx.noParent()) {
279  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
280  unsigned int idx = vtx.parentIndex();
281  SimTrackContainer::const_iterator it =
282  std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
283  if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
284  if (it->genpartIndex() != -1) {
285  std::vector<int>::const_iterator itIndex;
286  if (barcodesAreSorted) {
287  itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
288  } else {
289  itIndex = std::find(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
290  }
291 
292  // Check that I found something
293  // I need to check '*itIndex == it->genpartIndex()' because lower_bound just finds the right spot for an item in a sorted list, not the item
294  if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
295  // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
296  // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
297  unsigned int momidx = itIndex - genBarcodes->begin();
298  addGenParticle(*it,
299  *isimtrk,
300  momidx,
301  *simtracksSorted,
302  *simvertices,
303  *candsPtr,
304  ref,
305  *newGenBarcodes,
306  barcodesAreSorted);
307  }
308  }
309  }
310  }
311  }
312  }
313 
314  event.put(std::move(candsPtr));
315  event.put(std::move(newGenBarcodes));
316 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
size_t numberOfMothers() const override
number of mothers
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
const daughters & daughterRefVector() const
references to daughtes
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
void resetDaughters(const edm::ProductID &id)
set daughters product ID
int parentIndex() const
Definition: SimVertex.h:29
const mothers & motherRefVector() const
references to mothers
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
d
Definition: ztail.py:151
void addMother(const typename mothers::value_type &)
add a daughter via a reference
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
size_t numberOfDaughters() const override
number of daughters
edm::EDGetTokenT< reco::GenParticleCollection > gensToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
void addGenParticle(const SimTrack &stMom, const SimTrack &stDau, unsigned int momidx, const edm::SimTrackContainer &simtks, const edm::SimVertexContainer &simvtxs, reco::GenParticleCollection &mergedGens, const reco::GenParticleRefProd ref, std::vector< int > &genBarcodes, bool &barcodesAreSorted) const
Try to link the GEANT particle to the generator particle it came from.
bool noParent() const
Definition: SimVertex.h:30
void resetMothers(const edm::ProductID &id)
set mother product ID
std::vector< SimTrack > SimTrackContainer
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_

Member Data Documentation

std::unique_ptr<StrFilter> pat::GenPlusSimParticleProducer::filter_
private

Definition at line 56 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle(), and GenPlusSimParticleProducer().

bool pat::GenPlusSimParticleProducer::firstEvent_
private

Definition at line 48 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

edm::EDGetTokenT<std::vector<int> > pat::GenPlusSimParticleProducer::genBarcodesToken_
private

Definition at line 60 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

edm::EDGetTokenT<reco::GenParticleCollection> pat::GenPlusSimParticleProducer::gensToken_
private

Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of barcodes, aligned with them.

Definition at line 59 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

std::set<int> pat::GenPlusSimParticleProducer::pdgIds_
private

Definition at line 52 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

std::vector<PdtEntry> pat::GenPlusSimParticleProducer::pdts_
private

Definition at line 53 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().

int pat::GenPlusSimParticleProducer::setStatus_
private

Definition at line 51 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle().

edm::EDGetTokenT<edm::SimTrackContainer> pat::GenPlusSimParticleProducer::simtracksToken_
private

Definition at line 49 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

edm::EDGetTokenT<edm::SimVertexContainer> pat::GenPlusSimParticleProducer::simverticesToken_
private

Definition at line 50 of file GenPlusSimParticleProducer.cc.

Referenced by produce().