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 ()
 
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 ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef
StringCutObjectSelector
< reco::GenParticle
StrFilter
 

Private Member Functions

virtual void endJob ()
 
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::InputTag genParticles_
 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::InputTag src_
 
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
 
typedef WorkerT< EDProducerWorkerType
 
- 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::EDProducer
CurrentProcessingContext const * currentContext () const
 
- 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.9 2013/02/27 23:26:56 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 116 of file PATGenCandsFromSimTracksProducer.cc.

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

116  :
117  firstEvent_(true),
118  src_(cfg.getParameter<InputTag>("src")), // source sim tracks & vertices
119  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
120  makeMotherLink_(cfg.existsAs<bool>("makeMotherLink") ? cfg.getParameter<bool>("makeMotherLink") : false),
121  writeAncestors_(cfg.existsAs<bool>("writeAncestors") ? cfg.getParameter<bool>("writeAncestors") : false),
123 {
124  // Possibly allow a list of particle types
125  if (cfg.exists("particleTypes")) {
126  pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
127  }
128  if (cfg.exists("motherTypes")) {
129  motherPdts_ = cfg.getParameter<vector<PdtEntry> >("motherTypes");
130  }
131 
132  // Possibly allow a string cut
133  if (cfg.existsAs<string>("filter")) {
134  string filter = cfg.getParameter<string>("filter");
135  if (!filter.empty()) {
136  filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
137  }
138  }
139 
141  edm::LogWarning("Configuration") << "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:187
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::InputTag genParticles_
Collection of GenParticles I need to make refs to. It must also have its associated vector&lt;int&gt; of ba...
pat::PATGenCandsFromSimTracksProducer::~PATGenCandsFromSimTracksProducer ( )
inline

Definition at line 32 of file PATGenCandsFromSimTracksProducer.cc.

32 {}

Member Function Documentation

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

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 149 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().

149  {
150  assert(tk.genpartIndex() == -1); // MUST NOT be called with a PYTHIA track
151  if (!tk.noVertex()) {
152  const SimVertex &vtx = g.simvtxs[tk.vertIndex()];
153  if (!vtx.noParent()) {
154  unsigned int idx = vtx.parentIndex();
155  SimTrackContainer::const_iterator it = std::lower_bound(g.simtks.begin(), g.simtks.end(), idx, LessById());
156  if ((it != g.simtks.end()) && (it->trackId() == idx)) {
157  return &*it;
158  }
159  }
160  }
161  return 0;
162 }
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
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 165 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().

165  {
167  const SimTrack * simMother = findGeantMother(tk, g);
168 
170  if (simMother != 0) motherRef = findRef(*simMother,g);
171 
172  if (writeAncestors_) {
173  // If writing ancestors, I need to serialize myself, and then to return a ref to me
174  // But first check if I've already been serialized
175  std::map<unsigned int,int>::const_iterator it = g.simTksProcessed.find(tk.trackId());
176  if (it != g.simTksProcessed.end()) {
177  // just return a ref to it
178  assert(it->second > 0);
179  return edm::Ref<reco::GenParticleCollection>(g.refprod, (it->second) - 1);
180  } else {
181  // make genParticle, save, update the map, and return ref to myself
182  reco::GenParticle p = makeGenParticle_(tk, motherRef, g);
183  g.output.push_back(p);
184  g.simTksProcessed[tk.trackId()] = g.output.size();
185  return edm::Ref<reco::GenParticleCollection>(g.refprod, g.output.size()-1 );
186  }
187  } else {
188  // Otherwise, I just return a ref to my mum
189  return motherRef;
190  }
191 }
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.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
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:49
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().

194  {
195  assert(st.genpartIndex() != -1);
196  // Note that st.genpartIndex() is the barcode, not the index within GenParticleCollection, so I have to search the particle
197  std::vector<int>::const_iterator it;
198  if (g.barcodesAreSorted) {
199  it = std::lower_bound(g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
200  } else {
201  it = std::find( g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
202  }
203 
204  // Check that I found something
205  // 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
206  if ((it != g.genBarcodes->end()) && (*it == st.genpartIndex())) {
207  return reco::GenParticleRef(g.gens, it - g.genBarcodes->begin());
208  } else {
209  return reco::GenParticleRef();
210  }
211 }
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
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
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 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().

214  {
215  // Make up a GenParticleCandidate from the GEANT track info.
216  int charge = static_cast<int>(tk.charge());
218  Particle::Point vtx; // = (0,0,0) by default
219  if (!tk.noVertex()) vtx = g.simvtxs[tk.vertIndex()].position();
220  GenParticle gp(charge, p4, vtx, tk.type(), setStatus_, true);
221  if (mother.isNonnull()) gp.addMother(mother);
222  return gp;
223 }
float charge() const
charge
Definition: CoreSimTrack.cc:3
double charge(const std::vector< uint8_t > &Ampls)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
double p4[4]
Definition: TauolaWrapper.h:92
math::XYZPoint Point
point in the space
Definition: Particle.h:29
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:40
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
void PATGenCandsFromSimTracksProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 226 of file PATGenCandsFromSimTracksProducer.cc.

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

227  {
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.getByLabel(src_, simtracks);
252 
253  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
254  std::auto_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.getByLabel(src_, 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.getByLabel(genParticles_, gens);
274  event.getByLabel(genParticles_, genBarcodes);
275  if (gens->size() != genBarcodes->size()) throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
276  barcodesAreSorted = __gnu_cxx::is_sorted(genBarcodes->begin(), genBarcodes->end());
277  }
278 
279 
280  // make the output collection
281  auto_ptr<GenParticleCollection> cands(new 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();
287  isimtrk != simtracks->end(); ++isimtrk) {
288 
289  // Skip PYTHIA tracks.
290  if (isimtrk->genpartIndex() != -1) continue;
291 
292  // Maybe apply the PdgId filter
293  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
294  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end()) continue;
295  }
296 
298 
299  // Maybe apply filter on the particle
300  if (filter_.get() != 0) {
301  if (!(*filter_)(genp)) continue;
302  }
303 
304  if (!motherPdgIds_.empty()) {
305  const SimTrack *motherSimTk = findGeantMother(*isimtrk, globals);
306  if (motherSimTk == 0) continue;
307  if (motherPdgIds_.find(std::abs(motherSimTk->type())) == motherPdgIds_.end()) continue;
308  }
309 
311  Ref<GenParticleCollection> motherRef;
312  const SimTrack * mother = findGeantMother(*isimtrk, globals);
313  if (mother != 0) motherRef = findRef(*mother, globals);
314  if (motherRef.isNonnull()) genp.addMother(motherRef);
315  }
316 
317  cands->push_back(genp);
318  }
319 
320  // Write to the Event, and get back a handle (which can be useful for debugging)
321  edm::OrphanHandle<reco::GenParticleCollection> orphans = event.put(cands);
322 
323 #ifdef DEBUG_PATGenCandsFromSimTracksProducer
324  std::cout << "Produced a list of " << orphans->size() << " genParticles." << std::endl;
325  for (GenParticleCollection::const_iterator it = orphans->begin(), ed = orphans->end(); it != ed; ++it) {
326  std::cout << " ";
327  std::cout << "GenParticle #" << (it - orphans->begin()) << ": pdgId " << it->pdgId()
328  << ", pt = " << it->pt() << ", eta = " << it->eta() << ", phi = " << it->phi()
329  << ", rho = " << it->vertex().Rho() << ", z = " << it->vertex().Z() << std::endl;
330  edm::Ref<GenParticleCollection> mom = it->motherRef();
331  size_t depth = 2;
332  while (mom.isNonnull()) {
333  if (mom.id() == orphans.id()) {
334  // I need to re-make the ref because they are not working until this module returns.
335  mom = edm::Ref<GenParticleCollection>(orphans, mom.key());
336  }
337  for (size_t i = 0; i < depth; ++i) std::cout << " ";
338  std::cout << "GenParticleRef [" << mom.id() << "/" << mom.key() << "]: pdgId " << mom->pdgId() << ", status = " << mom->status()
339  << ", pt = " << mom->pt() << ", eta = " << mom->eta() << ", phi = " << mom->phi()
340  << ", rho = " << mom->vertex().Rho() << ", z = " << mom->vertex().Z() << std::endl;
341  if (mom.id() != orphans.id()) break;
342  if ((mom->motherRef().id() == mom.id()) && (mom->motherRef().key() == mom.key())) {
343  throw cms::Exception("Corrupt Data") << "A particle is it's own mother.\n";
344  }
345  mom = mom->motherRef();
346  depth++;
347  }
348  }
349  std::cout << std::endl;
350 #endif
351 
352 }
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::Ref< reco::GenParticleCollection > findRef(const SimTrack &tk, GlobalContext &g) const
#define abs(x)
Definition: mlp_lapack.h:159
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
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:40
edm::InputTag genParticles_
Collection of GenParticles I need to make refs to. It must also have its associated vector&lt;int&gt; of ba...
tuple cout
Definition: gather_cfg.py:121
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
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::InputTag pat::PATGenCandsFromSimTracksProducer::genParticles_
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::InputTag pat::PATGenCandsFromSimTracksProducer::src_
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().