CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
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::GenParticle
StrFilter
 

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::GenParticleCollection
gensToken_
 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::SimTrackContainer
simtracksToken_
 
edm::EDGetTokenT
< edm::SimVertexContainer
simverticesToken_
 
edm::ESGetToken
< HepPDT::ParticleDataTable,
edm::DefaultRecord
tableToken_
 

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 wantsInputProcessBlocks ()
 
static bool wantsProcessBlocks ()
 
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)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< B > consumes (edm::InputTag tag) noexcept
 
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<Transition Tr = Transition::Event>
constexpr auto esConsumes () noexcept
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag) noexcept
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

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 42 of file GenPlusSimParticleProducer.cc.

Member Typedef Documentation

Definition at line 58 of file GenPlusSimParticleProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 102 of file GenPlusSimParticleProducer.cc.

References edm::EDConsumerBase::esConsumes(), edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), alcazmumu_cfi::filter, filter_, edm::ParameterSet::getParameter(), pdts_, and tableToken_.

103  : firstEvent_(true),
104  simtracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
106  consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
107  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
108  gensToken_(consumes<GenParticleCollection>(
109  cfg.getParameter<InputTag>("genParticles"))), // get the genParticles to add GEANT particles to
110  genBarcodesToken_(consumes<std::vector<int>>(
111  cfg.getParameter<InputTag>("genParticles"))) // get the genParticles to add GEANT particles to
112 {
113  // Possibly allow a list of particle types
114  if (cfg.exists("particleTypes")) {
115  pdts_ = cfg.getParameter<vector<PdtEntry>>("particleTypes");
116  if (!pdts_.empty()) {
118  }
119  }
120 
121  // Possibly allow a string cut
122  if (cfg.existsAs<string>("filter")) {
123  string filter = cfg.getParameter<string>("filter");
124  if (!filter.empty()) {
125  filter_ = std::make_unique<StrFilter>(filter);
126  }
127  }
128  produces<GenParticleCollection>();
129  produces<vector<int>>();
130 }
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
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&lt;int&gt; of ba...
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_
pat::GenPlusSimParticleProducer::~GenPlusSimParticleProducer ( )
inlineoverride

Definition at line 45 of file GenPlusSimParticleProducer.cc.

45 {}

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 132 of file GenPlusSimParticleProducer.cc.

References reco::CompositeRefCandidateT< D >::addDaughter(), CoreSimTrack::charge(), RecoTauCleanerPlugins::charge, filter_, EgammaValidation_cff::genp, cuda_std::lower_bound(), CoreSimTrack::momentum(), SimVertex::noParent(), SimTrack::noVertex(), SimVertex::parentIndex(), setStatus_, CoreSimTrack::trackId(), CoreSimTrack::type(), and SimTrack::vertIndex().

Referenced by produce().

140  {
141  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
142  // Make up a GenParticleCandidate from the GEANT track info.
143  int charge = static_cast<int>(stDau.charge());
144  const Particle::LorentzVector &p4 = stDau.momentum();
145  Particle::Point vtx; // = (0,0,0) by default
146  if (!stDau.noVertex())
147  vtx = simvertices[stDau.vertIndex()].position();
148  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
149 
150  // Maybe apply filter on the particle
151  if (filter_.get() != nullptr) {
152  if (!(*filter_)(genp))
153  return;
154  }
155 
156  reco::GenParticleRef parentRef(ref, momidx);
157  genp.addMother(parentRef);
158  mergedGens.push_back(genp);
159  // get the index for the daughter just added
160  unsigned int dauidx = mergedGens.size() - 1;
161 
162  // update add daughter relationship
163  reco::GenParticle &cand = mergedGens[momidx];
164  cand.addDaughter(GenParticleRef(ref, dauidx));
165 
166  //look for simtrack daughters of stDau to see if we need to recur further down the chain
167 
168  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin(); isimtrk != simtracksSorted.end();
169  ++isimtrk) {
170  if (!isimtrk->noVertex()) {
171  // Pick the vertex (isimtrk.vertIndex() is really an index)
172  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
173 
174  // Check if the vertex has a parent track (otherwise, we're lost)
175  if (!vtx.noParent()) {
176  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
177  unsigned int idx = vtx.parentIndex();
178  SimTrackContainer::const_iterator it =
179  std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
180  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
181  if (it->trackId() == stDau.trackId()) {
182  //need the genparticle index of stDau which is dauidx
184  stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
185  }
186  }
187  }
188  }
189  }
190 }
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
bool noVertex() const
Definition: SimTrack.h:34
tuple genp
produce generated paricles in acceptance #
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:33
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
__host__ __device__ constexpr RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
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 49 of file GenPlusSimParticleProducer.cc.

49 {}
void GenPlusSimParticleProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 192 of file GenPlusSimParticleProducer.cc.

References funct::abs(), reco::CompositeRefCandidateT< D >::addDaughter(), addGenParticle(), reco::CompositeRefCandidateT< D >::addMother(), HLT_FULL_cff::cands, ztail::d, reco::CompositeRefCandidateT< D >::daughterRefVector(), spr::find(), firstEvent_, genBarcodesToken_, gensToken_, edm::EventSetup::getData(), mps_fire::i, edm::RefProd< T >::id(), submitPVResolutionJobs::key, cuda_std::lower_bound(), 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 tableToken_.

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

Member Data Documentation

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

Definition at line 59 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle(), and GenPlusSimParticleProducer().

bool pat::GenPlusSimParticleProducer::firstEvent_
private

Definition at line 51 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 63 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 62 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 55 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 56 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().

int pat::GenPlusSimParticleProducer::setStatus_
private

Definition at line 54 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle().

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

Definition at line 52 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 53 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> pat::GenPlusSimParticleProducer::tableToken_
private

Definition at line 65 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().