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::GenPlusSimParticleProducer Class Reference

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

Inheritance diagram for pat::GenPlusSimParticleProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Classes

struct  LessById
 

Public Member Functions

 GenPlusSimParticleProducer (const edm::ParameterSet &)
 
 ~GenPlusSimParticleProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Types

typedef
StringCutObjectSelector
< reco::GenParticle
StrFilter
 

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

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...
 
std::set< int > pdgIds_
 
std::vector< PdtEntrypdts_
 
int setStatus_
 
edm::InputTag src_
 

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
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

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 "famosSimHits" for FAMOS 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 39 of file GenPlusSimParticleProducer.cc.

Member Typedef Documentation

Definition at line 54 of file GenPlusSimParticleProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 97 of file GenPlusSimParticleProducer.cc.

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

97  :
98  firstEvent_(true),
99  src_(cfg.getParameter<InputTag>("src")), // source sim tracks & vertices
100  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
101  genParticles_(cfg.getParameter<InputTag>("genParticles")) // get the genParticles to add GEANT particles to
102 {
103  // Possibly allow a list of particle types
104  if (cfg.exists("particleTypes")) {
105  pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
106  }
107 
108  // Possibly allow a string cut
109  if (cfg.existsAs<string>("filter")) {
110  string filter = cfg.getParameter<string>("filter");
111  if (!filter.empty()) {
112  filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
113  }
114  }
115  produces<GenParticleCollection>();
116  produces<vector<int> >();
117 }
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
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::GenPlusSimParticleProducer::~GenPlusSimParticleProducer ( )
inline

Definition at line 42 of file GenPlusSimParticleProducer.cc.

42 {}

Member Function Documentation

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 119 of file GenPlusSimParticleProducer.cc.

References reco::CompositeRefCandidateT< D >::addDaughter(), DeDxDiscriminatorTools::charge(), CoreSimTrack::charge(), filter_, EgammaValidation_cff::genp, CoreSimTrack::momentum(), SimVertex::noParent(), SimTrack::noVertex(), p4, SimVertex::parentIndex(), setStatus_, CoreSimTrack::trackId(), CoreSimTrack::type(), and SimTrack::vertIndex().

Referenced by produce().

128 {
129  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
130  // Make up a GenParticleCandidate from the GEANT track info.
131  int charge = static_cast<int>(stDau.charge());
133  Particle::Point vtx; // = (0,0,0) by default
134  if (!stDau.noVertex()) vtx = simvertices[stDau.vertIndex()].position();
135  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
136 
137  // Maybe apply filter on the particle
138  if (filter_.get() != 0) {
139  if (!(*filter_)(genp)) return;
140  }
141 
142  reco::GenParticleRef parentRef(ref, momidx);
143  genp.addMother(parentRef);
144  mergedGens.push_back(genp);
145  // get the index for the daughter just added
146  unsigned int dauidx = mergedGens.size()-1;
147 
148  // update add daughter relationship
149  reco::GenParticle & cand = mergedGens[ momidx ];
150  cand.addDaughter( GenParticleRef( ref, dauidx) );
151 
152  //look for simtrack daughters of stDau to see if we need to recur further down the chain
153 
154  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin();
155  isimtrk != simtracksSorted.end(); ++isimtrk) {
156  if (!isimtrk->noVertex()) {
157  // Pick the vertex (isimtrk.vertIndex() is really an index)
158  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
159 
160  // Check if the vertex has a parent track (otherwise, we're lost)
161  if (!vtx.noParent()) {
162 
163  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
164  unsigned int idx = vtx.parentIndex();
165  SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
166  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
167  if (it->trackId()==stDau.trackId()) {
168  //need the genparticle index of stDau which is dauidx
169  addGenParticle(stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
170  }
171  }
172  }
173  }
174  }
175 }
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
float charge() const
charge
Definition: CoreSimTrack.cc:3
double charge(const std::vector< uint8_t > &Ampls)
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
int parentIndex() const
Definition: SimVertex.h:33
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
tuple genp
produce generated paricles in acceptance #
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:49
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
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.
bool noParent() const
Definition: SimVertex.h:34
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
virtual void pat::GenPlusSimParticleProducer::endJob ( void  )
inlineprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 46 of file GenPlusSimParticleProducer.cc.

