CMS 3D CMS Logo

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

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

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

Classes

struct  LessById
 

Public Member Functions

 GenPlusSimParticleProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Types

typedef StringCutObjectSelector< reco::GenParticleStrFilter
 

Private Member Functions

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. More...
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::unique_ptr< StrFilterfilter_
 
bool firstEvent_
 
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
 
edm::EDGetTokenT< reco::GenParticleCollectiongensToken_
 Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of barcodes, aligned with them. More...
 
std::set< int > pdgIds_
 
std::vector< PdtEntrypdts_
 
int setStatus_
 
edm::EDGetTokenT< edm::SimTrackContainersimtracksToken_
 
edm::EDGetTokenT< edm::SimVertexContainersimverticesToken_
 
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordtableToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Produces reco::GenParticle from SimTracks.

The GenPlusSimParticleProducer produces GenParticles from SimTracks, so they can be used for MC matching. A copy of the original genParticle list is made and the genParticles created from the GEANT/FAMOS particles are added to the list including all ancestors and correct mother/daughter references

Sample useage in cfg.py file: process.genParticlePlusGEANT = cms.EDProducer("GenPlusSimParticleProducer", src = cms.InputTag("g4SimHits"), # use "fastSimProducer" for FastSim setStatus = cms.int32(8), # set status = 8 for GEANT GPs particleTypes = cms.vstring("pi+"), # also picks pi- (optional) filter = cms.vstring("pt > 0.0"), # just for testing genParticles = cms.InputTag("genParticles") # original genParticle list )

Author
Jordan Tucker (original module), Keith Ulmer (generalization)

Definition at line 42 of file GenPlusSimParticleProducer.cc.

Member Typedef Documentation

◆ StrFilter

Definition at line 56 of file GenPlusSimParticleProducer.cc.

Constructor & Destructor Documentation

◆ GenPlusSimParticleProducer()

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

Definition at line 100 of file GenPlusSimParticleProducer.cc.

References looper::cfg, deDxTools::esConsumes(), ALCARECOTkAlBeamHalo_cff::filter, filter_, pdts_, and tableToken_.

101  : firstEvent_(true),
102  simtracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
104  consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
105  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
106  gensToken_(consumes<GenParticleCollection>(
107  cfg.getParameter<InputTag>("genParticles"))), // get the genParticles to add GEANT particles to
108  genBarcodesToken_(consumes<std::vector<int>>(
109  cfg.getParameter<InputTag>("genParticles"))) // get the genParticles to add GEANT particles to
110 {
111  // Possibly allow a list of particle types
112  if (cfg.exists("particleTypes")) {
113  pdts_ = cfg.getParameter<vector<PdtEntry>>("particleTypes");
114  if (!pdts_.empty()) {
116  }
117  }
118 
119  // Possibly allow a string cut
120  if (cfg.existsAs<string>("filter")) {
121  string filter = cfg.getParameter<string>("filter");
122  if (!filter.empty()) {
123  filter_ = std::make_unique<StrFilter>(filter);
124  }
125  }
126  produces<GenParticleCollection>();
127  produces<vector<int>>();
128 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
std::unique_ptr< StrFilter > filter_
edm::EDGetTokenT< reco::GenParticleCollection > gensToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_

Member Function Documentation

◆ addGenParticle()

void GenPlusSimParticleProducer::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
private

Try to link the GEANT particle to the generator particle it came from.

Arguments: – Specific – gp: GenParticle made from the GEANT particle st: The GEANT simTrack for which we create a genParticle

– GEANT related – simtks: A list of GEANT tracks, sorted by track id simvtxs: The list of GEANT vertices, in their natural order (skimtks have pointers into this vector!)

GenParticle related – gens : Handle to the GenParticles, to make the ref to genBarcodes : Barcodes for each GenParticle, in a vector aligned to the GenParticleCollection. barcodesAreSorted: true if the barcodes are sorted (which means I can use binary_search on them)

Definition at line 130 of file GenPlusSimParticleProducer.cc.

References CoreSimTrack::charge(), ALCARECOTkAlJpsiMuMu_cff::charge, filter_, EgammaValidation_cff::genp, heavyIonCSV_trainingSettings::idx, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, CoreSimTrack::momentum(), SimTrack::noVertex(), setStatus_, CoreSimTrack::trackId(), CoreSimTrack::type(), SimTrack::vertIndex(), and L1BJetProducer_cff::vtx.

Referenced by produce().

138  {
139  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
140  // Make up a GenParticleCandidate from the GEANT track info.
141  int charge = static_cast<int>(stDau.charge());
142  const Particle::LorentzVector &p4 = stDau.momentum();
143  Particle::Point vtx; // = (0,0,0) by default
144  if (!stDau.noVertex())
145  vtx = simvertices[stDau.vertIndex()].position();
146  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
147 
148  // Maybe apply filter on the particle
149  if (filter_.get() != nullptr) {
150  if (!(*filter_)(genp))
151  return;
152  }
153 
154  reco::GenParticleRef parentRef(ref, momidx);
155  genp.addMother(parentRef);
156  mergedGens.push_back(genp);
157  // get the index for the daughter just added
158  unsigned int dauidx = mergedGens.size() - 1;
159 
160  // update add daughter relationship
161  reco::GenParticle &cand = mergedGens[momidx];
162  cand.addDaughter(GenParticleRef(ref, dauidx));
163 
164  //look for simtrack daughters of stDau to see if we need to recur further down the chain
165 
166  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin(); isimtrk != simtracksSorted.end();
167  ++isimtrk) {
168  if (!isimtrk->noVertex()) {
169  // Pick the vertex (isimtrk.vertIndex() is really an index)
170  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
171 
172  // Check if the vertex has a parent track (otherwise, we're lost)
173  if (!vtx.noParent()) {
174  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
175  unsigned int idx = vtx.parentIndex();
176  SimTrackContainer::const_iterator it =
177  std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
178  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
179  if (it->trackId() == stDau.trackId()) {
180  //need the genparticle index of stDau which is dauidx
182  stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
183  }
184  }
185  }
186  }
187  }
188 }
genp
produce generated paricles in acceptance #
std::unique_ptr< StrFilter > filter_
float charge() const
charge
Definition: CoreSimTrack.cc:17
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:33
bool noVertex() const
Definition: SimTrack.h:34
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
unsigned int trackId() const
Definition: CoreSimTrack.h:31
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::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27

