CMS 3D CMS Logo

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

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

Inheritance diagram for pat::PATGenCandsFromSimTracksProducer:
edm::stream::EDProducer<>

Classes

struct  GlobalContext
 Global context for all recursive methods. More...
 
struct  LessById
 

Public Member Functions

 PATGenCandsFromSimTracksProducer (const edm::ParameterSet &)
 
 ~PATGenCandsFromSimTracksProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Private Types

typedef StringCutObjectSelector< reco::GenParticleStrFilter
 

Private Member Functions

const SimTrackfindGeantMother (const SimTrack &tk, const GlobalContext &g) const
 Find the mother of a given GEANT track (or NULL if it can't be found). More...
 
edm::Ref< reco::GenParticleCollectionfindRef (const SimTrack &tk, GlobalContext &g) const
 
edm::Ref< reco::GenParticleCollectiongeneratorRef_ (const SimTrack &tk, const GlobalContext &g) const
 Used by findRef if the track is a PYTHIA particle. More...
 
reco::GenParticle makeGenParticle_ (const SimTrack &tk, const edm::Ref< reco::GenParticleCollection > &mother, const GlobalContext &g) const
 Make a GenParticle for this SimTrack, with a given mother. More...
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const StrFilter filter_
 
bool firstEvent_
 
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
 
edm::EDGetTokenT< reco::GenParticleCollectiongenParticlesToken_
 Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of barcodes, aligned with them. More...
 
bool makeMotherLink_
 If true, I'll try to make a link from the GEANT particle to a GenParticle. More...
 
std::set< int > motherPdgIds_
 
std::vector< PdtEntrymotherPdts_
 
std::set< int > pdgIds_
 
std::vector< PdtEntrypdts_
 
int setStatus_
 
edm::EDGetTokenT< edm::SimTrackContainersimTracksToken_
 
edm::EDGetTokenT< edm::SimVertexContainersimVertexToken_
 
bool writeAncestors_
 If true, I'll save GenParticles corresponding to the ancestors of this GEANT particle. Common ancestors are only written once. More...
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Produces reco::GenParticle from SimTracks.

The PATGenCandsFromSimTracksProducer produces GenParticles from SimTracks, so they can be used for MC matching.

Author
Jordan Tucker (original module), Giovanni Petrucciani (PAT integration)
Version
Id
PATGenCandsFromSimTracksProducer.cc,v 1.8 2010/10/20 23:09:25 wmtan Exp

Definition at line 29 of file PATGenCandsFromSimTracksProducer.cc.

Member Typedef Documentation

Definition at line 46 of file PATGenCandsFromSimTracksProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 121 of file PATGenCandsFromSimTracksProducer.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), makeMotherLink_, motherPdts_, pdts_, and writeAncestors_.

122  : firstEvent_(true),
123  simTracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks
124  simVertexToken_(consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim vertices
125  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
126  filter_(cfg.existsAs<string>("filter") ? cfg.getParameter<string>("filter") : std::string("1 == 1")),
127  makeMotherLink_(cfg.existsAs<bool>("makeMotherLink") ? cfg.getParameter<bool>("makeMotherLink") : false),
128  writeAncestors_(cfg.existsAs<bool>("writeAncestors") ? cfg.getParameter<bool>("writeAncestors") : false),
129  genParticlesToken_(mayConsume<GenParticleCollection>(cfg.getParameter<InputTag>("genParticles"))),
130  genBarcodesToken_(mayConsume<std::vector<int> >(cfg.getParameter<InputTag>("genParticles"))) {
131  // Possibly allow a list of particle types
132  if (cfg.exists("particleTypes")) {
133  pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
134  }
135  if (cfg.exists("motherTypes")) {
136  motherPdts_ = cfg.getParameter<vector<PdtEntry> >("motherTypes");
137  }
138 
140  edm::LogWarning("Configuration")
141  << "PATGenCandsFromSimTracksProducer: "
142  << "you have set 'writeAncestors' to 'true' and 'makeMotherLink' to false;"
143  << "GEANT particles with generator level (e.g.PYHIA) mothers won't have mother links.\n";
144  }
145  produces<GenParticleCollection>();
146 }
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< reco::GenParticleCollection > genParticlesToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool writeAncestors_
If true, I&#39;ll save GenParticles corresponding to the ancestors of this GEANT particle. Common ancestors are only written once.
bool makeMotherLink_
If true, I&#39;ll try to make a link from the GEANT particle to a GenParticle.
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken_
edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
pat::PATGenCandsFromSimTracksProducer::~PATGenCandsFromSimTracksProducer ( )
inlineoverride

