CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

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

Public Member Functions

 PATGenCandsFromSimTracksProducer (const edm::ParameterSet &)
 
 ~PATGenCandsFromSimTracksProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef
StringCutObjectSelector
< reco::GenParticle
StrFilter
 

Private Member Functions

virtual void endJob () override
 
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::GenParticleCollection
findRef (const SimTrack &tk, GlobalContext &g) const
 
edm::Ref
< reco::GenParticleCollection
generatorRef_ (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...
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::auto_ptr< StrFilterfilter_
 
bool firstEvent_
 
edm::EDGetTokenT< std::vector
< int > > 
genBarcodesToken_
 
edm::EDGetTokenT
< reco::GenParticleCollection
genParticlesToken_
 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::SimTrackContainer
simTracksToken_
 
edm::EDGetTokenT
< edm::SimVertexContainer
simVertexToken_
 
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::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- 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 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 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 47 of file PATGenCandsFromSimTracksProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 118 of file PATGenCandsFromSimTracksProducer.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), alcazmumu_cfi::filter, filter_, edm::ParameterSet::getParameter(), makeMotherLink_, motherPdts_, pdts_, and writeAncestors_.

118  :
119  firstEvent_(true),
120  simTracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks
121  simVertexToken_(consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim vertices
122  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
123  makeMotherLink_(cfg.existsAs<bool>("makeMotherLink") ? cfg.getParameter<bool>("makeMotherLink") : false),
124  writeAncestors_(cfg.existsAs<bool>("writeAncestors") ? cfg.getParameter<bool>("writeAncestors") : false),
125  genParticlesToken_(mayConsume<GenParticleCollection>(cfg.getParameter<InputTag>("genParticles"))),
126  genBarcodesToken_(mayConsume<std::vector<int> >(cfg.getParameter<InputTag>("genParticles")))
127 {
128  // Possibly allow a list of particle types
129  if (cfg.exists("particleTypes")) {
130  pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
131  }
132  if (cfg.exists("motherTypes")) {
133  motherPdts_ = cfg.getParameter<vector<PdtEntry> >("motherTypes");
134  }
135 
136  // Possibly allow a string cut
137  if (cfg.existsAs<string>("filter")) {
138  string filter = cfg.getParameter<string>("filter");
139  if (!filter.empty()) {
140  filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
141  }
142  }
143 
145  edm::LogWarning("Configuration") << "PATGenCandsFromSimTracksProducer: " <<
146  "you have set 'writeAncestors' to 'true' and 'makeMotherLink' to false;" <<
147  "GEANT particles with generator level (e.g.PYHIA) mothers won't have mother links.\n";
148  }
149  produces<GenParticleCollection>();
150 }
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:184
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector&lt;int&gt; 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.
StringCutObjectSelector< reco::GenParticle > StrFilter
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken_
edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
pat::PATGenCandsFromSimTracksProducer::~PATGenCandsFromSimTracksProducer ( )
inline

Definition at line 32 of file PATGenCandsFromSimTracksProducer.cc.

32 {}

Member Function Documentation

virtual void pat::PATGenCandsFromSimTracksProducer::endJob ( void  )
inlineoverrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 36 of file PATGenCandsFromSimTracksProducer.cc.

36 {}
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 153 of file PATGenCandsFromSimTracksProducer.cc.

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

Referenced by findRef(), and produce().

153  {
154  assert(tk.genpartIndex() == -1); // MUST NOT be called with a PYTHIA track
155  if (!tk.noVertex()) {
156  const SimVertex &vtx = g.simvtxs[tk.vertIndex()];
157  if (!vtx.noParent()) {
158  unsigned int idx = vtx.parentIndex();
159  SimTrackContainer::const_iterator it = std::lower_bound(g.simtks.begin(), g.simtks.end(), idx, LessById());
160  if ((it != g.simtks.end()) && (it->trackId() == idx)) {
161  return &*it;
162  }
163  }
164  }
165  return 0;
166 }
int parentIndex() const
Definition: SimVertex.h:33
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
bool noVertex() const
Definition: SimTrack.h:30
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
bool noParent() const
Definition: SimVertex.h:34
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 169 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().