◆ produce()

void GenPlusSimParticleProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 190 of file GenPlusSimParticleProducer.cc.

References funct::abs(), addGenParticle(), HLT_2023v12_cff::cands, ztail::d, spr::find(), firstEvent_, genBarcodesToken_, gensToken_, edm::EventSetup::getData(), mps_fire::i, edm::RefProd< C >::id(), heavyIonCSV_trainingSettings::idx, crabWrapper::key, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, visualization-live-secondInstance_cfg::m, eostools::move(), pdgIds_, pdts_, simtracksToken_, simverticesToken_, jetUpdater_cfi::sort, tableToken_, and L1BJetProducer_cff::vtx.

190  {
191  if (firstEvent_) {
192  if (!pdts_.empty()) {
193  auto const &pdt = iSetup.getData(tableToken_);
194  pdgIds_.clear();
195  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
196  itp->setup(pdt); // decode string->pdgId and vice-versa
197  pdgIds_.insert(std::abs(itp->pdgId()));
198  }
199  pdts_.clear();
200  }
201  firstEvent_ = false;
202  }
203 
204  // Simulated tracks (i.e. GEANT particles).
205  Handle<SimTrackContainer> simtracks;
206  event.getByToken(simtracksToken_, simtracks);
207 
208  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
209  std::unique_ptr<SimTrackContainer> simtracksTmp;
210  const SimTrackContainer *simtracksSorted = &*simtracks;
211  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
212  simtracksTmp = std::make_unique<SimTrackContainer>(*simtracks);
213  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
214  simtracksSorted = &*simtracksTmp;
215  }
216 
217  // Get the associated vertices
218  Handle<SimVertexContainer> simvertices;
219  event.getByToken(simverticesToken_, simvertices);
220 
221  // Get the GenParticles and barcodes, if needed to set mother and daughter links
223  Handle<std::vector<int>> genBarcodes;
224  bool barcodesAreSorted = true;
225  event.getByToken(gensToken_, gens);
226  event.getByToken(genBarcodesToken_, genBarcodes);
227  if (gens->size() != genBarcodes->size())
228  throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
229 
230  // make the output collection
231  auto candsPtr = std::make_unique<GenParticleCollection>();
232  GenParticleCollection &cands = *candsPtr;
233 
234  const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
235 
236  // add the original genParticles to the merged output list
237  for (size_t i = 0; i < gens->size(); ++i) {
238  reco::GenParticle cand((*gens)[i]);
239  cands.push_back(cand);
240  }
241 
242  // make new barcodes vector and fill it with the original list
243  auto newGenBarcodes = std::make_unique<vector<int>>();
244  for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
245  newGenBarcodes->push_back((*genBarcodes)[i]);
246  }
247  barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
248 
249  for (size_t i = 0; i < cands.size(); ++i) {
251  size_t nDaus = cand.numberOfDaughters();
252  GenParticleRefVector daus = cand.daughterRefVector();
253  cand.resetDaughters(ref.id());
254  for (size_t d = 0; d < nDaus; ++d) {
255  cand.addDaughter(GenParticleRef(ref, daus[d].key()));
256  }
257 
258  size_t nMoms = cand.numberOfMothers();
259  GenParticleRefVector moms = cand.motherRefVector();
260  cand.resetMothers(ref.id());
261  for (size_t m = 0; m < nMoms; ++m) {
262  cand.addMother(GenParticleRef(ref, moms[m].key()));
263  }
264  }
265 
266  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin(); isimtrk != simtracks->end(); ++isimtrk) {
267  // Skip PYTHIA tracks.
268  if (isimtrk->genpartIndex() != -1)
269  continue;
270 
271  // Maybe apply the PdgId filter
272  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
273  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end())
274  continue;
275  }
276 
277  // find simtrack that has a genParticle match to its parent
278  // Look at the production vertex. If there is no vertex, I can do nothing...
279  if (!isimtrk->noVertex()) {
280  // Pick the vertex (isimtrk.vertIndex() is really an index)
281  const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
282 
283  // Check if the vertex has a parent track (otherwise, we're lost)
284  if (!vtx.noParent()) {
285  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
286  unsigned int idx = vtx.parentIndex();
287  SimTrackContainer::const_iterator it =
288  std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
289  if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
290  if (it->genpartIndex() != -1) {
291  std::vector<int>::const_iterator itIndex;
292  if (barcodesAreSorted) {
293  itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
294  } else {
295  itIndex = std::find(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
296  }
297 
298  // Check that I found something
299  // 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
300  if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
301  // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
302  // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
303  unsigned int momidx = itIndex - genBarcodes->begin();
304  addGenParticle(*it,
305  *isimtrk,
306  momidx,
307  *simtracksSorted,
308  *simvertices,
309  *candsPtr,
310  ref,
311  *newGenBarcodes,
312  barcodesAreSorted);
313  }
314  }
315  }
316  }
317  }
318  }
319 
320  event.put(std::move(candsPtr));
321  event.put(std::move(newGenBarcodes));
322 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
d
Definition: ztail.py:151
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
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<int> of ba...
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.
std::vector< SimTrack > SimTrackContainer
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_

