CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATGenCandsFromSimTracksProducer.cc
Go to the documentation of this file.
1 
17 
22 
25 
26 #include <ext/algorithm>
27 
28 namespace pat {
30 public:
33 
34 private:
35  virtual void produce(edm::Event&, const edm::EventSetup&) override;
36  virtual void endJob() override {}
37 
42  std::set<int> pdgIds_; // these are the ones we really use
43  std::vector<PdtEntry> pdts_; // these are needed before we get the EventSetup
44  std::set<int> motherPdgIds_; // these are the ones we really use
45  std::vector<PdtEntry> motherPdts_; // these are needed before we get the EventSetup
46 
48  std::auto_ptr<StrFilter> filter_;
49 
54 
58 
60  struct GlobalContext {
62  const edm::SimVertexContainer &simvtxs1,
64  const edm::Handle<std::vector<int> > &genBarcodes1,
65  bool barcodesAreSorted1,
68  simtks(simtks1), simvtxs(simvtxs1),
69  gens(gens1), genBarcodes(genBarcodes1), barcodesAreSorted(barcodesAreSorted1),
70  output(output1), refprod(refprod1), simTksProcessed() {}
71  // GEANT info
74  // PYTHIA info
78  // MY OUTPUT info
81  // BOOK-KEEPING
82  std::map<unsigned int,int> simTksProcessed; // key = sim track id;
83  // val = 0: not processed;
84  // i>0: (index+1) in my output
85  // i<0: -(index+1) in pythia [NOT USED]
86  };
87 
89  const SimTrack * findGeantMother(const SimTrack &tk, const GlobalContext &g) const ;
96 
101 
102 
103 
104  struct LessById {
105  bool operator()(const SimTrack &tk1, const SimTrack &tk2) const { return tk1.trackId() < tk2.trackId(); }
106  bool operator()(const SimTrack &tk1, unsigned int id ) const { return tk1.trackId() < id; }
107  bool operator()(unsigned int id, const SimTrack &tk2) const { return id < tk2.trackId(); }
108  };
109 
110 };
111 }
112 
113 using namespace std;
114 using namespace edm;
115 using namespace reco;
117 
118 PATGenCandsFromSimTracksProducer::PATGenCandsFromSimTracksProducer(const ParameterSet& cfg) :
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 }
151 
152 const SimTrack *
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 }
167 
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();
190  }
191  } else {
192  // Otherwise, I just return a ref to my mum
193  return motherRef;
194  }
195 }
196 
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 }
216 
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 }
228 
229 
231  const EventSetup& iSetup) {
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 }
357 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
const edm::Handle< reco::GenParticleCollection > & gens
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...
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
bool exists(std::string const &parameterName) const
checks if a parameter exists
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:266
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
float charge() const
charge
Definition: CoreSimTrack.cc:18
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
GlobalContext(const edm::SimTrackContainer &simtks1, const edm::SimVertexContainer &simvtxs1, const edm::Handle< reco::GenParticleCollection > &gens1, const edm::Handle< std::vector< int > > &genBarcodes1, bool barcodesAreSorted1, reco::GenParticleCollection &output1, const edm::RefProd< reco::GenParticleCollection > &refprod1)
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
Produces reco::GenParticle from SimTracks.
ProductID id() const
int parentIndex() const
Definition: SimVertex.h:33
double p4[4]
Definition: TauolaWrapper.h:92
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool noVertex() const
Definition: SimTrack.h:30
bool makeMotherLink_
If true, I&#39;ll try to make a link from the GEANT particle to a GenParticle.
const edm::RefProd< reco::GenParticleCollection > & refprod
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
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 #
StringCutObjectSelector< reco::GenParticle > StrFilter
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
unsigned int trackId() const
Definition: CoreSimTrack.h:34
bool operator()(const SimTrack &tk1, unsigned int id) const
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< SimVertex > SimVertexContainer
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
edm::Ref< reco::GenParticleCollection > generatorRef_(const SimTrack &tk, const GlobalContext &g) const
Used by findRef if the track is a PYTHIA particle.
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
bool operator()(unsigned int id, const SimTrack &tk2) const
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
tuple cout
Definition: gather_cfg.py:121
bool operator()(const SimTrack &tk1, const SimTrack &tk2) const
bool noParent() const
Definition: SimVertex.h:34
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken_
volatile std::atomic< bool > shutdown_flag false
edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
virtual void produce(edm::Event &, const edm::EventSetup &) override
std::vector< SimTrack > SimTrackContainer