169  {
171  const SimTrack * simMother = findGeantMother(tk, g);
172 
174  if (simMother != 0) motherRef = findRef(*simMother,g);
175 
176  if (writeAncestors_) {
177  // If writing ancestors, I need to serialize myself, and then to return a ref to me
178  // But first check if I've already been serialized
179  std::map<unsigned int,int>::const_iterator it = g.simTksProcessed.find(tk.trackId());
180  if (it != g.simTksProcessed.end()) {
181  // just return a ref to it
182  assert(it->second > 0);
183  return edm::Ref<reco::GenParticleCollection>(g.refprod, (it->second) - 1);
184  } else {
185  // make genParticle, save, update the map, and return ref to myself
186  reco::GenParticle p = makeGenParticle_(tk, motherRef, g);
187  g.output.push_back(p);
188  g.simTksProcessed[tk.trackId()] = g.output.size();
189  return edm::Ref<reco::GenParticleCollection>(g.refprod, g.output.size()-1 );
190  }
191  } else {
192  // Otherwise, I just return a ref to my mum
193  return motherRef;
194  }
195 }
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:33
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:34
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 198 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().

198  {
199  assert(st.genpartIndex() != -1);
200  // Note that st.genpartIndex() is the barcode, not the index within GenParticleCollection, so I have to search the particle
201  std::vector<int>::const_iterator it;
202  if (g.barcodesAreSorted) {
203  it = std::lower_bound(g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
204  } else {
205  it = std::find( g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
206  }
207 
208  // Check that I found something
209  // 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
210  if ((it != g.genBarcodes->end()) && (*it == st.genpartIndex())) {
211  return reco::GenParticleRef(g.gens, it - g.genBarcodes->begin());
212  } else {
213  return reco::GenParticleRef();
214  }
215 }
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:7
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 218 of file PATGenCandsFromSimTracksProducer.cc.

References reco::CompositeRefCandidateT< D >::addMother(), DeDxDiscriminatorTools::charge(), CoreSimTrack::charge(), edm::Ref< C, T, F >::isNonnull(), CoreSimTrack::momentum(), SimTrack::noVertex(), p4, setStatus_, pat::PATGenCandsFromSimTracksProducer::GlobalContext::simvtxs, CoreSimTrack::type(), and SimTrack::vertIndex().

Referenced by findRef(), and produce().

218  {
219  // Make up a GenParticleCandidate from the GEANT track info.
220  int charge = static_cast<int>(tk.charge());
222  Particle::Point vtx; // = (0,0,0) by default
223  if (!tk.noVertex()) vtx = g.simvtxs[tk.vertIndex()].position();
224  GenParticle gp(charge, p4, vtx, tk.type(), setStatus_, true);
225  if (mother.isNonnull()) gp.addMother(mother);
226  return gp;
227 }
float charge() const
charge
Definition: CoreSimTrack.cc:18
double charge(const std::vector< uint8_t > &Ampls)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
double p4[4]
Definition: TauolaWrapper.h:92
bool noVertex() const
Definition: SimTrack.h:30
void addMother(const typename mothers::value_type &)
add a daughter via a reference
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:28
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:32
void PATGenCandsFromSimTracksProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 230 of file PATGenCandsFromSimTracksProducer.cc.

References funct::abs(), reco::CompositeRefCandidateT< D >::addMother(), gather_cfg::cout, edm::hlt::Exception, filter_, findGeantMother(), findRef(), firstEvent_, genBarcodesToken_, EgammaValidation_cff::genp, genParticlesToken_, 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_, pdgIds_, pdts_, simTracksToken_, simVertexToken_, python.multivaluedict::sort(), CoreSimTrack::type(), and writeAncestors_.

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

Member Data Documentation

std::auto_ptr<StrFilter> pat::PATGenCandsFromSimTracksProducer::filter_
private
bool pat::PATGenCandsFromSimTracksProducer::firstEvent_
private

Definition at line 38 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

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

Definition at line 57 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 56 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 51 of file PATGenCandsFromSimTracksProducer.cc.

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

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

Definition at line 44 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

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

Definition at line 42 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

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

Definition at line 41 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by makeGenParticle_().

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

Definition at line 39 of file PATGenCandsFromSimTracksProducer.cc.

Referenced by produce().

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

Definition at line 40 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 53 of file PATGenCandsFromSimTracksProducer.cc.

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