CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenPlusSimParticleProducer.cc
Go to the documentation of this file.
1 
27 
32 
35 
36 #include <ext/algorithm>
37 
38 namespace pat {
40 public:
43 
44 private:
45  virtual void produce(edm::Event&, const edm::EventSetup&) override;
46  virtual void endJob() override {}
47 
52  std::set<int> pdgIds_; // these are the ones we really use
53  std::vector<PdtEntry> pdts_; // these are needed before we get the EventSetup
54 
56  std::auto_ptr<StrFilter> filter_;
57 
61 
63 
76  void addGenParticle(const SimTrack &stMom,
77  const SimTrack &stDau,
78  unsigned int momidx,
79  const edm::SimTrackContainer &simtks,
80  const edm::SimVertexContainer &simvtxs,
81  reco::GenParticleCollection &mergedGens,
82  const reco::GenParticleRefProd ref,
83  std::vector<int> &genBarcodes,
84  bool &barcodesAreSorted ) const ;
85  struct LessById {
86  bool operator()(const SimTrack &tk1, const SimTrack &tk2) const { return tk1.trackId() < tk2.trackId(); }
87  bool operator()(const SimTrack &tk1, unsigned int id ) const { return tk1.trackId() < id; }
88  bool operator()(unsigned int id, const SimTrack &tk2) const { return id < tk2.trackId(); }
89  };
90 
91 };
92 }
93 
94 using namespace std;
95 using namespace edm;
96 using namespace reco;
98 
99 GenPlusSimParticleProducer::GenPlusSimParticleProducer(const ParameterSet& cfg) :
100  firstEvent_(true),
101  simtracksToken_(consumes<SimTrackContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
102  simverticesToken_(consumes<SimVertexContainer>(cfg.getParameter<InputTag>("src"))), // source sim tracks & vertices
103  setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
104  gensToken_(consumes<GenParticleCollection>(cfg.getParameter<InputTag>("genParticles"))), // get the genParticles to add GEANT particles to
105  genBarcodesToken_(consumes<std::vector<int> >(cfg.getParameter<InputTag>("genParticles"))) // get the genParticles to add GEANT particles to
106 {
107  // Possibly allow a list of particle types
108  if (cfg.exists("particleTypes")) {
109  pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
110  }
111 
112  // Possibly allow a string cut
113  if (cfg.existsAs<string>("filter")) {
114  string filter = cfg.getParameter<string>("filter");
115  if (!filter.empty()) {
116  filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
117  }
118  }
119  produces<GenParticleCollection>();
120  produces<vector<int> >();
121 }
122 
124  const SimTrack &stDau,
125  unsigned int momidx,
126  const SimTrackContainer &simtracksSorted,
127  const SimVertexContainer &simvertices,
128  reco::GenParticleCollection &mergedGens,
129  const GenParticleRefProd ref,
130  std::vector<int> &genBarcodes,
131  bool &barcodesAreSorted) const
132 {
133  // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
134  // Make up a GenParticleCandidate from the GEANT track info.
135  int charge = static_cast<int>(stDau.charge());
137  Particle::Point vtx; // = (0,0,0) by default
138  if (!stDau.noVertex()) vtx = simvertices[stDau.vertIndex()].position();
139  GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
140 
141  // Maybe apply filter on the particle
142  if (filter_.get() != 0) {
143  if (!(*filter_)(genp)) return;
144  }
145 
146  reco::GenParticleRef parentRef(ref, momidx);
147  genp.addMother(parentRef);
148  mergedGens.push_back(genp);
149  // get the index for the daughter just added
150  unsigned int dauidx = mergedGens.size()-1;
151 
152  // update add daughter relationship
153  reco::GenParticle & cand = mergedGens[ momidx ];
154  cand.addDaughter( GenParticleRef( ref, dauidx) );
155 
156  //look for simtrack daughters of stDau to see if we need to recur further down the chain
157 
158  for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin();
159  isimtrk != simtracksSorted.end(); ++isimtrk) {
160  if (!isimtrk->noVertex()) {
161  // Pick the vertex (isimtrk.vertIndex() is really an index)
162  const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
163 
164  // Check if the vertex has a parent track (otherwise, we're lost)
165  if (!vtx.noParent()) {
166 
167  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
168  unsigned int idx = vtx.parentIndex();
169  SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
170  if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
171  if (it->trackId()==stDau.trackId()) {
172  //need the genparticle index of stDau which is dauidx
173  addGenParticle(stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
174  }
175  }
176  }
177  }
178  }
179 }
180 
182  const EventSetup& iSetup) {
183  if (firstEvent_){
184  if (!pdts_.empty()) {
185  pdgIds_.clear();
186  for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
187  itp->setup(iSetup); // decode string->pdgId and vice-versa
188  pdgIds_.insert(std::abs(itp->pdgId()));
189  }
190  pdts_.clear();
191  }
192  firstEvent_ = false;
193  }
194 
195  // Simulated tracks (i.e. GEANT particles).
196  Handle<SimTrackContainer> simtracks;
197  event.getByToken(simtracksToken_, simtracks);
198 
199  // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
200  std::auto_ptr<SimTrackContainer> simtracksTmp;
201  const SimTrackContainer * simtracksSorted = &* simtracks;
202  if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
203  simtracksTmp.reset(new SimTrackContainer(*simtracks));
204  std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
205  simtracksSorted = &* simtracksTmp;
206  }
207 
208  // Get the associated vertices
209  Handle<SimVertexContainer> simvertices;
210  event.getByToken(simverticesToken_, simvertices);
211 
212  // Get the GenParticles and barcodes, if needed to set mother and daughter links
214  Handle<std::vector<int> > genBarcodes;
215  bool barcodesAreSorted = true;
216  event.getByToken(gensToken_, gens);
217  event.getByToken(genBarcodesToken_, genBarcodes);
218  if (gens->size() != genBarcodes->size()) throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
219 
220  // make the output collection
221  auto_ptr<GenParticleCollection> candsPtr(new GenParticleCollection);
222  GenParticleCollection & cands = * candsPtr;
223 
224  const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
225 
226  // add the original genParticles to the merged output list
227  for( size_t i = 0; i < gens->size(); ++ i ) {
228  reco::GenParticle cand((*gens)[i]);
229  cands.push_back(cand);
230  }
231 
232  // make new barcodes vector and fill it with the original list
233  auto_ptr<vector<int> > newGenBarcodes( new vector<int>() );
234  for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
235  newGenBarcodes->push_back((*genBarcodes)[i]);
236  }
237  barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
238 
239  for( size_t i = 0; i < cands.size(); ++ i ) {
240  reco::GenParticle & cand = cands[ i ];
241  size_t nDaus = cand.numberOfDaughters();
243  cand.resetDaughters( ref.id() );
244  for ( size_t d = 0; d < nDaus; ++d) {
245  cand.addDaughter( GenParticleRef( ref, daus[d].key() ) );
246  }
247 
248  size_t nMoms = cand.numberOfMothers();
250  cand.resetMothers( ref.id() );
251  for ( size_t m = 0; m < nMoms; ++m) {
252  cand.addMother( GenParticleRef( ref, moms[m].key() ) );
253  }
254  }
255 
256  for (SimTrackContainer::const_iterator isimtrk = simtracks->begin();
257  isimtrk != simtracks->end(); ++isimtrk) {
258 
259  // Skip PYTHIA tracks.
260  if (isimtrk->genpartIndex() != -1) continue;
261 
262  // Maybe apply the PdgId filter
263  if (!pdgIds_.empty()) { // if we have a filter on pdg ids
264  if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end()) continue;
265  }
266 
267  // find simtrack that has a genParticle match to its parent
268  // Look at the production vertex. If there is no vertex, I can do nothing...
269  if (!isimtrk->noVertex()) {
270 
271  // Pick the vertex (isimtrk.vertIndex() is really an index)
272  const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
273 
274  // Check if the vertex has a parent track (otherwise, we're lost)
275  if (!vtx.noParent()) {
276 
277  // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it
278  unsigned int idx = vtx.parentIndex();
279  SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
280  if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
281  if (it->genpartIndex() != -1) {
282  std::vector<int>::const_iterator itIndex;
283  if (barcodesAreSorted) {
284  itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
285  } else {
286  itIndex = std::find( genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
287  }
288 
289  // Check that I found something
290  // 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
291  if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
292  // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
293  // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
294  unsigned int momidx = itIndex - genBarcodes->begin();
295  addGenParticle(*it, *isimtrk, momidx, *simtracksSorted, *simvertices, *candsPtr, ref, *newGenBarcodes, barcodesAreSorted);
296  }
297  }
298  }
299  }
300  }
301  }
302 
303  event.put(candsPtr);
304  event.put(newGenBarcodes);
305 }
306 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool operator()(const SimTrack &tk1, const SimTrack &tk2) const
tuple cfg
Definition: looper.py:293
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
#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
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
float charge() const
charge
Definition: CoreSimTrack.cc:18
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
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:33
const mothers & motherRefVector() const
references to mothers
virtual size_t numberOfMothers() const
number of mothers
double p4[4]
Definition: TauolaWrapper.h:92
virtual size_t numberOfDaughters() const
number of daughters
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool noVertex() const
Definition: SimTrack.h:30
bool operator()(unsigned int id, const SimTrack &tk2) const
GenPlusSimParticleProducer(const edm::ParameterSet &)
Produces reco::GenParticle from SimTracks.
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
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
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:29
bool operator()(const SimTrack &tk1, unsigned int id) const
unsigned int trackId() const
Definition: CoreSimTrack.h:34
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
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
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
virtual void produce(edm::Event &, const edm::EventSetup &) override
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
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:41
bool noParent() const
Definition: SimVertex.h:34
void resetMothers(const edm::ProductID &id)
set mother product ID
std::vector< SimTrack > SimTrackContainer
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_