Definition at line 32 of file PATGenCandsFromSimTracksProducer.cc.

References produce().

32 {}

Member Function Documentation

const SimTrack * PATGenCandsFromSimTracksProducer::findGeantMother ( const SimTrack tk,
const GlobalContext g 
) const
private

Find the mother of a given GEANT track (or NULL if it can't be found).

Definition at line 148 of file PATGenCandsFromSimTracksProducer.cc.

References SimTrack::genpartIndex(), training_settings::idx, SimVertex::noParent(), SimTrack::noVertex(), SimVertex::parentIndex(), pat::PATGenCandsFromSimTracksProducer::GlobalContext::simtks, pat::PATGenCandsFromSimTracksProducer::GlobalContext::simvtxs, SimTrack::vertIndex(), and badGlobalMuonTaggersAOD_cff::vtx.

Referenced by findRef(), and produce().

148  {
149  assert(tk.genpartIndex() == -1); // MUST NOT be called with a PYTHIA track
150  if (!tk.noVertex()) {
151  const SimVertex &vtx = g.simvtxs[tk.vertIndex()];
152  if (!vtx.noParent()) {
153  unsigned int idx = vtx.parentIndex();
154  SimTrackContainer::const_iterator it = std::lower_bound(g.simtks.begin(), g.simtks.end(), idx, LessById());
155  if ((it != g.simtks.end()) && (it->trackId() == idx)) {
156  return &*it;
157  }
158  }
159  }
160  return nullptr;
161 }
int parentIndex() const
Definition: SimVertex.h:29
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:34
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
bool noParent() const
Definition: SimVertex.h:30
edm::Ref< reco::GenParticleCollection > PATGenCandsFromSimTracksProducer::findRef ( const SimTrack tk,
GlobalContext g 
) const
private

Find the GenParticle reference for a given GEANT or PYTHIA track.

  • if the track corresponds to a PYTHIA particle, return a ref to that particle
  • otherwise, if this simtrack has no mother simtrack, return a null ref
  • otherwise, if writeAncestors is true, make a GenParticle for it and return a ref to it
  • otherwise, if writeAncestors is false, return the ref to the GEANT mother of this track

Definition at line 163 of file PATGenCandsFromSimTracksProducer.cc.

References findGeantMother(), generatorRef_(), SimTrack::genpartIndex(), makeGenParticle_(), makeMotherLink_, pat::PATGenCandsFromSimTracksProducer::GlobalContext::output, AlCaHLTBitMon_ParallelJobs::p, pat::PATGenCandsFromSimTracksProducer::GlobalContext::refprod, pat::PATGenCandsFromSimTracksProducer::GlobalContext::simTksProcessed, CoreSimTrack::trackId(), and writeAncestors_.

Referenced by produce().

164  {
165  if (tk.genpartIndex() != -1)
167  const SimTrack *simMother = findGeantMother(tk, g);
168 
170  if (simMother != nullptr)
171  motherRef = findRef(*simMother, g);
172 
173  if (writeAncestors_) {
174  // If writing ancestors, I need to serialize myself, and then to return a ref to me
175  // But first check if I've already been serialized
176  std::map<unsigned int, int>::const_iterator it = g.simTksProcessed.find(tk.trackId());
177  if (it != g.simTksProcessed.end()) {
178  // just return a ref to it
179  assert(it->second > 0);
180  return edm::Ref<reco::GenParticleCollection>(g.refprod, (it->second) - 1);
181  } else {
182  // make genParticle, save, update the map, and return ref to myself
183  reco::GenParticle p = makeGenParticle_(tk, motherRef, g);
184  g.output.push_back(p);
185  g.simTksProcessed[tk.trackId()] = g.output.size();
186  return edm::Ref<reco::GenParticleCollection>(g.refprod, g.output.size() - 1);
187  }
188  } else {
189  // Otherwise, I just return a ref to my mum
190  return motherRef;
191  }
192 }
const SimTrack * findGeantMother(const SimTrack &tk, const GlobalContext &g) const
Find the mother of a given GEANT track (or NULL if it can&#39;t be found).
edm::Ref< reco::GenParticleCollection > findRef(const SimTrack &tk, GlobalContext &g) const
bool writeAncestors_
If true, I&#39;ll save GenParticles corresponding to the ancestors of this GEANT particle. Common ancestors are only written once.
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:34
bool makeMotherLink_
If true, I&#39;ll try to make a link from the GEANT particle to a GenParticle.
reco::GenParticle makeGenParticle_(const SimTrack &tk, const edm::Ref< reco::GenParticleCollection > &mother, const GlobalContext &g) const
Make a GenParticle for this SimTrack, with a given mother.
unsigned int trackId() const
Definition: CoreSimTrack.h:31
edm::Ref< reco::GenParticleCollection > generatorRef_(const SimTrack &tk, const GlobalContext &g) const
Used by findRef if the track is a PYTHIA particle.
edm::Ref< reco::GenParticleCollection > PATGenCandsFromSimTracksProducer::generatorRef_ ( const SimTrack tk,
const GlobalContext g 
) const
private

Used by findRef if the track is a PYTHIA particle.

Definition at line 194 of file PATGenCandsFromSimTracksProducer.cc.

References pat::PATGenCandsFromSimTracksProducer::GlobalContext::barcodesAreSorted, spr::find(), pat::PATGenCandsFromSimTracksProducer::GlobalContext::genBarcodes, SimTrack::genpartIndex(), and pat::PATGenCandsFromSimTracksProducer::GlobalContext::gens.

Referenced by findRef().

195  {
196  assert(st.genpartIndex() != -1);
197  // Note that st.genpartIndex() is the barcode, not the index within GenParticleCollection, so I have to search the particle
198  std::vector<int>::const_iterator it;
199  if (g.barcodesAreSorted) {
200  it = std::lower_bound(g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
201  } else {
202  it = std::find(g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
203  }
204 
205  // Check that I found something
206  // I need to check '*it == st.genpartIndex()' because lower_bound just finds the right spot for an item in a sorted list, not the item
207  if ((it != g.genBarcodes->end()) && (*it == st.genpartIndex())) {
208  return reco::GenParticleRef(g.gens, it - g.genBarcodes->begin());
209  } else {
210  return reco::GenParticleRef();
211  }
212 }
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
reco::GenParticle PATGenCandsFromSimTracksProducer::makeGenParticle_ ( const SimTrack tk,
const edm::Ref< reco::GenParticleCollection > &  mother,
const GlobalContext g 
) const
private

Make a GenParticle for this SimTrack, with a given mother.

Definition at line 214 of file PATGenCandsFromSimTracksProducer.cc.

References CoreSimTrack::charge(), ALCARECOTkAlJpsiMuMu_cff::charge, runTauDisplay::gp, edm::Ref< C, T, F >::isNonnull(), CoreSimTrack::momentum(), SimTrack::noVertex(), p4, setStatus_, pat::PATGenCandsFromSimTracksProducer::GlobalContext::simvtxs, CoreSimTrack::type(), SimTrack::vertIndex(), and badGlobalMuonTaggersAOD_cff::vtx.

Referenced by findRef(), and produce().

215  {
216  // Make up a GenParticleCandidate from the GEANT track info.
217  int charge = static_cast<int>(tk.charge());
218  const Particle::LorentzVector &p4 = tk.momentum();
219  Particle::Point vtx; // = (0,0,0) by default
220  if (!tk.noVertex())
221  vtx = g.simvtxs[tk.vertIndex()].position();
222  GenParticle gp(charge, p4, vtx, tk.type(), setStatus_, true);
223  if (mother.isNonnull())
224  gp.addMother(mother);
225  return gp;
226 }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
float charge() const
charge
Definition: CoreSimTrack.cc:17
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
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
void PATGenCandsFromSimTracksProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 228 of file PATGenCandsFromSimTracksProducer.cc.

References funct::abs(), reco::CompositeRefCandidateT< D >::addMother(), HLT_2018_cff::cands, gather_cfg::cout, DEFINE_FWK_MODULE, LEDCalibrationChannels::depth, Exception, filter_, findGeantMother(), findRef(), firstEvent_, genBarcodesToken_, EgammaValidation_cff::genp, genParticlesToken_, mps_fire::i, edm::OrphanHandleBase::id(), edm::Ref< C, T, F >::id(), edm::Ref< C, T, F >::isNonnull(), edm::Ref< C, T, F >::key(), makeGenParticle_(), makeMotherLink_, motherPdgIds_, motherPdts_, eostools::move(), pdgIds_, pdts_, simTracksToken_, simVertexToken_, CoreSimTrack::type(), and writeAncestors_.

Referenced by ~PATGenCandsFromSimTracksProducer().

228  {
229  if (firstEvent_) {
230  if (!pdts_.empty()) {
231  pdgIds_.clear();
232  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
233  itp->setup(iSetup); // decode string->pdgId and vice-versa
234  pdgIds_.insert(std::abs(itp->pdgId()));
235  }
236  pdts_.clear();
237  }
238  if (!motherPdts_.empty()) {
239  motherPdgIds_.clear();
240  for (vector<PdtEntry>::iterator itp = motherPdts_.begin(), edp = motherPdts_.end(); itp != edp; ++itp) {
241  itp->setup(iSetup); // decode string->pdgId and vice-versa
242  motherPdgIds_.insert(std::abs(itp->pdgId()));
243  }
244  motherPdts_.clear();
245  }
246  firstEvent_ = false;
247  }
248 
249  // Simulated tracks (i.e. GEANT particles).
250  Handle<SimTrackContainer> simtracks;
251  event.getByToken(simTracksToken_, simtracks);
252 
253  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
254  std::unique_ptr<SimTrackContainer> simtracksTmp;
255  const SimTrackContainer *simtracksSorted = &*simtracks;
257  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
258  simtracksTmp.reset(new SimTrackContainer(*simtracks));
259  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
260  simtracksSorted = &*simtracksTmp;
261  }
262  }
263 
264  // Get the associated vertices
265  Handle<SimVertexContainer> simvertices;
266  event.getByToken(simVertexToken_, simvertices);
267 
268  // Get the GenParticles and barcodes, if needed to set mother links
270  Handle<std::vector<int> > genBarcodes;
271  bool barcodesAreSorted = true;
272  if (makeMotherLink_) {
273  event.getByToken(genParticlesToken_, gens);
274  event.getByToken(genBarcodesToken_, genBarcodes);
275  if (gens->size() != genBarcodes->size())
276  throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
277  barcodesAreSorted = __gnu_cxx::is_sorted(genBarcodes->begin(), genBarcodes->end());
278  }
279 
280  // make the output collection
281  auto cands = std::make_unique<GenParticleCollection>();
282  edm::RefProd<GenParticleCollection> refprod = event.getRefBeforePut<GenParticleCollection>();
283 
284  GlobalContext globals(*simtracksSorted, *simvertices, gens, genBarcodes, barcodesAreSorted, *cands, refprod);
285 
286  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin(); isimtrk != simtracks->end(); ++isimtrk) {
287  // Skip PYTHIA tracks.
288  if (isimtrk->genpartIndex() != -1)
289  continue;
290 
291  // Maybe apply the PdgId filter
292  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
293  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end())
294  continue;
295  }
296 
298 
299  // Maybe apply filter on the particle
300  if (!(filter_(genp)))
301  continue;
302 
303  if (!motherPdgIds_.empty()) {
304  const SimTrack *motherSimTk = findGeantMother(*isimtrk, globals);
305  if (motherSimTk == nullptr)
306  continue;
307  if (motherPdgIds_.find(std::abs(motherSimTk->type())) == motherPdgIds_.end())
308  continue;
309  }
310 
312  Ref<GenParticleCollection> motherRef;
313  const SimTrack *mother = findGeantMother(*isimtrk, globals);
314  if (mother != nullptr)
315  motherRef = findRef(*mother, globals);
316  if (motherRef.isNonnull())
317  genp.addMother(motherRef);
318  }
319 
320  cands->push_back(genp);
321  }
322 
323  // Write to the Event, and get back a handle (which can be useful for debugging)
324  edm::OrphanHandle<reco::GenParticleCollection> orphans = event.put(std::move(cands));
325 
326 #ifdef DEBUG_PATGenCandsFromSimTracksProducer
327  std::cout << "Produced a list of " << orphans->size() << " genParticles." << std::endl;
328  for (GenParticleCollection::const_iterator it = orphans->begin(), ed = orphans->end(); it != ed; ++it) {
329  std::cout << " ";
330  std::cout << "GenParticle #" << (it - orphans->begin()) << ": pdgId " << it->pdgId() << ", pt = " << it->pt()
331  << ", eta = " << it->eta() << ", phi = " << it->phi() << ", rho = " << it->vertex().Rho()
332  << ", z = " << it->vertex().Z() << std::endl;
333  edm::Ref<GenParticleCollection> mom = it->motherRef();
334  size_t depth = 2;
335  while (mom.isNonnull()) {
336  if (mom.id() == orphans.id()) {
337  // I need to re-make the ref because they are not working until this module returns.
338  mom = edm::Ref<GenParticleCollection>(orphans, mom.key());
339  }
340  for (size_t i = 0; i < depth; ++i)
341  std::cout << " ";
342  std::cout << "GenParticleRef [" << mom.id() << "/" << mom.key() << "]: pdgId " << mom->pdgId()
343  << ", status = " << mom->status() << ", pt = " << mom->pt() << ", eta = " << mom->eta()
344  << ", phi = " << mom->phi() << ", rho = " << mom->vertex().Rho() << ", z = " << mom->vertex().Z()
345  << std::endl;
346  if (mom.id() != orphans.id())
347  break;
348  if ((mom->motherRef().id() == mom.id()) && (mom->motherRef().key() == mom.key())) {
349  throw cms::Exception("Corrupt Data") << "A particle is it's own mother.\n";
350  }
351  mom = mom->motherRef();
352  depth++;
353  }
354  }
355  std::cout << std::endl;
356 #endif
357 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
genp
produce generated paricles in acceptance #
const SimTrack * findGeantMother(const SimTrack &tk, const GlobalContext &g) const
Find the mother of a given GEANT track (or NULL if it can&#39;t be found).
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
edm::Ref< reco::GenParticleCollection > findRef(const SimTrack &tk, GlobalContext &g) const
bool writeAncestors_
If true, I&#39;ll save GenParticles corresponding to the ancestors of this GEANT particle. Common ancestors are only written once.
key_type key() const
Accessor for product key.
Definition: Ref.h:250
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
ProductID id() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool makeMotherLink_
If true, I&#39;ll try to make a link from the GEANT particle to a GenParticle.
reco::GenParticle makeGenParticle_(const SimTrack &tk, const edm::Ref< reco::GenParticleCollection > &mother, const GlobalContext &g) const
Make a GenParticle for this SimTrack, with a given mother.
void addMother(const typename mothers::value_type &)
add a daughter via a reference
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken_
edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
std::vector< SimTrack > SimTrackContainer
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

const StrFilter pat::PATGenCandsFromSimTracksProducer::filter_
private

Definition at line 47 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

bool pat::PATGenCandsFromSimTracksProducer::firstEvent_
private

Definition at line 37 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

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

Definition at line 56 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

edm::EDGetTokenT<reco::GenParticleCollection> pat::PATGenCandsFromSimTracksProducer::genParticlesToken_
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 55 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

bool pat::PATGenCandsFromSimTracksProducer::makeMotherLink_
private

If true, I'll try to make a link from the GEANT particle to a GenParticle.

Definition at line 50 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by findRef(), PATGenCandsFromSimTracksProducer(), and produce().

std::set<int> pat::PATGenCandsFromSimTracksProducer::motherPdgIds_
private

Definition at line 43 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

std::vector<PdtEntry> pat::PATGenCandsFromSimTracksProducer::motherPdts_
private
std::set<int> pat::PATGenCandsFromSimTracksProducer::pdgIds_
private

Definition at line 41 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

std::vector<PdtEntry> pat::PATGenCandsFromSimTracksProducer::pdts_
private
int pat::PATGenCandsFromSimTracksProducer::setStatus_
private

Definition at line 40 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by makeGenParticle_().

edm::EDGetTokenT<edm::SimTrackContainer> pat::PATGenCandsFromSimTracksProducer::simTracksToken_
private

Definition at line 38 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

edm::EDGetTokenT<edm::SimVertexContainer> pat::PATGenCandsFromSimTracksProducer::simVertexToken_
private

Definition at line 39 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

bool pat::PATGenCandsFromSimTracksProducer::writeAncestors_
private

If true, I'll save GenParticles corresponding to the ancestors of this GEANT particle. Common ancestors are only written once.

Definition at line 52 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by findRef(), PATGenCandsFromSimTracksProducer(), and produce().