CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GenPlusSimParticleProducer.cc
Go to the documentation of this file.
1 
28 
33 
37 
38 #include <ext/algorithm>
39 #include <memory>
40 
41 namespace pat {
43  public:
46 
47  private:
48  void produce(edm::Event &, const edm::EventSetup &) override;
49  void endJob() override {}
50 
55  std::set<int> pdgIds_; // these are the ones we really use
56  std::vector<PdtEntry> pdts_; // these are needed before we get the EventSetup
57 
59  std::unique_ptr<StrFilter> filter_;
60 
64 
67 
80  void addGenParticle(const SimTrack &stMom,
81  const SimTrack &stDau,
82  unsigned int momidx,
83  const edm::SimTrackContainer &simtks,
84  const edm::SimVertexContainer &simvtxs,
85  reco::GenParticleCollection &mergedGens,
86  const reco::GenParticleRefProd ref,
87  std::vector<int> &genBarcodes,
88  bool &barcodesAreSorted) const;
89  struct LessById {
90  bool operator()(const SimTrack &tk1, const SimTrack &tk2) const { return tk1.trackId() < tk2.trackId(); }
91  bool operator()(const SimTrack &tk1, unsigned int id) const { return tk1.trackId() < id; }
92  bool operator()(unsigned int id, const SimTrack &tk2) const { return id < tk2.trackId(); }
93  };
94  };
95 } // namespace pat
96 
97 using namespace std;
98 using namespace edm;
99 using namespace reco;
101 
102 GenPlusSimParticleProducer::GenPlusSimParticleProducer(const ParameterSet &cfg)
103  : firstEvent_(true),
104  simtracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
105  simverticesToken_(
106  consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
107  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
108  gensToken_(consumes<GenParticleCollection>(
109  cfg.getParameter<InputTag>("genParticles"))), // get the genParticles to add GEANT particles to
110  genBarcodesToken_(consumes<std::vector<int>>(
111  cfg.getParameter<InputTag>("genParticles"))) // get the genParticles to add GEANT particles to
112 {
113  // Possibly allow a list of particle types
114  if (cfg.exists("particleTypes")) {
115  pdts_ = cfg.getParameter<vector<PdtEntry>>("particleTypes");
116  if (!pdts_.empty()) {
118  }
119  }
120 
121  // Possibly allow a string cut
122  if (cfg.existsAs<string>("filter")) {
123  string filter = cfg.getParameter<string>("filter");
124  if (!filter.empty()) {
125  filter_ = std::make_unique<StrFilter>(filter);
126  }
127  }
128  produces<GenParticleCollection>();
129  produces<vector<int>>();
130 }
131 
133  const SimTrack &stDau,
134  unsigned int momidx,
135  const SimTrackContainer &simtracksSorted,
136  const SimVertexContainer &simvertices,
137  reco::GenParticleCollection &mergedGens,
138  const GenParticleRefProd ref,
139  std::vector<int> &genBarcodes,
140  bool &barcodesAreSorted) const {
141  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
142  // Make up a GenParticleCandidate from the GEANT track info.
143  int charge = static_cast<int>(stDau.charge());
144  const Particle::LorentzVector &p4 = stDau.momentum();
145  Particle::Point vtx; // = (0,0,0) by default
146  if (!stDau.noVertex())
147  vtx = simvertices[stDau.vertIndex()].position();
148  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
149 
150  // Maybe apply filter on the particle
151  if (filter_.get() != nullptr) {
152  if (!(*filter_)(genp))
153  return;
154  }
155 
156  reco::GenParticleRef parentRef(ref, momidx);
157  genp.addMother(parentRef);
158  mergedGens.push_back(genp);
159  // get the index for the daughter just added
160  unsigned int dauidx = mergedGens.size() - 1;
161 
162  // update add daughter relationship
163  reco::GenParticle &cand = mergedGens[momidx];
164  cand.addDaughter(GenParticleRef(ref, dauidx));
165 
166  //look for simtrack daughters of stDau to see if we need to recur further down the chain
167 
168  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin(); isimtrk != simtracksSorted.end();
169  ++isimtrk) {
170  if (!isimtrk->noVertex()) {
171  // Pick the vertex (isimtrk.vertIndex() is really an index)
172  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
173 
174  // Check if the vertex has a parent track (otherwise, we're lost)
175  if (!vtx.noParent()) {
176  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
177  unsigned int idx = vtx.parentIndex();
178  SimTrackContainer::const_iterator it =
179  std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
180  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
181  if (it->trackId() == stDau.trackId()) {
182  //need the genparticle index of stDau which is dauidx
184  stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
185  }
186  }
187  }
188  }
189  }
190 }
191 
193  if (firstEvent_) {
194  if (!pdts_.empty()) {
195  auto const &pdt = iSetup.getData(tableToken_);
196  pdgIds_.clear();
197  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
198  itp->setup(pdt); // decode string->pdgId and vice-versa
199  pdgIds_.insert(std::abs(itp->pdgId()));
200  }
201  pdts_.clear();
202  }
203  firstEvent_ = false;
204  }
205 
206  // Simulated tracks (i.e. GEANT particles).
207  Handle<SimTrackContainer> simtracks;
208  event.getByToken(simtracksToken_, simtracks);
209 
210  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
211  std::unique_ptr<SimTrackContainer> simtracksTmp;
212  const SimTrackContainer *simtracksSorted = &*simtracks;
213  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
214  simtracksTmp = std::make_unique<SimTrackContainer>(*simtracks);
215  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
216  simtracksSorted = &*simtracksTmp;
217  }
218 
219  // Get the associated vertices
220  Handle<SimVertexContainer> simvertices;
221  event.getByToken(simverticesToken_, simvertices);
222 
223  // Get the GenParticles and barcodes, if needed to set mother and daughter links
225  Handle<std::vector<int>> genBarcodes;
226  bool barcodesAreSorted = true;
227  event.getByToken(gensToken_, gens);
228  event.getByToken(genBarcodesToken_, genBarcodes);
229  if (gens->size() != genBarcodes->size())
230  throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
231 
232  // make the output collection
233  auto candsPtr = std::make_unique<GenParticleCollection>();
234  GenParticleCollection &cands = *candsPtr;
235 
236  const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
237 
238  // add the original genParticles to the merged output list
239  for (size_t i = 0; i < gens->size(); ++i) {
240  reco::GenParticle cand((*gens)[i]);
241  cands.push_back(cand);
242  }
243 
244  // make new barcodes vector and fill it with the original list
245  auto newGenBarcodes = std::make_unique<vector<int>>();
246  for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
247  newGenBarcodes->push_back((*genBarcodes)[i]);
248  }
249  barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
250 
251  for (size_t i = 0; i < cands.size(); ++i) {
252  reco::GenParticle &cand = cands[i];
253  size_t nDaus = cand.numberOfDaughters();
255  cand.resetDaughters(ref.id());
256  for (size_t d = 0; d < nDaus; ++d) {
257  cand.addDaughter(GenParticleRef(ref, daus[d].key()));
258  }
259 
260  size_t nMoms = cand.numberOfMothers();
262  cand.resetMothers(ref.id());
263  for (size_t m = 0; m < nMoms; ++m) {
264  cand.addMother(GenParticleRef(ref, moms[m].key()));
265  }
266  }
267 
268  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin(); isimtrk != simtracks->end(); ++isimtrk) {
269  // Skip PYTHIA tracks.
270  if (isimtrk->genpartIndex() != -1)
271  continue;
272 
273  // Maybe apply the PdgId filter
274  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
275  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end())
276  continue;
277  }
278 
279  // find simtrack that has a genParticle match to its parent
280  // Look at the production vertex. If there is no vertex, I can do nothing...
281  if (!isimtrk->noVertex()) {
282  // Pick the vertex (isimtrk.vertIndex() is really an index)
283  const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
284 
285  // Check if the vertex has a parent track (otherwise, we're lost)
286  if (!vtx.noParent()) {
287  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
288  unsigned int idx = vtx.parentIndex();
289  SimTrackContainer::const_iterator it =
290  std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
291  if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
292  if (it->genpartIndex() != -1) {
293  std::vector<int>::const_iterator itIndex;
294  if (barcodesAreSorted) {
295  itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
296  } else {
297  itIndex = std::find(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
298  }
299 
300  // Check that I found something
301  // I need to check '*itIndex == it->genpartIndex()' because lower_bound just finds the right spot for an item in a sorted list, not the item
302  if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
303  // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
304  // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
305  unsigned int momidx = itIndex - genBarcodes->begin();
306  addGenParticle(*it,
307  *isimtrk,
308  momidx,
309  *simtracksSorted,
310  *simvertices,
311  *candsPtr,
312  ref,
313  *newGenBarcodes,
314  barcodesAreSorted);
315  }
316  }
317  }
318  }
319  }
320  }
321 
322  event.put(std::move(candsPtr));
323  event.put(std::move(newGenBarcodes));
324 }
325 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool operator()(const SimTrack &tk1, const SimTrack &tk2) const
tuple cfg
Definition: looper.py:296
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
uint16_t *__restrict__ id
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
size_t numberOfDaughters() const override
number of daughters
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
bool exists(std::string const &parameterName) const
checks if a parameter exists
const daughters & daughterRefVector() const
references to daughtes
std::unique_ptr< StrFilter > filter_
size_t numberOfMothers() const override
number of mothers
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
float charge() const
charge
Definition: CoreSimTrack.cc:17
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
bool getData(T &iHolder) const
Definition: EventSetup.h:128
tuple d
Definition: ztail.py:151
void resetDaughters(const edm::ProductID &id)
set daughters product ID
StringCutObjectSelector< reco::GenParticle > StrFilter
int parentIndex() const
Definition: SimVertex.h:29
const mothers & motherRefVector() const
references to mothers
def move
Definition: eostools.py:511
tuple key
prepare the HTCondor submission files and eventually submit them
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool noVertex() const
Definition: SimTrack.h:34
bool operator()(unsigned int id, const SimTrack &tk2) const
GenPlusSimParticleProducer(const edm::ParameterSet &)
Produces reco::GenParticle from SimTracks.
tuple genp
produce generated paricles in acceptance #
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:33
bool operator()(const SimTrack &tk1, unsigned int id) const
unsigned int trackId() const
Definition: CoreSimTrack.h:31
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
edm::EDGetTokenT< reco::GenParticleCollection > gensToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector&lt;int&gt; of ba...
std::vector< SimVertex > SimVertexContainer
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
void produce(edm::Event &, const edm::EventSetup &) override
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
__host__ __device__ constexpr RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
void addGenParticle(const SimTrack &stMom, const SimTrack &stDau, unsigned int momidx, const edm::SimTrackContainer &simtks, const edm::SimVertexContainer &simvtxs, reco::GenParticleCollection &mergedGens, const reco::GenParticleRefProd ref, std::vector< int > &genBarcodes, bool &barcodesAreSorted) const
Try to link the GEANT particle to the generator particle it came from.
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
bool noParent() const
Definition: SimVertex.h:30
void resetMothers(const edm::ProductID &id)
set mother product ID
std::vector< SimTrack > SimTrackContainer
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_