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
23 
25 //#include "FastSimulation/Utilities/interface/Histos.h"
26 
27 using namespace HepPDT;
28 
29 // system include
30 #include <iostream>
31 #include <iomanip>
32 #include <map>
33 #include <string>
34 
36  :
37  nSimTracks(0),
38  nSimVertices(0),
39  nGenParticles(0),
40  nChargedParticleTracks(0),
41  initialSize(5000),
42  random(0)
43 {
44 
46  theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
47 
48  // Initialize the vectors of particles and vertices
49  theGenParticles = new std::vector<HepMC::GenParticle*>();
50  theSimTracks = new std::vector<FSimTrack>;
51  theSimVertices = new std::vector<FSimVertex>;
52  theChargedTracks = new std::vector<unsigned>();
54 
55  // Reserve some size to avoid mutiple copies
56  /* */
57  theSimTracks->resize(initialSize);
58  theSimVertices->resize(initialSize);
61  theFSimVerticesType->resize(initialSize);
66  /* */
67 
68  // Initialize the Particle filter
69  myFilter = new KineParticleFilter(kine);
70 
71 }
72 
74  const edm::ParameterSet& kine,
75  const RandomEngine* engine)
76  :
77  nSimTracks(0),
78  nSimVertices(0),
79  nGenParticles(0),
80  nChargedParticleTracks(0),
81  initialSize(5000),
82  theVertexGenerator(0),
83  random(engine)
84 {
85 
86  // Initialize the vertex generator
87  std::string vtxType = vtx.getParameter<std::string>("type");
88  if ( vtxType == "Gaussian" )
90  else if ( vtxType == "Flat" )
92  else if ( vtxType == "BetaFunc" )
94  else
96  // Initialize the beam spot, if not read from the DataBase
97  theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
98 
99  // Initialize the distance from (0,0,0) after which *generated* particles are
100  // no longer considered - because the mother could have interacted before.
101  // unit : cm x cm
102  lateVertexPosition = 2.5*2.5;
103 
104  // Initialize the vectors of particles and vertices
105  theGenParticles = new std::vector<HepMC::GenParticle*>();
106  theSimTracks = new std::vector<FSimTrack>;
107  theSimVertices = new std::vector<FSimVertex>;
108  theChargedTracks = new std::vector<unsigned>();
110 
111  // Reserve some size to avoid mutiple copies
112  /* */
113  theSimTracks->resize(initialSize);
114  theSimVertices->resize(initialSize);
115  theGenParticles->resize(initialSize);
116  theChargedTracks->resize(initialSize);
117  theFSimVerticesType->resize(initialSize);
122  /* */
123 
124  // Initialize the Particle filter
125  myFilter = new KineParticleFilter(kine);
126 
127  // Get the Famos Histos pointer
128  // myHistos = Histos::instance();
129 
130  // Initialize a few histograms
131  /*
132  myHistos->book("hvtx",100,-0.1,0.1);
133  myHistos->book("hvty",100,-0.1,0.1);
134  myHistos->book("hvtz",100,-500.,500.);
135  */
136 }
137 
139 
140  // Clear the vectors
141  theGenParticles->clear();
142  theSimTracks->clear();
143  theSimVertices->clear();
144  theChargedTracks->clear();
145  theFSimVerticesType->clear();
146 
147  // Delete
148  delete theGenParticles;
149  delete theSimTracks;
150  delete theSimVertices;
151  delete theChargedTracks;
152  delete theFSimVerticesType;
153  delete myFilter;
154 
155  //Write the histograms
156  // myHistos->put("histos.root");
157  // delete myHistos;
158 }
159 
160 void
162 
163  pdt = aPdt;
164 
165 }
166 
167 /*
168 const HepPDT::ParticleDataTable*
169 FBaseSimEvent::theTable() const {
170  return pdt;
171 }
172 */
173 
174 void
175 FBaseSimEvent::fill(const HepMC::GenEvent& myGenEvent) {
176 
177  // Clear old vectors
178  clear();
179 
180  // Add the particles in the FSimEvent
181  addParticles(myGenEvent);
182 
183  /*
184  std::cout << "The MC truth! " << std::endl;
185  printMCTruth(myGenEvent);
186 
187  std::cout << std::endl << "The FAMOS event! " << std::endl;
188  print();
189  */
190 
191 }
192 
193 void
195 
196  // Clear old vectors
197  clear();
198 
199  // Add the particles in the FSimEvent
200  addParticles(myGenParticles);
201 
202 }
203 
204 void
205 FBaseSimEvent::fill(const std::vector<SimTrack>& simTracks,
206  const std::vector<SimVertex>& simVertices) {
207 
208  // Watch out there ! A SimVertex is in mm (stupid),
209  // while a FSimVertex is in cm (clever).
210 
211  clear();
212 
213  unsigned nVtx = simVertices.size();
214  unsigned nTks = simTracks.size();
215 
216  // Empty event, do nothin'
217  if ( nVtx == 0 ) return;
218 
219  // Two arrays for internal use.
220  std::vector<int> myVertices(nVtx,-1);
221  std::vector<int> myTracks(nTks,-1);
222 
223  // create a map associating geant particle id and position in the
224  // event SimTrack vector
225 
226  std::map<unsigned, unsigned> geantToIndex;
227  for( unsigned it=0; it<simTracks.size(); ++it ) {
228  geantToIndex[ simTracks[it].trackId() ] = it;
229  }
230 
231  // Create also a map associating a SimTrack with its endVertex
232  /*
233  std::map<unsigned, unsigned> endVertex;
234  for ( unsigned iv=0; iv<simVertices.size(); ++iv ) {
235  endVertex[ simVertices[iv].parentIndex() ] = iv;
236  }
237  */
238 
239  // Set the main vertex for the kine particle filter
240  // SimVertices were in mm until 110_pre2
241  // HepLorentzVector primaryVertex = simVertices[0].position()/10.;
242  // SImVertices are now in cm
243  // Also : position is copied until SimVertex switches to Mathcore.
244  // XYZTLorentzVector primaryVertex = simVertices[0].position();
245  // The next 5 lines to be then replaced by the previous line
246  XYZTLorentzVector primaryVertex(simVertices[0].position().x(),
247  simVertices[0].position().y(),
248  simVertices[0].position().z(),
249  simVertices[0].position().t());
250  //
251  myFilter->setMainVertex(primaryVertex);
252  // Add the main vertex to the list.
254  myVertices[0] = 0;
255 
256  for( unsigned trackId=0; trackId<nTks; ++trackId ) {
257 
258  // The track
259  const SimTrack& track = simTracks[trackId];
260  // std::cout << std::endl << "SimTrack " << trackId << " " << track << std::endl;
261 
262  // The origin vertex
263  int vertexId = track.vertIndex();
264  const SimVertex& vertex = simVertices[vertexId];
265  //std::cout << "Origin vertex " << vertexId << " " << vertex << std::endl;
266 
267  // The mother track
268  int motherId = -1;
269  if( !vertex.noParent() ) { // there is a parent to this vertex
270  // geant id of the mother
271  unsigned motherGeantId = vertex.parentIndex();
272  std::map<unsigned, unsigned >::iterator association
273  = geantToIndex.find( motherGeantId );
274  if(association != geantToIndex.end() )
275  motherId = association->second;
276  }
277  int originId = motherId == - 1 ? -1 : myTracks[motherId];
278  //std::cout << "Origin id " << originId << std::endl;
279 
280  /*
281  if ( endVertex.find(trackId) != endVertex.end() )
282  std::cout << "End vertex id = " << endVertex[trackId] << std::endl;
283  else
284  std::cout << "No endVertex !!! " << std::endl;
285  std::cout << "Tracker surface position " << track.trackerSurfacePosition() << std::endl;
286  */
287 
288  // Add the vertex (if it does not already exist!)
289  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
290  vertex.position().pz(),vertex.position().e());
291  if ( myVertices[vertexId] == -1 )
292  // Momentum and position are copied until SimTrack and SimVertex
293  // switch to Mathcore.
294  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
295  // The next line to be then replaced by the previous line
296  myVertices[vertexId] = addSimVertex(position,originId);
297 
298  // Add the track (with protection for brem'ing electrons and muons)
299  int motherType = motherId == -1 ? 0 : simTracks[motherId].type();
300 
301  bool notBremInDetector =
302  (abs(motherType) != 11 && abs(motherType) != 13) ||
303  motherType != track.type() ||
304  position.Perp2() < lateVertexPosition ;
305 
306  if ( notBremInDetector ) {
307  // Momentum and position are copied until SimTrack and SimVertex
308  // switch to Mathcore.
309  // RawParticle part(track.momentum(), vertex.position());
310  // The next 3 lines to be then replaced by the previous line
311  XYZTLorentzVector momentum(track.momentum().px(),track.momentum().py(),
312  track.momentum().pz(),track.momentum().e());
313  RawParticle part(momentum,position);
314  //
315  part.setID(track.type());
316  //std::cout << "Ctau = " << part.PDGcTau() << std::endl;
317  // Don't save tracks that have decayed immediately but for which no daughters
318  // were saved (probably due to cuts on E, pT and eta)
319  // if ( part.PDGcTau() > 0.1 || endVertex.find(trackId) != endVertex.end() )
320  myTracks[trackId] = addSimTrack(&part,myVertices[vertexId],track.genpartIndex());
321  if ( myTracks[trackId] >= 0 ) {
322  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
323  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
324  }
325  } else {
326 
327  myTracks[trackId] = myTracks[motherId];
328  if ( myTracks[trackId] >= 0 ) {
329  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
330  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
331  }
332  }
333 
334  }
335 
336  // Now loop over the remaining end vertices !
337  for( unsigned vertexId=0; vertexId<nVtx; ++vertexId ) {
338 
339  // if the vertex is already saved, just ignore.
340  if ( myVertices[vertexId] != -1 ) continue;
341 
342  // The yet unused vertex
343  const SimVertex& vertex = simVertices[vertexId];
344 
345  // The mother track
346  int motherId = -1;
347  if( !vertex.noParent() ) { // there is a parent to this vertex
348 
349  // geant id of the mother
350  unsigned motherGeantId = vertex.parentIndex();
351  std::map<unsigned, unsigned >::iterator association
352  = geantToIndex.find( motherGeantId );
353  if(association != geantToIndex.end() )
354  motherId = association->second;
355  }
356  int originId = motherId == - 1 ? -1 : myTracks[motherId];
357 
358  // Add the vertex
359  // Momentum and position are copied until SimTrack and SimVertex
360  // switch to Mathcore.
361  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
362  // The next 3 lines to be then replaced by the previous line
363  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
364  vertex.position().pz(),vertex.position().e());
365  myVertices[vertexId] = addSimVertex(position,originId);
366  }
367 
368  // Finally, propagate all particles to the calorimeters
369  BaseParticlePropagator myPart;
370  XYZTLorentzVector mom;
372 
373  // Loop over the tracks
374  for( int fsimi=0; fsimi < (int)nTracks() ; ++fsimi) {
375 
376 
377  FSimTrack& myTrack = track(fsimi);
378  double trackerSurfaceTime = myTrack.vertex().position().t()
379  + myTrack.momentum().e()/myTrack.momentum().pz()
380  * ( myTrack.trackerSurfacePosition().z()
381  - myTrack.vertex().position().z() );
382  pos = XYZTLorentzVector(myTrack.trackerSurfacePosition().x(),
383  myTrack.trackerSurfacePosition().y(),
384  myTrack.trackerSurfacePosition().z(),
385  trackerSurfaceTime);
386  mom = XYZTLorentzVector(myTrack.trackerSurfaceMomentum().x(),
387  myTrack.trackerSurfaceMomentum().y(),
388  myTrack.trackerSurfaceMomentum().z(),
389  myTrack.trackerSurfaceMomentum().t());
390 
391  if ( mom.T() > 0. ) {
392  // The particle to be propagated
393  myPart = BaseParticlePropagator(RawParticle(mom,pos),0.,0.,4.);
394  myPart.setCharge(myTrack.charge());
395 
396  // Propagate to Preshower layer 1
397  myPart.propagateToPreshowerLayer1(false);
398  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
399  myTrack.setLayer1(myPart,myPart.getSuccess());
400 
401  // Propagate to Preshower Layer 2
402  myPart.propagateToPreshowerLayer2(false);
403  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
404  myTrack.setLayer2(myPart,myPart.getSuccess());
405 
406  // Propagate to Ecal Endcap
407  myPart.propagateToEcalEntrance(false);
408  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
409  myTrack.setEcal(myPart,myPart.getSuccess());
410 
411  // Propagate to HCAL entrance
412  myPart.propagateToHcalEntrance(false);
413  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
414  myTrack.setHcal(myPart,myPart.getSuccess());
415 
416  // Propagate to VFCAL entrance
417  myPart.propagateToVFcalEntrance(false);
418  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
419  myTrack.setVFcal(myPart,myPart.getSuccess());
420  }
421 
422  }
423 
424 }
425 
426 
427 void
428 FBaseSimEvent::addParticles(const HepMC::GenEvent& myGenEvent) {
429 
431  int genEventSize = myGenEvent.particles_size();
432  std::vector<int> myGenVertices(genEventSize, static_cast<int>(0));
433 
434  // If no particles, no work to be done !
435  if ( myGenEvent.particles_empty() ) return;
436 
437  // Are there particles in the FSimEvent already ?
438  int offset = nGenParts();
439 
440  // Primary vertex (already smeared by the SmearedVtx module)
441  HepMC::GenVertex* primaryVertex = *(myGenEvent.vertices_begin());
442 
443  // Beginning of workaround a bug in pythia particle gun
444  unsigned primaryMother = primaryVertex->particles_in_size();
445  if ( primaryMother ) {
446  unsigned partId = (*(primaryVertex->particles_in_const_begin()))->pdg_id();
447  if ( abs(partId) == 2212 ) primaryMother = 0;
448  }
449  // End of workaround a bug in pythia particle gun
450 
451  XYZTLorentzVector primaryVertexPosition(primaryVertex->position().x()/10.,
452  primaryVertex->position().y()/10.,
453  primaryVertex->position().z()/10.,
454  primaryVertex->position().t()/10.);
455  // Actually this is the true end of the workaround
456  primaryVertexPosition *= (1-primaryMother);
457  // THE END.
458 
459  // Smear the main vertex if needed
460  // Now takes the origin from the database
461  XYZTLorentzVector smearedVertex;
462  if ( primaryVertexPosition.Vect().Mag2() < 1E-16 ) {
464  smearedVertex = XYZTLorentzVector(
468  0.);
469  }
470 
471  // Set the main vertex
472  myFilter->setMainVertex(primaryVertexPosition+smearedVertex);
473 
474  // This is the smeared main vertex
475  int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
476 
477  HepMC::GenEvent::particle_const_iterator piter;
478  HepMC::GenEvent::particle_const_iterator pbegin = myGenEvent.particles_begin();
479  HepMC::GenEvent::particle_const_iterator pend = myGenEvent.particles_end();
480 
481  int initialBarcode = 0;
482  if ( pbegin != pend ) initialBarcode = (*pbegin)->barcode();
483  // Loop on the particles of the generated event
484  for ( piter = pbegin; piter != pend; ++piter ) {
485 
486  // This is the generated particle pointer - for the signal event only
487  HepMC::GenParticle* p = *piter;
488 
489  if ( !offset ) {
490  (*theGenParticles)[nGenParticles++] = p;
492  theGenSize *= 2;
493  theGenParticles->resize(theGenSize);
494  }
495 
496  }
497 
498  // Reject particles with late origin vertex (i.e., coming from late decays)
499  // This should not happen, but one never knows what users may be up to!
500  // For example exotic particles might decay late - keep the decay products in the case.
501  XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
502  HepMC::GenVertex* productionVertex = p->production_vertex();
503  if ( productionVertex ) {
504  unsigned productionMother = productionVertex->particles_in_size();
505  if ( productionMother ) {
506  unsigned motherId = (*(productionVertex->particles_in_const_begin()))->pdg_id();
507  if ( abs(motherId) < 1000000 )
508  productionVertexPosition =
509  XYZTLorentzVector(productionVertex->position().x()/10.,
510  productionVertex->position().y()/10.,
511  productionVertex->position().z()/10.,
512  productionVertex->position().t()/10.) + smearedVertex;
513  }
514  }
515  if ( !myFilter->accept(productionVertexPosition) ) continue;
516 
517  int abspdgId = abs(p->pdg_id());
518  HepMC::GenVertex* endVertex = p->end_vertex();
519 
520  // Keep only:
521  // 1) Stable particles (watch out! New status code = 1001!)
522  bool testStable = p->status()%1000==1;
523  // Declare stable standard particles that decay after a macroscopic path length
524  // (except if exotic)
525  if ( p->status() == 2 && abspdgId < 1000000) {
526  if ( endVertex ) {
527  XYZTLorentzVector decayPosition =
528  XYZTLorentzVector(endVertex->position().x()/10.,
529  endVertex->position().y()/10.,
530  endVertex->position().z()/10.,
531  endVertex->position().t()/10.) + smearedVertex;
532  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
533  if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
534  }
535  }
536 
537  // 2) or particles with stable daughters (watch out! New status code = 1001!)
538  bool testDaugh = false;
539  if ( !testStable &&
540  p->status() == 2 &&
541  endVertex &&
542  endVertex->particles_out_size() ) {
543  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
544  endVertex->particles_out_const_begin();
545  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
546  endVertex->particles_out_const_end();
547  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
548  HepMC::GenParticle* daugh = *firstDaughterIt;
549  if ( daugh->status()%1000==1 ) {
550  // Check that it is not a "prompt electron or muon brem":
551  if (abspdgId == 11 || abspdgId == 13) {
552  if ( endVertex ) {
553  XYZTLorentzVector endVertexPosition = XYZTLorentzVector(endVertex->position().x()/10.,
554  endVertex->position().y()/10.,
555  endVertex->position().z()/10.,
556  endVertex->position().t()/10.);
557  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
558  if ( endVertexPosition.Perp2() < lateVertexPosition ) {
559  break;
560  }
561  }
562  }
563  testDaugh=true;
564  break;
565  }
566  }
567  }
568 
569  // 3) or particles that fly more than one micron.
570  double dist = 0.;
571  if ( !testStable && !testDaugh && p->production_vertex() ) {
573  productionVertexPosition(p->production_vertex()->position().x()/10.,
574  p->production_vertex()->position().y()/10.,
575  p->production_vertex()->position().z()/10.,
576  p->production_vertex()->position().t()/10.);
577  dist = (primaryVertexPosition-productionVertexPosition).Vect().Mag2();
578  }
579  bool testDecay = ( dist > 1e-8 ) ? true : false;
580 
581  // Save the corresponding particle and vertices
582  if ( testStable || testDaugh || testDecay ) {
583 
584  /*
585  const HepMC::GenParticle* mother = p->production_vertex() ?
586  *(p->production_vertex()->particles_in_const_begin()) : 0;
587  */
588 
589  int motherBarcode = p->production_vertex() &&
590  p->production_vertex()->particles_in_const_begin() !=
591  p->production_vertex()->particles_in_const_end() ?
592  (*(p->production_vertex()->particles_in_const_begin()))->barcode() : 0;
593 
594  int originVertex =
595  motherBarcode && myGenVertices[motherBarcode-initialBarcode] ?
596  myGenVertices[motherBarcode-initialBarcode] : mainVertex;
597 
598  XYZTLorentzVector momentum(p->momentum().px(),
599  p->momentum().py(),
600  p->momentum().pz(),
601  p->momentum().e());
602  RawParticle part(momentum, vertex(originVertex).position());
603  part.setID(p->pdg_id());
604 
605  // Add the particle to the event and to the various lists
606 
607  int theTrack = testStable && p->end_vertex() ?
608  // The particle is scheduled to decay
609  addSimTrack(&part,originVertex, nGenParts()-offset,p->end_vertex()) :
610  // The particle is not scheduled to decay
611  addSimTrack(&part,originVertex, nGenParts()-offset);
612 
613  if (
614  // This one deals with particles with no end vertex
615  !p->end_vertex() ||
616  // This one deals with particles that have a pre-defined
617  // decay proper time, but have not decayed yet
618  ( testStable && p->end_vertex() && !p->end_vertex()->particles_out_size() )
619  // In both case, just don't add a end vertex in the FSimEvent
620  ) continue;
621 
622  // Add the vertex to the event and to the various lists
623  XYZTLorentzVector decayVertex =
624  XYZTLorentzVector(p->end_vertex()->position().x()/10.,
625  p->end_vertex()->position().y()/10.,
626  p->end_vertex()->position().z()/10.,
627  p->end_vertex()->position().t()/10.) +
628  smearedVertex;
629  // vertex(mainVertex).position();
630  int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
631 
632  if ( theVertex != -1 ) myGenVertices[p->barcode()-initialBarcode] = theVertex;
633 
634  // There we are !
635  }
636  }
637 
638 }
639 
640 void
642 
643  // If no particles, no work to be done !
644  unsigned int nParticles = myGenParticles.size();
646 
647  if ( !nParticles ) return;
648 
650  std::map<const reco::Candidate*,int> myGenVertices;
651 
652  // Are there particles in the FSimEvent already ?
653  int offset = nTracks();
654 
655  // Skip the incoming protons
656  nGenParticles = 0;
657  unsigned int ip = 0;
658  if ( nParticles > 1 &&
659  myGenParticles[0].pdgId() == 2212 &&
660  myGenParticles[1].pdgId() == 2212 ) {
661  ip = 2;
662  nGenParticles = 2;
663  }
664 
665  // Primary vertex (already smeared by the SmearedVtx module)
666  XYZTLorentzVector primaryVertex (myGenParticles[ip].vx(),
667  myGenParticles[ip].vy(),
668  myGenParticles[ip].vz(),
669  0.);
670 
671  // Smear the main vertex if needed
672  XYZTLorentzVector smearedVertex;
673  if ( primaryVertex.mag() < 1E-8 ) {
675  smearedVertex = XYZTLorentzVector(
679  0.);
680  }
681 
682  // Set the main vertex
683  myFilter->setMainVertex(primaryVertex+smearedVertex);
684 
685  // This is the smeared main vertex
686  int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
687 
688  // Loop on the particles of the generated event
689  for ( ; ip<nParticles; ++ip ) {
690 
691  // nGenParticles = ip;
692 
693  nGenParticles++;
694  const reco::GenParticle& p = myGenParticles[ip];
695 
696  // Reject particles with late origin vertex (i.e., coming from late decays)
697  // This should not happen, but one never knows what users may be up to!
698  // For example exotic particles might decay late - keep the decay products in the case.
699  XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
700  const reco::Candidate* productionMother = p.numberOfMothers() ? p.mother(0) : 0;
701  if ( productionMother ) {
702  unsigned motherId = productionMother->pdgId();
703  if ( abs(motherId) < 1000000 )
704  productionVertexPosition = XYZTLorentzVector(p.vx(), p.vy(), p.vz(), 0.) + smearedVertex;
705  }
706  if ( !myFilter->accept(productionVertexPosition) ) continue;
707 
708  // Keep only:
709  // 1) Stable particles
710  bool testStable = p.status()%1000==1;
711  // Declare stable standard particles that decay after a macroscopic path length
712  // (except if exotic particle)
713  if ( p.status() == 2 && abs(p.pdgId()) < 1000000 ) {
714  unsigned int nDaughters = p.numberOfDaughters();
715  if ( nDaughters ) {
716  const reco::Candidate* daughter = p.daughter(0);
717  XYZTLorentzVector decayPosition =
718  XYZTLorentzVector(daughter->vx(), daughter->vy(), daughter->vz(), 0.) + smearedVertex;
719  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
720  if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
721  }
722  }
723 
724  // 2) or particles with stable daughters
725  bool testDaugh = false;
726  unsigned int nDaughters = p.numberOfDaughters();
727  if ( !testStable &&
728  // p.status() == 2 &&
729  nDaughters ) {
730  for ( unsigned iDaughter=0; iDaughter<nDaughters; ++iDaughter ) {
731  const reco::Candidate* daughter = p.daughter(iDaughter);
732  if ( daughter->status()%1000==1 ) {
733  testDaugh=true;
734  break;
735  }
736  }
737  }
738 
739  // 3) or particles that fly more than one micron.
740  double dist = 0.;
741  if ( !testStable && !testDaugh ) {
742  XYZTLorentzVector productionVertex(p.vx(),p.vy(),p.vz(),0.);
743  dist = (primaryVertex-productionVertex).Vect().Mag2();
744  }
745  bool testDecay = ( dist > 1e-8 ) ? true : false;
746 
747  // Save the corresponding particle and vertices
748  if ( testStable || testDaugh || testDecay ) {
749 
750  const reco::Candidate* mother = p.numberOfMothers() ? p.mother(0) : 0;
751 
752  int originVertex =
753  mother &&
754  myGenVertices.find(mother) != myGenVertices.end() ?
755  myGenVertices[mother] : mainVertex;
756 
757  XYZTLorentzVector momentum(p.px(),p.py(),p.pz(),p.energy());
758  RawParticle part(momentum, vertex(originVertex).position());
759  part.setID(p.pdgId());
760 
761  // Add the particle to the event and to the various lists
762  int theTrack = addSimTrack(&part,originVertex, nGenParts()-offset);
763 
764  // It there an end vertex ?
765  if ( !nDaughters ) continue;
766  const reco::Candidate* daughter = p.daughter(0);
767 
768  // Add the vertex to the event and to the various lists
769  XYZTLorentzVector decayVertex =
770  XYZTLorentzVector(daughter->vx(), daughter->vy(),
771  daughter->vz(), 0.) + smearedVertex;
772  int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
773 
774  if ( theVertex != -1 ) myGenVertices[&p] = theVertex;
775 
776  // There we are !
777  }
778  }
779 
780  // There is no GenParticle's in that case...
781  // nGenParticles=0;
782 
783 }
784 
785 int
786 FBaseSimEvent::addSimTrack(const RawParticle* p, int iv, int ig,
787  const HepMC::GenVertex* ev) {
788 
789  // Check that the particle is in the Famos "acceptance"
790  // Keep all primaries of pile-up events, though
791  if ( !myFilter->accept(p) && ig >= -1 ) return -1;
792 
793  // The new track index
794  int trackId = nSimTracks++;
796  theTrackSize *= 2;
797  theSimTracks->resize(theTrackSize);
798  }
799 
800  // Attach the particle to the origin vertex, and to the mother
801  vertex(iv).addDaughter(trackId);
802  if ( !vertex(iv).noParent() ) {
803  track(vertex(iv).parent().id()).addDaughter(trackId);
804 
805  if ( ig == -1 ) {
806  int motherId = track(vertex(iv).parent().id()).genpartIndex();
807  if ( motherId < -1 ) ig = motherId;
808  }
809  }
810 
811  // Some transient information for FAMOS internal use
812  (*theSimTracks)[trackId] = ev ?
813  // A proper decay time is scheduled
814  FSimTrack(p,iv,ig,trackId,this,
815  ev->position().t()/10.
816  * p->PDGmass()
817  / std::sqrt(p->momentum().Vect().Mag2())) :
818  // No proper decay time is scheduled
819  FSimTrack(p,iv,ig,trackId,this);
820 
821  return trackId;
822 
823 }
824 
825 int
827 
828  // Check that the vertex is in the Famos "acceptance"
829  if ( !myFilter->accept(v) ) return -1;
830 
831  // The number of vertices
832  int vertexId = nSimVertices++;
834  theVertexSize *= 2;
835  theSimVertices->resize(theVertexSize);
837  }
838 
839  // Attach the end vertex to the particle (if accepted)
840  if ( im !=-1 ) track(im).setEndVertex(vertexId);
841 
842  // Some transient information for FAMOS internal use
843  (*theSimVertices)[vertexId] = FSimVertex(v,im,vertexId,this);
844 
845  (*theFSimVerticesType)[vertexId] = FSimVertexType(type);
846 
847  return vertexId;
848 
849 }
850 
851 void
852 FBaseSimEvent::printMCTruth(const HepMC::GenEvent& myGenEvent) {
853 
854  std::cout << "Id Gen Name eta phi pT E Vtx1 "
855  << " x y z "
856  << "Moth Vtx2 eta phi R Z Da1 Da2 Ecal?" << std::endl;
857 
858  for ( HepMC::GenEvent::particle_const_iterator
859  piter = myGenEvent.particles_begin();
860  piter != myGenEvent.particles_end();
861  ++piter ) {
862 
863  HepMC::GenParticle* p = *piter;
864  /* */
865  int partId = p->pdg_id();
866  std::string name;
867 
868  if ( pdt->particle(ParticleID(partId)) !=0 ) {
869  name = (pdt->particle(ParticleID(partId)))->name();
870  } else {
871  name = "none";
872  }
873 
874  XYZTLorentzVector momentum1(p->momentum().px(),
875  p->momentum().py(),
876  p->momentum().pz(),
877  p->momentum().e());
878 
879  int vertexId1 = 0;
880 
881  if ( !p->production_vertex() ) continue;
882 
883  XYZVector vertex1 (p->production_vertex()->position().x()/10.,
884  p->production_vertex()->position().y()/10.,
885  p->production_vertex()->position().z()/10.);
886  vertexId1 = p->production_vertex()->barcode();
887 
888  std::cout.setf(std::ios::fixed, std::ios::floatfield);
889  std::cout.setf(std::ios::right, std::ios::adjustfield);
890 
891  std::cout << std::setw(4) << p->barcode() << " "
892  << name;
893 
894  for(unsigned int k=0;k<11-name.length() && k<12; k++) std::cout << " ";
895 
896  double eta = momentum1.eta();
897  if ( eta > +10. ) eta = +10.;
898  if ( eta < -10. ) eta = -10.;
899  std::cout << std::setw(6) << std::setprecision(2) << eta << " "
900  << std::setw(6) << std::setprecision(2) << momentum1.phi() << " "
901  << std::setw(7) << std::setprecision(2) << momentum1.pt() << " "
902  << std::setw(7) << std::setprecision(2) << momentum1.e() << " "
903  << std::setw(4) << vertexId1 << " "
904  << std::setw(6) << std::setprecision(1) << vertex1.x() << " "
905  << std::setw(6) << std::setprecision(1) << vertex1.y() << " "
906  << std::setw(6) << std::setprecision(1) << vertex1.z() << " ";
907 
908  const HepMC::GenParticle* mother =
909  *(p->production_vertex()->particles_in_const_begin());
910 
911  if ( mother )
912  std::cout << std::setw(4) << mother->barcode() << " ";
913  else
914  std::cout << " " ;
915 
916  if ( p->end_vertex() ) {
917  XYZTLorentzVector vertex2(p->end_vertex()->position().x()/10.,
918  p->end_vertex()->position().y()/10.,
919  p->end_vertex()->position().z()/10.,
920  p->end_vertex()->position().t()/10.);
921  int vertexId2 = p->end_vertex()->barcode();
922 
923  std::vector<const HepMC::GenParticle*> children;
924  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
925  p->end_vertex()->particles_out_const_begin();
926  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
927  p->end_vertex()->particles_out_const_end();
928  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
929  children.push_back(*firstDaughterIt);
930  }
931 
932  std::cout << std::setw(4) << vertexId2 << " "
933  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
934  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
935  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
936  << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
937  for ( unsigned id=0; id<children.size(); ++id )
938  std::cout << std::setw(4) << children[id]->barcode() << " ";
939  }
940  std::cout << std::endl;
941 
942  }
943 
944 }
945 
946 void
948 
949  std::cout << " Id Gen Name eta phi pT E Vtx1 "
950  << " x y z "
951  << "Moth Vtx2 eta phi R Z Daughters Ecal?" << std::endl;
952 
953  for( int i=0; i<(int)nTracks(); i++ )
954  std::cout << track(i) << std::endl;
955 
956  for( int i=0; i<(int)nVertices(); i++ )
957  std::cout << "i = " << i << " " << vertexType(i) << std::endl;
958 
959 
960 
961 }
962 
963 void
965 
966  nSimTracks = 0;
967  nSimVertices = 0;
968  nGenParticles = 0;
970 
971 }
972 
973 void
975  (*theChargedTracks)[nChargedParticleTracks++] = id;
978  theChargedSize *= 2;
980  }
981 }
982 
983 int
985  if (id>=0 && id<(int)nChargedParticleTracks)
986  return (*theChargedTracks)[id];
987  else
988  return -1;
989 }
990 
991 /*
992 const SimTrack &
993 FBaseSimEvent::embdTrack(int i) const {
994  return (*theSimTracks)[i].simTrack();
995 }
996 
997 const SimVertex &
998 FBaseSimEvent::embdVertex(int i) const {
999  return (*theSimVertices)[i].simVertex();
1000 }
1001 */
1002 
1003 const HepMC::GenParticle*
1005  return (*theGenParticles)[i];
1006 }
1007 
1008 /*
1009 FSimTrack&
1010 FBaseSimEvent::track(int id) const {
1011  return (*theSimTracks)[id];
1012 }
1013 
1014 
1015 FSimVertex&
1016 FBaseSimEvent::vertex(int id) const {
1017  return (*theSimVertices)[id];
1018 }
1019 */
const ParticleDataTable * pdt
double lateVertexPosition
PrimaryVertexGenerator * theVertexGenerator
type
Definition: HCALResponse.h:22
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:138
T getParameter(std::string const &) const
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:254
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:168
unsigned int theTrackSize
unsigned int theVertexSize
virtual int status() const =0
status word
KineParticleFilter * myFilter
The particle filter.
math::XYZPoint theBeamSpot
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 ?
#define abs(x)
Definition: mlp_lapack.h:159
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
TRandom random
Definition: MVATrainer.cc:138
unsigned int nSimTracks
virtual double vy() const
y coordinate of vertex position
unsigned int theGenSize
T eta() const
bool notYetToEndVertex(const XYZTLorentzVector &pos) const
Compare the end vertex position with another position.
Definition: FSimTrack.cc:33
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
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.
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:285
double t() const
vertex time
Definition: RawParticle.h:272
math::XYZVector XYZVector
int parentIndex() const
Definition: SimVertex.h:33
Definition: DDAxes.h:10
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:28
unsigned int nTracks() const
Number of tracks.
Definition: FBaseSimEvent.h:95
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) ...
A FSimVertexType hold the information on the vertex origine.
#define nParticles
Definition: Histograms.cc:25
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)
int k[5][pyjets_maxn]
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:284
virtual void generate()=0
Generation process (to be implemented)
const RandomEngine * random
virtual double px() const
x coordinate of momentum vector
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
unsigned int theChargedSize
part
Definition: HCALResponse.h:21
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.
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:40
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:132
unsigned int nGenParts() const
Number of generator particles.
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
const math::XYZPoint & beamSpot() const
Return x0, y0, z0.
void addDaughter(int i)
Add a RecHit for a track on a layer.
Definition: FSimTrack.h:159
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:41
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
std::vector< HepMC::GenParticle * > * theGenParticles
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) ...
mathSSE::Vec4< T > v
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