CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FBaseSimEvent.cc
Go to the documentation of this file.
1 //HepMC Headers
2 #include "HepMC/GenEvent.h"
3 #include "HepMC/GenVertex.h"
4 #include "HepMC/GenParticle.h"
5 
6 //Framework Headers
8 
9 //CMSSW Data Formats
12 
13 //FAMOS Headers
19 
21 
22 using namespace HepPDT;
23 
24 // system include
25 #include <iostream>
26 #include <iomanip>
27 #include <map>
28 #include <string>
29 
31  :
32  nSimTracks(0),
33  nSimVertices(0),
34  nGenParticles(0),
35  nChargedParticleTracks(0),
36  initialSize(5000)
37 {
38 
39  // Initialize the vectors of particles and vertices
40  theGenParticles = new std::vector<HepMC::GenParticle*>();
41  theSimTracks = new std::vector<FSimTrack>;
42  theSimVertices = new std::vector<FSimVertex>;
43  theChargedTracks = new std::vector<unsigned>();
45 
46  // Reserve some size to avoid mutiple copies
47  /* */
48  theSimTracks->resize(initialSize);
49  theSimVertices->resize(initialSize);
52  theFSimVerticesType->resize(initialSize);
57  /* */
58 
59  // Initialize the Particle filter
60  myFilter = new KineParticleFilter(kine);
61 
62 }
63 
65 
66  // Clear the vectors
67  theGenParticles->clear();
68  theSimTracks->clear();
69  theSimVertices->clear();
70  theChargedTracks->clear();
71  theFSimVerticesType->clear();
72 
73  // Delete
74  delete theGenParticles;
75  delete theSimTracks;
76  delete theSimVertices;
77  delete theChargedTracks;
78  delete theFSimVerticesType;
79  delete myFilter;
80 
81 }
82 
83 void
85 
86  pdt = aPdt;
87 
88 }
89 
90 void
91 FBaseSimEvent::fill(const HepMC::GenEvent& myGenEvent) {
92 
93  // Clear old vectors
94  clear();
95 
96  // Add the particles in the FSimEvent
97  addParticles(myGenEvent);
98 
99  /*
100  std::cout << "The MC truth! " << std::endl;
101  printMCTruth(myGenEvent);
102 
103  std::cout << std::endl << "The FAMOS event! " << std::endl;
104  print();
105  */
106 
107 }
108 
109 void
111 
112  // Clear old vectors
113  clear();
114 
115  // Add the particles in the FSimEvent
116  addParticles(myGenParticles);
117 
118 }
119 
120 void
121 FBaseSimEvent::fill(const std::vector<SimTrack>& simTracks,
122  const std::vector<SimVertex>& simVertices) {
123 
124  // Watch out there ! A SimVertex is in mm (stupid),
125  // while a FSimVertex is in cm (clever).
126 
127  clear();
128 
129  unsigned nVtx = simVertices.size();
130  unsigned nTks = simTracks.size();
131 
132  // Empty event, do nothin'
133  if ( nVtx == 0 ) return;
134 
135  // Two arrays for internal use.
136  std::vector<int> myVertices(nVtx,-1);
137  std::vector<int> myTracks(nTks,-1);
138 
139  // create a map associating geant particle id and position in the
140  // event SimTrack vector
141 
142  std::map<unsigned, unsigned> geantToIndex;
143  for( unsigned it=0; it<simTracks.size(); ++it ) {
144  geantToIndex[ simTracks[it].trackId() ] = it;
145  }
146 
147  // Create also a map associating a SimTrack with its endVertex
148  /*
149  std::map<unsigned, unsigned> endVertex;
150  for ( unsigned iv=0; iv<simVertices.size(); ++iv ) {
151  endVertex[ simVertices[iv].parentIndex() ] = iv;
152  }
153  */
154 
155  // Set the main vertex for the kine particle filter
156  // SimVertices were in mm until 110_pre2
157  // HepLorentzVector primaryVertex = simVertices[0].position()/10.;
158  // SImVertices are now in cm
159  // Also : position is copied until SimVertex switches to Mathcore.
160  // XYZTLorentzVector primaryVertex = simVertices[0].position();
161  // The next 5 lines to be then replaced by the previous line
162  XYZTLorentzVector primaryVertex(simVertices[0].position().x(),
163  simVertices[0].position().y(),
164  simVertices[0].position().z(),
165  simVertices[0].position().t());
166  //
167  myFilter->setMainVertex(primaryVertex);
168  // Add the main vertex to the list.
170  myVertices[0] = 0;
171 
172  for( unsigned trackId=0; trackId<nTks; ++trackId ) {
173 
174  // The track
175  const SimTrack& track = simTracks[trackId];
176  // std::cout << std::endl << "SimTrack " << trackId << " " << track << std::endl;
177 
178  // The origin vertex
179  int vertexId = track.vertIndex();
180  const SimVertex& vertex = simVertices[vertexId];
181  //std::cout << "Origin vertex " << vertexId << " " << vertex << std::endl;
182 
183  // The mother track
184  int motherId = -1;
185  if( !vertex.noParent() ) { // there is a parent to this vertex
186  // geant id of the mother
187  unsigned motherGeantId = vertex.parentIndex();
188  std::map<unsigned, unsigned >::iterator association
189  = geantToIndex.find( motherGeantId );
190  if(association != geantToIndex.end() )
191  motherId = association->second;
192  }
193  int originId = motherId == - 1 ? -1 : myTracks[motherId];
194  //std::cout << "Origin id " << originId << std::endl;
195 
196  /*
197  if ( endVertex.find(trackId) != endVertex.end() )
198  std::cout << "End vertex id = " << endVertex[trackId] << std::endl;
199  else
200  std::cout << "No endVertex !!! " << std::endl;
201  std::cout << "Tracker surface position " << track.trackerSurfacePosition() << std::endl;
202  */
203 
204  // Add the vertex (if it does not already exist!)
205  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
206  vertex.position().pz(),vertex.position().e());
207  if ( myVertices[vertexId] == -1 )
208  // Momentum and position are copied until SimTrack and SimVertex
209  // switch to Mathcore.
210  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
211  // The next line to be then replaced by the previous line
212  myVertices[vertexId] = addSimVertex(position,originId);
213 
214  // Add the track (with protection for brem'ing electrons and muons)
215  int motherType = motherId == -1 ? 0 : simTracks[motherId].type();
216 
217  bool notBremInDetector =
218  (abs(motherType) != 11 && abs(motherType) != 13) ||
219  motherType != track.type() ||
220  position.Perp2() < lateVertexPosition ;
221 
222  if ( notBremInDetector ) {
223  // Momentum and position are copied until SimTrack and SimVertex
224  // switch to Mathcore.
225  // RawParticle part(track.momentum(), vertex.position());
226  // The next 3 lines to be then replaced by the previous line
227  XYZTLorentzVector momentum(track.momentum().px(),track.momentum().py(),
228  track.momentum().pz(),track.momentum().e());
229  RawParticle part(momentum,position);
230  //
231  part.setID(track.type());
232  //std::cout << "Ctau = " << part.PDGcTau() << std::endl;
233  // Don't save tracks that have decayed immediately but for which no daughters
234  // were saved (probably due to cuts on E, pT and eta)
235  // if ( part.PDGcTau() > 0.1 || endVertex.find(trackId) != endVertex.end() )
236  myTracks[trackId] = addSimTrack(&part,myVertices[vertexId],track.genpartIndex());
237  if ( myTracks[trackId] >= 0 ) {
238  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
239  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
240  }
241  } else {
242 
243  myTracks[trackId] = myTracks[motherId];
244  if ( myTracks[trackId] >= 0 ) {
245  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
246  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
247  }
248  }
249 
250  }
251 
252  // Now loop over the remaining end vertices !
253  for( unsigned vertexId=0; vertexId<nVtx; ++vertexId ) {
254 
255  // if the vertex is already saved, just ignore.
256  if ( myVertices[vertexId] != -1 ) continue;
257 
258  // The yet unused vertex
259  const SimVertex& vertex = simVertices[vertexId];
260 
261  // The mother track
262  int motherId = -1;
263  if( !vertex.noParent() ) { // there is a parent to this vertex
264 
265  // geant id of the mother
266  unsigned motherGeantId = vertex.parentIndex();
267  std::map<unsigned, unsigned >::iterator association
268  = geantToIndex.find( motherGeantId );
269  if(association != geantToIndex.end() )
270  motherId = association->second;
271  }
272  int originId = motherId == - 1 ? -1 : myTracks[motherId];
273 
274  // Add the vertex
275  // Momentum and position are copied until SimTrack and SimVertex
276  // switch to Mathcore.
277  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
278  // The next 3 lines to be then replaced by the previous line
279  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
280  vertex.position().pz(),vertex.position().e());
281  myVertices[vertexId] = addSimVertex(position,originId);
282  }
283 
284  // Finally, propagate all particles to the calorimeters
285  BaseParticlePropagator myPart;
286  XYZTLorentzVector mom;
287  XYZTLorentzVector pos;
288 
289 
290  // Loop over the tracks
291  for( int fsimi=0; fsimi < (int)nTracks() ; ++fsimi) {
292 
293 
294  FSimTrack& myTrack = track(fsimi);
295  double trackerSurfaceTime = myTrack.vertex().position().t()
296  + myTrack.momentum().e()/myTrack.momentum().pz()
297  * ( myTrack.trackerSurfacePosition().z()
298  - myTrack.vertex().position().z() );
299  pos = XYZTLorentzVector(myTrack.trackerSurfacePosition().x(),
300  myTrack.trackerSurfacePosition().y(),
301  myTrack.trackerSurfacePosition().z(),
302  trackerSurfaceTime);
303  mom = XYZTLorentzVector(myTrack.trackerSurfaceMomentum().x(),
304  myTrack.trackerSurfaceMomentum().y(),
305  myTrack.trackerSurfaceMomentum().z(),
306  myTrack.trackerSurfaceMomentum().t());
307 
308  if ( mom.T() > 0. ) {
309  // The particle to be propagated
310  myPart = BaseParticlePropagator(RawParticle(mom,pos),0.,0.,4.);
311  myPart.setCharge(myTrack.charge());
312 
313  // Propagate to Preshower layer 1
314  myPart.propagateToPreshowerLayer1(false);
315  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
316  myTrack.setLayer1(myPart,myPart.getSuccess());
317 
318  // Propagate to Preshower Layer 2
319  myPart.propagateToPreshowerLayer2(false);
320  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
321  myTrack.setLayer2(myPart,myPart.getSuccess());
322 
323  // Propagate to Ecal Endcap
324  myPart.propagateToEcalEntrance(false);
325  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
326  myTrack.setEcal(myPart,myPart.getSuccess());
327 
328  // Propagate to HCAL entrance
329  myPart.propagateToHcalEntrance(false);
330  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
331  myTrack.setHcal(myPart,myPart.getSuccess());
332 
333  // Attempt propagation to HF for low pt and high eta
334  if ( myPart.cos2ThetaV()>0.8 || mom.T() < 3. ) {
335  // Propagate to VFCAL entrance
336  myPart.propagateToVFcalEntrance(false);
337  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
338  myTrack.setVFcal(myPart,myPart.getSuccess());
339 
340  // Otherwise propagate to the HCAL exit and HO.
341  } else {
342  // Propagate to HCAL exit
343  myPart.propagateToHcalExit(false);
344  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
345  myTrack.setHcalExit(myPart,myPart.getSuccess());
346  // Propagate to HOLayer entrance
347  myPart.setMagneticField(0);
348  myPart.propagateToHOLayer(false);
349  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
350  myTrack.setHO(myPart,myPart.getSuccess());
351  }
352  }
353  }
354 }
355 
356 
357 void
358 FBaseSimEvent::addParticles(const HepMC::GenEvent& myGenEvent) {
359 
361  int genEventSize = myGenEvent.particles_size();
362  std::vector<int> myGenVertices(genEventSize, static_cast<int>(0));
363 
364  // If no particles, no work to be done !
365  if ( myGenEvent.particles_empty() ) return;
366 
367  // Are there particles in the FSimEvent already ?
368  int offset = nGenParts();
369 
370  // Primary vertex
371  HepMC::GenVertex* primaryVertex = *(myGenEvent.vertices_begin());
372 
373  // Beginning of workaround a bug in pythia particle gun
374  // ??? what is happening here ??? may be becuase pythia guns can have some bogus particles at the start of the event ???
375  unsigned primaryMother = primaryVertex->particles_in_size();
376  if ( primaryMother ) {
377  unsigned partId = (*(primaryVertex->particles_in_const_begin()))->pdg_id();
378  if ( abs(partId) == 2212 ) primaryMother = 0;
379  }
380  // End of workaround a bug in pythia particle gun
381 
382  XYZTLorentzVector primaryVertexPosition(primaryVertex->position().x()/10.,
383  primaryVertex->position().y()/10.,
384  primaryVertex->position().z()/10.,
385  primaryVertex->position().t()/10.);
386  // Actually this is the true end of the workaround
387  primaryVertexPosition *= (1-primaryMother);
388  // THE END.
389 
390 
391  // Set the main vertex
392  myFilter->setMainVertex(primaryVertexPosition);
393 
394  // This is the main vertex index
395  int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
396 
397  HepMC::GenEvent::particle_const_iterator piter;
398  HepMC::GenEvent::particle_const_iterator pbegin = myGenEvent.particles_begin();
399  HepMC::GenEvent::particle_const_iterator pend = myGenEvent.particles_end();
400 
401  int initialBarcode = 0;
402  if ( pbegin != pend ) initialBarcode = (*pbegin)->barcode();
403  // Loop on the particles of the generated event
404  for ( piter = pbegin; piter != pend; ++piter ) {
405 
406  // This is the generated particle pointer - for the signal event only
407  HepMC::GenParticle* p = *piter;
408 
409  if ( !offset ) {
410  (*theGenParticles)[nGenParticles++] = p;
412  theGenSize *= 2;
413  theGenParticles->resize(theGenSize);
414  }
415 
416  }
417 
418  // Reject particles with late origin vertex (i.e., coming from late decays)
419  // This should not happen, but one never knows what users may be up to!
420  // For example exotic particles might decay late - keep the decay products in the case.
421  XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
422  HepMC::GenVertex* productionVertex = p->production_vertex();
423  if ( productionVertex ) {
424  unsigned productionMother = productionVertex->particles_in_size();
425  if ( productionMother ) {
426  unsigned motherId = (*(productionVertex->particles_in_const_begin()))->pdg_id();
427  if ( abs(motherId) < 1000000 )
428  productionVertexPosition =
429  XYZTLorentzVector(productionVertex->position().x()/10.,
430  productionVertex->position().y()/10.,
431  productionVertex->position().z()/10.,
432  productionVertex->position().t()/10.);
433  }
434  }
435  if ( !myFilter->accept(productionVertexPosition) ) continue;
436 
437  int abspdgId = abs(p->pdg_id());
438  HepMC::GenVertex* endVertex = p->end_vertex();
439 
440  // Keep only:
441  // 1) Stable particles (watch out! New status code = 1001!)
442  bool testStable = p->status()%1000==1;
443  // Declare stable standard particles that decay after a macroscopic path length
444  // (except if exotic)
445  if ( p->status() == 2 && abspdgId < 1000000) {
446  if ( endVertex ) {
447  XYZTLorentzVector decayPosition =
448  XYZTLorentzVector(endVertex->position().x()/10.,
449  endVertex->position().y()/10.,
450  endVertex->position().z()/10.,
451  endVertex->position().t()/10.);
452  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
453  if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
454  }
455  }
456 
457  // 2) or particles with stable daughters (watch out! New status code = 1001!)
458  bool testDaugh = false;
459  if ( !testStable &&
460  p->status() == 2 &&
461  endVertex &&
462  endVertex->particles_out_size() ) {
463  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
464  endVertex->particles_out_const_begin();
465  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
466  endVertex->particles_out_const_end();
467  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
468  HepMC::GenParticle* daugh = *firstDaughterIt;
469  if ( daugh->status()%1000==1 ) {
470  // Check that it is not a "prompt electron or muon brem":
471  if (abspdgId == 11 || abspdgId == 13) {
472  if ( endVertex ) {
473  XYZTLorentzVector endVertexPosition = XYZTLorentzVector(endVertex->position().x()/10.,
474  endVertex->position().y()/10.,
475  endVertex->position().z()/10.,
476  endVertex->position().t()/10.);
477  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
478  if ( endVertexPosition.Perp2() < lateVertexPosition ) {
479  break;
480  }
481  }
482  }
483  testDaugh=true;
484  break;
485  }
486  }
487  }
488 
489  // 3) or particles that fly more than one micron.
490  double dist = 0.;
491  if ( !testStable && !testDaugh && p->production_vertex() ) {
493  productionVertexPosition(p->production_vertex()->position().x()/10.,
494  p->production_vertex()->position().y()/10.,
495  p->production_vertex()->position().z()/10.,
496  p->production_vertex()->position().t()/10.);
497  dist = (primaryVertexPosition-productionVertexPosition).Vect().Mag2();
498  }
499  bool testDecay = ( dist > 1e-8 ) ? true : false;
500 
501  // Save the corresponding particle and vertices
502  if ( testStable || testDaugh || testDecay ) {
503 
504  /*
505  const HepMC::GenParticle* mother = p->production_vertex() ?
506  *(p->production_vertex()->particles_in_const_begin()) : 0;
507  */
508 
509  int motherBarcode = p->production_vertex() &&
510  p->production_vertex()->particles_in_const_begin() !=
511  p->production_vertex()->particles_in_const_end() ?
512  (*(p->production_vertex()->particles_in_const_begin()))->barcode() : 0;
513 
514  int originVertex =
515  motherBarcode && myGenVertices[motherBarcode-initialBarcode] ?
516  myGenVertices[motherBarcode-initialBarcode] : mainVertex;
517 
518  XYZTLorentzVector momentum(p->momentum().px(),
519  p->momentum().py(),
520  p->momentum().pz(),
521  p->momentum().e());
522  RawParticle part(momentum, vertex(originVertex).position());
523  part.setID(p->pdg_id());
524 
525  // Add the particle to the event and to the various lists
526 
527  int theTrack = testStable && p->end_vertex() ?
528  // The particle is scheduled to decay
529  addSimTrack(&part,originVertex, nGenParts()-offset,p->end_vertex()) :
530  // The particle is not scheduled to decay
531  addSimTrack(&part,originVertex, nGenParts()-offset);
532 
533  if (
534  // This one deals with particles with no end vertex
535  !p->end_vertex() ||
536  // This one deals with particles that have a pre-defined
537  // decay proper time, but have not decayed yet
538  ( testStable && p->end_vertex() && !p->end_vertex()->particles_out_size() )
539  // In both case, just don't add a end vertex in the FSimEvent
540  ) continue;
541 
542  // Add the vertex to the event and to the various lists
543  XYZTLorentzVector decayVertex =
544  XYZTLorentzVector(p->end_vertex()->position().x()/10.,
545  p->end_vertex()->position().y()/10.,
546  p->end_vertex()->position().z()/10.,
547  p->end_vertex()->position().t()/10.);
548  // vertex(mainVertex).position();
549  int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
550 
551  if ( theVertex != -1 ) myGenVertices[p->barcode()-initialBarcode] = theVertex;
552 
553  // There we are !
554  }
555  }
556 
557 }
558 
559 void
561 
562  // If no particles, no work to be done !
563  unsigned int nParticles = myGenParticles.size();
564  nGenParticles = nParticles;
565 
566  if ( !nParticles ) return;
567 
569  std::map<const reco::Candidate*,int> myGenVertices;
570 
571  // Are there particles in the FSimEvent already ?
572  int offset = nTracks();
573 
574  // Skip the incoming protons
575  nGenParticles = 0;
576  unsigned int ip = 0;
577  if ( nParticles > 1 &&
578  myGenParticles[0].pdgId() == 2212 &&
579  myGenParticles[1].pdgId() == 2212 ) {
580  ip = 2;
581  nGenParticles = 2;
582  }
583 
584  // Primary vertex
585  XYZTLorentzVector primaryVertex (myGenParticles[ip].vx(),
586  myGenParticles[ip].vy(),
587  myGenParticles[ip].vz(),
588  0.);
589 
590  // Set the main vertex
591  myFilter->setMainVertex(primaryVertex);
592 
593  // This is the main vertex index
594  int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
595 
596  // Loop on the particles of the generated event
597  for ( ; ip<nParticles; ++ip ) {
598 
599  // nGenParticles = ip;
600 
601  nGenParticles++;
602  const reco::GenParticle& p = myGenParticles[ip];
603 
604  // Reject particles with late origin vertex (i.e., coming from late decays)
605  // This should not happen, but one never knows what users may be up to!
606  // For example exotic particles might decay late - keep the decay products in the case.
607  XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
608  const reco::Candidate* productionMother = p.numberOfMothers() ? p.mother(0) : 0;
609  if ( productionMother ) {
610  unsigned motherId = productionMother->pdgId();
611  if ( abs(motherId) < 1000000 )
612  productionVertexPosition = XYZTLorentzVector(p.vx(), p.vy(), p.vz(), 0.);
613  }
614  if ( !myFilter->accept(productionVertexPosition) ) continue;
615 
616  // Keep only:
617  // 1) Stable particles
618  bool testStable = p.status()%1000==1;
619  // Declare stable standard particles that decay after a macroscopic path length
620  // (except if exotic particle)
621  if ( p.status() == 2 && abs(p.pdgId()) < 1000000 ) {
622  unsigned int nDaughters = p.numberOfDaughters();
623  if ( nDaughters ) {
624  const reco::Candidate* daughter = p.daughter(0);
625  XYZTLorentzVector decayPosition =
626  XYZTLorentzVector(daughter->vx(), daughter->vy(), daughter->vz(), 0.);
627  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
628  if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
629  }
630  }
631 
632  // 2) or particles with stable daughters
633  bool testDaugh = false;
634  unsigned int nDaughters = p.numberOfDaughters();
635  if ( !testStable &&
636  // p.status() == 2 &&
637  nDaughters ) {
638  for ( unsigned iDaughter=0; iDaughter<nDaughters; ++iDaughter ) {
639  const reco::Candidate* daughter = p.daughter(iDaughter);
640  if ( daughter->status()%1000==1 ) {
641  testDaugh=true;
642  break;
643  }
644  }
645  }
646 
647  // 3) or particles that fly more than one micron.
648  double dist = 0.;
649  if ( !testStable && !testDaugh ) {
650  XYZTLorentzVector productionVertex(p.vx(),p.vy(),p.vz(),0.);
651  dist = (primaryVertex-productionVertex).Vect().Mag2();
652  }
653  bool testDecay = ( dist > 1e-8 ) ? true : false;
654 
655  // Save the corresponding particle and vertices
656  if ( testStable || testDaugh || testDecay ) {
657 
658  const reco::Candidate* mother = p.numberOfMothers() ? p.mother(0) : 0;
659 
660  int originVertex =
661  mother &&
662  myGenVertices.find(mother) != myGenVertices.end() ?
663  myGenVertices[mother] : mainVertex;
664 
665  XYZTLorentzVector momentum(p.px(),p.py(),p.pz(),p.energy());
666  RawParticle part(momentum, vertex(originVertex).position());
667  part.setID(p.pdgId());
668 
669  // Add the particle to the event and to the various lists
670  int theTrack = addSimTrack(&part,originVertex, nGenParts()-offset);
671 
672  // It there an end vertex ?
673  if ( !nDaughters ) continue;
674  const reco::Candidate* daughter = p.daughter(0);
675 
676  // Add the vertex to the event and to the various lists
677  XYZTLorentzVector decayVertex =
678  XYZTLorentzVector(daughter->vx(), daughter->vy(),
679  daughter->vz(), 0.);
680  int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
681 
682  if ( theVertex != -1 ) myGenVertices[&p] = theVertex;
683 
684  // There we are !
685  }
686  }
687 
688  // There is no GenParticle's in that case...
689  // nGenParticles=0;
690 
691 }
692 
693 int
694 FBaseSimEvent::addSimTrack(const RawParticle* p, int iv, int ig,
695  const HepMC::GenVertex* ev) {
696 
697  // Check that the particle is in the Famos "acceptance"
698  // Keep all primaries of pile-up events, though
699  if ( !myFilter->accept(p) && ig >= -1 ) return -1;
700 
701  // The new track index
702  int trackId = nSimTracks++;
704  theTrackSize *= 2;
705  theSimTracks->resize(theTrackSize);
706  }
707 
708  // Attach the particle to the origin vertex, and to the mother
709  vertex(iv).addDaughter(trackId);
710  if ( !vertex(iv).noParent() ) {
711  track(vertex(iv).parent().id()).addDaughter(trackId);
712 
713  if ( ig == -1 ) {
714  int motherId = track(vertex(iv).parent().id()).genpartIndex();
715  if ( motherId < -1 ) ig = motherId;
716  }
717  }
718 
719  // Some transient information for FAMOS internal use
720  (*theSimTracks)[trackId] = ev ?
721  // A proper decay time is scheduled
722  FSimTrack(p,iv,ig,trackId,this,
723  ev->position().t()/10.
724  * p->PDGmass()
725  / std::sqrt(p->momentum().Vect().Mag2())) :
726  // No proper decay time is scheduled
727  FSimTrack(p,iv,ig,trackId,this);
728 
729  return trackId;
730 
731 }
732 
733 int
735 
736  // Check that the vertex is in the Famos "acceptance"
737  if ( !myFilter->accept(v) ) return -1;
738 
739  // The number of vertices
740  int vertexId = nSimVertices++;
742  theVertexSize *= 2;
743  theSimVertices->resize(theVertexSize);
745  }
746 
747  // Attach the end vertex to the particle (if accepted)
748  if ( im !=-1 ) track(im).setEndVertex(vertexId);
749 
750  // Some transient information for FAMOS internal use
751  (*theSimVertices)[vertexId] = FSimVertex(v,im,vertexId,this);
752 
753  (*theFSimVerticesType)[vertexId] = FSimVertexType(type);
754 
755  return vertexId;
756 
757 }
758 
759 void
760 FBaseSimEvent::printMCTruth(const HepMC::GenEvent& myGenEvent) {
761 
762  std::cout << "Id Gen Name eta phi pT E Vtx1 "
763  << " x y z "
764  << "Moth Vtx2 eta phi R Z Da1 Da2 Ecal?" << std::endl;
765 
766  for ( HepMC::GenEvent::particle_const_iterator
767  piter = myGenEvent.particles_begin();
768  piter != myGenEvent.particles_end();
769  ++piter ) {
770 
771  HepMC::GenParticle* p = *piter;
772  /* */
773  int partId = p->pdg_id();
775 
776  if ( pdt->particle(ParticleID(partId)) !=0 ) {
777  name = (pdt->particle(ParticleID(partId)))->name();
778  } else {
779  name = "none";
780  }
781 
782  XYZTLorentzVector momentum1(p->momentum().px(),
783  p->momentum().py(),
784  p->momentum().pz(),
785  p->momentum().e());
786 
787  int vertexId1 = 0;
788 
789  if ( !p->production_vertex() ) continue;
790 
791  XYZVector vertex1 (p->production_vertex()->position().x()/10.,
792  p->production_vertex()->position().y()/10.,
793  p->production_vertex()->position().z()/10.);
794  vertexId1 = p->production_vertex()->barcode();
795 
796  std::cout.setf(std::ios::fixed, std::ios::floatfield);
797  std::cout.setf(std::ios::right, std::ios::adjustfield);
798 
799  std::cout << std::setw(4) << p->barcode() << " "
800  << name;
801 
802  for(unsigned int k=0;k<11-name.length() && k<12; k++) std::cout << " ";
803 
804  double eta = momentum1.eta();
805  if ( eta > +10. ) eta = +10.;
806  if ( eta < -10. ) eta = -10.;
807  std::cout << std::setw(6) << std::setprecision(2) << eta << " "
808  << std::setw(6) << std::setprecision(2) << momentum1.phi() << " "
809  << std::setw(7) << std::setprecision(2) << momentum1.pt() << " "
810  << std::setw(7) << std::setprecision(2) << momentum1.e() << " "
811  << std::setw(4) << vertexId1 << " "
812  << std::setw(6) << std::setprecision(1) << vertex1.x() << " "
813  << std::setw(6) << std::setprecision(1) << vertex1.y() << " "
814  << std::setw(6) << std::setprecision(1) << vertex1.z() << " ";
815 
816  const HepMC::GenParticle* mother =
817  *(p->production_vertex()->particles_in_const_begin());
818 
819  if ( mother )
820  std::cout << std::setw(4) << mother->barcode() << " ";
821  else
822  std::cout << " " ;
823 
824  if ( p->end_vertex() ) {
825  XYZTLorentzVector vertex2(p->end_vertex()->position().x()/10.,
826  p->end_vertex()->position().y()/10.,
827  p->end_vertex()->position().z()/10.,
828  p->end_vertex()->position().t()/10.);
829  int vertexId2 = p->end_vertex()->barcode();
830 
831  std::vector<const HepMC::GenParticle*> children;
832  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
833  p->end_vertex()->particles_out_const_begin();
834  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
835  p->end_vertex()->particles_out_const_end();
836  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
837  children.push_back(*firstDaughterIt);
838  }
839 
840  std::cout << std::setw(4) << vertexId2 << " "
841  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
842  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
843  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
844  << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
845  for ( unsigned id=0; id<children.size(); ++id )
846  std::cout << std::setw(4) << children[id]->barcode() << " ";
847  }
848  std::cout << std::endl;
849 
850  }
851 
852 }
853 
854 void
856 
857  std::cout << " Id Gen Name eta phi pT E Vtx1 "
858  << " x y z "
859  << "Moth Vtx2 eta phi R Z Daughters Ecal?" << std::endl;
860 
861  for( int i=0; i<(int)nTracks(); i++ )
862  std::cout << track(i) << std::endl;
863 
864  for( int i=0; i<(int)nVertices(); i++ )
865  std::cout << "i = " << i << " " << vertexType(i) << std::endl;
866 
867 
868 
869 }
870 
871 void
873 
874  nSimTracks = 0;
875  nSimVertices = 0;
876  nGenParticles = 0;
878 
879 }
880 
881 void
883  (*theChargedTracks)[nChargedParticleTracks++] = id;
886  theChargedSize *= 2;
888  }
889 }
890 
891 int
893  if (id>=0 && id<(int)nChargedParticleTracks)
894  return (*theChargedTracks)[id];
895  else
896  return -1;
897 }
898 
899 /*
900 const SimTrack &
901 FBaseSimEvent::embdTrack(int i) const {
902  return (*theSimTracks)[i].simTrack();
903 }
904 
905 const SimVertex &
906 FBaseSimEvent::embdVertex(int i) const {
907  return (*theSimVertices)[i].simVertex();
908 }
909 */
910 
911 const HepMC::GenParticle*
913  return (*theGenParticles)[i];
914 }
915 
916 /*
917 FSimTrack&
918 FBaseSimEvent::track(int id) const {
919  return (*theSimTracks)[id];
920 }
921 
922 
923 FSimVertex&
924 FBaseSimEvent::vertex(int id) const {
925  return (*theSimVertices)[id];
926 }
927 */
const ParticleDataTable * pdt
double lateVertexPosition
type
Definition: HCALResponse.h:21
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:36
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:139
int addSimVertex(const XYZTLorentzVector &decayVertex, int im=-1, FSimVertexType::VertexType type=FSimVertexType::ANY)
Add a new vertex to the Event and to the various lists.
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
float charge() const
charge
Definition: FSimTrack.h:47
int addSimTrack(const RawParticle *p, int iv, int ig=-1, const HepMC::GenVertex *ev=0)
Add a new track to the Event and to the various lists.
bool propagateToPreshowerLayer1(bool first=true)
const HepMC::GenParticle * embdGenpart(int i) const
return MC track with a given id
list parent
Definition: dbtoconf.py:74
HepPDT::ParticleDataTable ParticleDataTable
std::vector< FSimTrack > * theSimTracks
virtual int status() const
status word
double PDGmass() const
get the THEORETICAL mass
Definition: RawParticle.cc:255
std::vector< unsigned > * theChargedTracks
FBaseSimEvent(const edm::ParameterSet &kine)
Default constructor.
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:190
unsigned int theTrackSize
unsigned int theVertexSize
virtual int status() const =0
status word
KineParticleFilter * myFilter
The particle filter.
virtual double vx() const =0
x coordinate of vertex position
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:69
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
void addParticles(const HepMC::GenEvent &hev)
Add the particles and their vertices to the list.
void setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:62
void clear()
clear the FBaseSimEvent content before the next event
bool ev
void setMagneticField(double b)
Set the magnetic field.
unsigned int nSimTracks
T eta() const
virtual double vy() const
y coordinate of vertex position
unsigned int theGenSize
float float float z
bool notYetToEndVertex(const XYZTLorentzVector &pos) const
Compare the end vertex position with another position.
Definition: FSimTrack.cc:33
virtual double vy() const =0
y coordinate of vertex position
FSimVertex & vertex(int id) const
Return vertex with given Id.
unsigned int nVertices() const
Number of vertices.
Definition: FBaseSimEvent.h:93
FSimVertexType & vertexType(int id) const
Return vertex with given Id.
bool propagateToVFcalEntrance(bool first=true)
virtual double energy() const
energy
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:55
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
double t() const
vertex time
Definition: RawParticle.h:273
math::XYZVector XYZVector
int parentIndex() const
Definition: SimVertex.h:33
virtual const_iterator end() const =0
last daughter const_iterator
~FBaseSimEvent()
usual virtual destructor
virtual size_t numberOfMothers() const
number of mothers
const math::XYZTLorentzVector & position() const
Temporary (until CMSSW moves to Mathcore) - No ! Actually very useful.
Definition: FSimVertex.h:49
T sqrt(T t)
Definition: SSEVec.h:48
unsigned int nTracks() const
Number of tracks.
Definition: FBaseSimEvent.h:88
bool propagateToHcalExit(bool first=true)
virtual size_t numberOfDaughters() const
number of daughters
unsigned int nSimVertices
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
double cos2ThetaV() const
Definition: RawParticle.h:269
A FSimVertexType hold the information on the vertex origine.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
std::vector< FSimVertexType > FSimVertexTypeCollection
collection of FSimVertexType objects
const FSimVertex & vertex() const
Origin vertex.
unsigned int nGenParticles
unsigned int initialSize
std::vector< FSimVertex > * theSimVertices
unsigned int offset(bool)
virtual double vz() const
z coordinate of vertex position
bool propagateToEcalEntrance(bool first=true)
void fill(const HepMC::GenEvent &hev)
fill the FBaseSimEvent from the current HepMC::GenEvent
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
virtual int pdgId() const =0
PDG identifier.
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
virtual double px() const
x coordinate of momentum vector
unsigned int theChargedSize
part
Definition: HCALResponse.h:20
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:38
void addChargedTrack(int id)
Add an id in the vector of charged tracks id&#39;s.
void setHcalExit(const RawParticle &pp, int success)
Set the hcal exit variables.
Definition: FSimTrack.cc:90
virtual double pz() const
z coordinate of momentum vector
virtual double vz() const =0
z coordinate of vertex position
unsigned int nChargedParticleTracks
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
bool propagateToHcalEntrance(bool first=true)
void initializePdt(const HepPDT::ParticleDataTable *aPdt)
Initialize the particle data table.
void setEndVertex(int endv)
Set the end vertex.
Definition: FSimTrack.h:148
unsigned int nGenParts() const
Number of generator particles.
Definition: FBaseSimEvent.h:98
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
static int position[264][3]
Definition: ReadPGInfo.cc:509
void addDaughter(int i)
Add a RecHit for a track on a layer.
Definition: FSimTrack.h:181
void setHO(const RawParticle &pp, int success)
Set the ho variables.
Definition: FSimTrack.cc:96
int chargedTrack(int id) const
return &quot;reconstructed&quot; charged tracks index.
virtual double vx() const
x coordinate of vertex position
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:83
tuple cout
Definition: gather_cfg.py:121
void addDaughter(int i)
Definition: FSimVertex.h:46
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:76
bool noParent() const
Definition: SimVertex.h:34
void printMCTruth(const HepMC::GenEvent &hev)
print the original MCTruth event
Definition: DDAxes.h:10
std::vector< HepMC::GenParticle * > * theGenParticles
bool propagateToHOLayer(bool first=true)
bool propagateToPreshowerLayer2(bool first=true)
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
virtual double py() const
y coordinate of momentum vector
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void print() const
print the FBaseSimEvent in an intelligible way
FSimTrack & track(int id) const
Return track with given Id.
FSimVertexTypeCollection * theFSimVerticesType