Member Data Documentation

◆ filter_

std::unique_ptr<StrFilter> pat::GenPlusSimParticleProducer::filter_
private

Definition at line 57 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle(), and GenPlusSimParticleProducer().

◆ firstEvent_

bool pat::GenPlusSimParticleProducer::firstEvent_
private

Definition at line 49 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ genBarcodesToken_

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

Definition at line 61 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ gensToken_

edm::EDGetTokenT<reco::GenParticleCollection> pat::GenPlusSimParticleProducer::gensToken_
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 60 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ pdgIds_

std::set<int> pat::GenPlusSimParticleProducer::pdgIds_
private

Definition at line 53 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ pdts_

std::vector<PdtEntry> pat::GenPlusSimParticleProducer::pdts_
private

Definition at line 54 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().

◆ setStatus_

int pat::GenPlusSimParticleProducer::setStatus_
private

Definition at line 52 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle().

◆ simtracksToken_

edm::EDGetTokenT<edm::SimTrackContainer> pat::GenPlusSimParticleProducer::simtracksToken_
private

Definition at line 50 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ simverticesToken_

edm::EDGetTokenT<edm::SimVertexContainer> pat::GenPlusSimParticleProducer::simverticesToken_
private

Definition at line 51 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

◆ tableToken_

edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> pat::GenPlusSimParticleProducer::tableToken_
private

Definition at line 63 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().