46 {}
void GenPlusSimParticleProducer::produce ( edm::Event event,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 177 of file GenPlusSimParticleProducer.cc.

References abs, reco::CompositeRefCandidateT< D >::addDaughter(), addGenParticle(), reco::CompositeRefCandidateT< D >::addMother(), reco::CompositeRefCandidateT< D >::daughterRefVector(), spr::find(), firstEvent_, genParticles_, i, edm::RefProd< T >::id(), combine::key, m, reco::CompositeRefCandidateT< D >::motherRefVector(), SimVertex::noParent(), reco::CompositeRefCandidateT< D >::numberOfDaughters(), reco::CompositeRefCandidateT< D >::numberOfMothers(), SimVertex::parentIndex(), pdgIds_, pdts_, reco::CompositeRefCandidateT< D >::resetDaughters(), reco::CompositeRefCandidateT< D >::resetMothers(), python.multivaluedict::sort(), and src_.

178  {
179  if (firstEvent_){
180  if (!pdts_.empty()) {
181  pdgIds_.clear();
182  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
183  itp->setup(iSetup); // decode string->pdgId and vice-versa
184  pdgIds_.insert(std::abs(itp->pdgId()));
185  }
186  pdts_.clear();
187  }
188  firstEvent_ = false;
189  }
190 
191  // Simulated tracks (i.e. GEANT particles).
192  Handle<SimTrackContainer> simtracks;
193  event.getByLabel(src_, simtracks);
194 
195  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
196  std::auto_ptr<SimTrackContainer> simtracksTmp;
197  const SimTrackContainer * simtracksSorted = &* simtracks;
198  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
199  simtracksTmp.reset(new SimTrackContainer(*simtracks));
200  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
201  simtracksSorted = &* simtracksTmp;
202  }
203 
204  // Get the associated vertices
205  Handle<SimVertexContainer> simvertices;
206  event.getByLabel(src_, simvertices);
207 
208  // Get the GenParticles and barcodes, if needed to set mother and daughter links
210  Handle<std::vector<int> > genBarcodes;
211  bool barcodesAreSorted = true;
212  event.getByLabel(genParticles_, gens);
213  event.getByLabel(genParticles_, genBarcodes);
214  if (gens->size() != genBarcodes->size()) throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
215 
216  // make the output collection
217  auto_ptr<GenParticleCollection> candsPtr(new GenParticleCollection);
218  GenParticleCollection & cands = * candsPtr;
219 
220  const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
221 
222  // add the original genParticles to the merged output list
223  for( size_t i = 0; i < gens->size(); ++ i ) {
224  reco::GenParticle cand((*gens)[i]);
225  cands.push_back(cand);
226  }
227 
228  // make new barcodes vector and fill it with the original list
229  auto_ptr<vector<int> > newGenBarcodes( new vector<int>() );
230  for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
231  newGenBarcodes->push_back((*genBarcodes)[i]);
232  }
233  barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
234 
235  for( size_t i = 0; i < cands.size(); ++ i ) {
236  reco::GenParticle & cand = cands[ i ];
237  size_t nDaus = cand.numberOfDaughters();
239  cand.resetDaughters( ref.id() );
240  for ( size_t d = 0; d < nDaus; ++d) {
241  cand.addDaughter( GenParticleRef( ref, daus[d].key() ) );
242  }
243 
244  size_t nMoms = cand.numberOfMothers();
246  cand.resetMothers( ref.id() );
247  for ( size_t m = 0; m < nMoms; ++m) {
248  cand.addMother( GenParticleRef( ref, moms[m].key() ) );
249  }
250  }
251 
252  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin();
253  isimtrk != simtracks->end(); ++isimtrk) {
254 
255  // Skip PYTHIA tracks.
256  if (isimtrk->genpartIndex() != -1) continue;
257 
258  // Maybe apply the PdgId filter
259  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
260  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end()) continue;
261  }
262 
263  // find simtrack that has a genParticle match to its parent
264  // Look at the production vertex. If there is no vertex, I can do nothing...
265  if (!isimtrk->noVertex()) {
266 
267  // Pick the vertex (isimtrk.vertIndex() is really an index)
268  const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
269 
270  // Check if the vertex has a parent track (otherwise, we're lost)
271  if (!vtx.noParent()) {
272 
273  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
274  unsigned int idx = vtx.parentIndex();
275  SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
276  if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
277  if (it->genpartIndex() != -1) {
278  std::vector<int>::const_iterator itIndex;
279  if (barcodesAreSorted) {
280  itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
281  } else {
282  itIndex = std::find( genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
283  }
284 
285  // Check that I found something
286  // 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
287  if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
288  // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
289  // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
290  unsigned int momidx = itIndex - genBarcodes->begin();
291  addGenParticle(*it, *isimtrk, momidx, *simtracksSorted, *simvertices, *candsPtr, ref, *newGenBarcodes, barcodesAreSorted);
292  }
293  }
294  }
295  }
296  }
297  }
298 
299  event.put(candsPtr);
300  event.put(newGenBarcodes);
301 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int i
Definition: DBlmapReader.cc:9
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
#define abs(x)
Definition: mlp_lapack.h:159
const daughters & daughterRefVector() const
references to daughtes
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
void resetDaughters(const edm::ProductID &id)
set daughters product ID
int parentIndex() const
Definition: SimVertex.h:33
const mothers & motherRefVector() const
references to mothers
virtual size_t numberOfMothers() const
number of mothers
virtual size_t numberOfDaughters() const
number of daughters
void addMother(const typename mothers::value_type &)
add a daughter via a reference
edm::InputTag genParticles_
Collection of GenParticles I need to make refs to. It must also have its associated vector&lt;int&gt; of ba...
list key
Definition: combine.py:13
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:140
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.
bool noParent() const
Definition: SimVertex.h:34
void resetMothers(const edm::ProductID &id)
set mother product ID
std::vector< SimTrack > SimTrackContainer

Member Data Documentation

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

Definition at line 55 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle(), and GenPlusSimParticleProducer().

bool pat::GenPlusSimParticleProducer::firstEvent_
private

Definition at line 48 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

edm::InputTag pat::GenPlusSimParticleProducer::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 58 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 51 of file GenPlusSimParticleProducer.cc.

Referenced by produce().

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

Definition at line 52 of file GenPlusSimParticleProducer.cc.

Referenced by GenPlusSimParticleProducer(), and produce().

int pat::GenPlusSimParticleProducer::setStatus_
private

Definition at line 50 of file GenPlusSimParticleProducer.cc.

Referenced by addGenParticle().

edm::InputTag pat::GenPlusSimParticleProducer::src_
private

Definition at line 49 of file GenPlusSimParticleProducer.cc.

Referenced by produce().