test
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  // Initialize the distance from (0,0,0) after which *generated* particles are
63  // no longer considered - because the mother could have interacted before.
64  // unit : cm x cm
65  lateVertexPosition = 2.5*2.5;
66 }
67 
69 
70  // Clear the vectors
71  theGenParticles->clear();
72  theSimTracks->clear();
73  theSimVertices->clear();
74  theChargedTracks->clear();
75  theFSimVerticesType->clear();
76 
77  // Delete
78  delete theGenParticles;
79  delete theSimTracks;
80  delete theSimVertices;
81  delete theChargedTracks;
82  delete theFSimVerticesType;
83  delete myFilter;
84 
85 }
86 
87 void
89 
90  pdt = aPdt;
91 
92 }
93 
94 void
95 FBaseSimEvent::fill(const HepMC::GenEvent& myGenEvent) {
96 
97  // Clear old vectors
98  clear();
99 
100  // Add the particles in the FSimEvent
101  addParticles(myGenEvent);
102 
103 }
104 
105 void
106 FBaseSimEvent::fill(const std::vector<SimTrack>& simTracks,
107  const std::vector<SimVertex>& simVertices) {
108 
109  // Watch out there ! A SimVertex is in mm (stupid),
110  // while a FSimVertex is in cm (clever).
111 
112  clear();
113 
114  unsigned nVtx = simVertices.size();
115  unsigned nTks = simTracks.size();
116 
117  // Empty event, do nothin'
118  if ( nVtx == 0 ) return;
119 
120  // Two arrays for internal use.
121  std::vector<int> myVertices(nVtx,-1);
122  std::vector<int> myTracks(nTks,-1);
123 
124  // create a map associating geant particle id and position in the
125  // event SimTrack vector
126 
127  std::map<unsigned, unsigned> geantToIndex;
128  for( unsigned it=0; it<simTracks.size(); ++it ) {
129  geantToIndex[ simTracks[it].trackId() ] = it;
130  }
131 
132  // Create also a map associating a SimTrack with its endVertex
133  /*
134  std::map<unsigned, unsigned> endVertex;
135  for ( unsigned iv=0; iv<simVertices.size(); ++iv ) {
136  endVertex[ simVertices[iv].parentIndex() ] = iv;
137  }
138  */
139 
140  // Set the main vertex for the kine particle filter
141  // SimVertices were in mm until 110_pre2
142  // HepLorentzVector primaryVertex = simVertices[0].position()/10.;
143  // SImVertices are now in cm
144  // Also : position is copied until SimVertex switches to Mathcore.
145  // XYZTLorentzVector primaryVertex = simVertices[0].position();
146  // The next 5 lines to be then replaced by the previous line
147  XYZTLorentzVector primaryVertex(simVertices[0].position().x(),
148  simVertices[0].position().y(),
149  simVertices[0].position().z(),
150  simVertices[0].position().t());
151  //
152  myFilter->setMainVertex(primaryVertex);
153  // Add the main vertex to the list.
155  myVertices[0] = 0;
156 
157  for( unsigned trackId=0; trackId<nTks; ++trackId ) {
158 
159  // The track
160  const SimTrack& track = simTracks[trackId];
161  // std::cout << std::endl << "SimTrack " << trackId << " " << track << std::endl;
162 
163  // The origin vertex
164  int vertexId = track.vertIndex();
165  const SimVertex& vertex = simVertices[vertexId];
166  //std::cout << "Origin vertex " << vertexId << " " << vertex << std::endl;
167 
168  // The mother track
169  int motherId = -1;
170  if( !vertex.noParent() ) { // there is a parent to this vertex
171  // geant id of the mother
172  unsigned motherGeantId = vertex.parentIndex();
173  std::map<unsigned, unsigned >::iterator association
174  = geantToIndex.find( motherGeantId );
175  if(association != geantToIndex.end() )
176  motherId = association->second;
177  }
178  int originId = motherId == - 1 ? -1 : myTracks[motherId];
179  //std::cout << "Origin id " << originId << std::endl;
180 
181  /*
182  if ( endVertex.find(trackId) != endVertex.end() )
183  std::cout << "End vertex id = " << endVertex[trackId] << std::endl;
184  else
185  std::cout << "No endVertex !!! " << std::endl;
186  std::cout << "Tracker surface position " << track.trackerSurfacePosition() << std::endl;
187  */
188 
189  // Add the vertex (if it does not already exist!)
190  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
191  vertex.position().pz(),vertex.position().e());
192  if ( myVertices[vertexId] == -1 )
193  // Momentum and position are copied until SimTrack and SimVertex
194  // switch to Mathcore.
195  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
196  // The next line to be then replaced by the previous line
197  myVertices[vertexId] = addSimVertex(position,originId);
198 
199  // Add the track (with protection for brem'ing electrons and muons)
200  int motherType = motherId == -1 ? 0 : simTracks[motherId].type();
201 
202  bool notBremInDetector =
203  (abs(motherType) != 11 && abs(motherType) != 13) ||
204  motherType != track.type() ||
205  position.Perp2() < lateVertexPosition ;
206 
207  if ( notBremInDetector ) {
208  // Momentum and position are copied until SimTrack and SimVertex
209  // switch to Mathcore.
210  // RawParticle part(track.momentum(), vertex.position());
211  // The next 3 lines to be then replaced by the previous line
212  XYZTLorentzVector momentum(track.momentum().px(),track.momentum().py(),
213  track.momentum().pz(),track.momentum().e());
214  RawParticle part(momentum,position);
215  //
216  part.setID(track.type());
217  //std::cout << "Ctau = " << part.PDGcTau() << std::endl;
218  // Don't save tracks that have decayed immediately but for which no daughters
219  // were saved (probably due to cuts on E, pT and eta)
220  // if ( part.PDGcTau() > 0.1 || endVertex.find(trackId) != endVertex.end() )
221  myTracks[trackId] = addSimTrack(&part,myVertices[vertexId],track.genpartIndex());
222  if ( myTracks[trackId] >= 0 ) {
223  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
224  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
225  }
226  } else {
227 
228  myTracks[trackId] = myTracks[motherId];
229  if ( myTracks[trackId] >= 0 ) {
230  (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
231  (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
232  }
233  }
234 
235  }
236 
237  // Now loop over the remaining end vertices !
238  for( unsigned vertexId=0; vertexId<nVtx; ++vertexId ) {
239 
240  // if the vertex is already saved, just ignore.
241  if ( myVertices[vertexId] != -1 ) continue;
242 
243  // The yet unused vertex
244  const SimVertex& vertex = simVertices[vertexId];
245 
246  // The mother track
247  int motherId = -1;
248  if( !vertex.noParent() ) { // there is a parent to this vertex
249 
250  // geant id of the mother
251  unsigned motherGeantId = vertex.parentIndex();
252  std::map<unsigned, unsigned >::iterator association
253  = geantToIndex.find( motherGeantId );
254  if(association != geantToIndex.end() )
255  motherId = association->second;
256  }
257  int originId = motherId == - 1 ? -1 : myTracks[motherId];
258 
259  // Add the vertex
260  // Momentum and position are copied until SimTrack and SimVertex
261  // switch to Mathcore.
262  // myVertices[vertexId] = addSimVertex(vertex.position(),originId);
263  // The next 3 lines to be then replaced by the previous line
264  XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
265  vertex.position().pz(),vertex.position().e());
266  myVertices[vertexId] = addSimVertex(position,originId);
267  }
268 
269  // Finally, propagate all particles to the calorimeters
270  BaseParticlePropagator myPart;
271  XYZTLorentzVector mom;
272  XYZTLorentzVector pos;
273 
274 
275  // Loop over the tracks
276  for( int fsimi=0; fsimi < (int)nTracks() ; ++fsimi) {
277 
278 
279  FSimTrack& myTrack = track(fsimi);
280  double trackerSurfaceTime = myTrack.vertex().position().t()
281  + myTrack.momentum().e()/myTrack.momentum().pz()
282  * ( myTrack.trackerSurfacePosition().z()
283  - myTrack.vertex().position().z() );
284  pos = XYZTLorentzVector(myTrack.trackerSurfacePosition().x(),
285  myTrack.trackerSurfacePosition().y(),
286  myTrack.trackerSurfacePosition().z(),
287  trackerSurfaceTime);
288  mom = XYZTLorentzVector(myTrack.trackerSurfaceMomentum().x(),
289  myTrack.trackerSurfaceMomentum().y(),
290  myTrack.trackerSurfaceMomentum().z(),
291  myTrack.trackerSurfaceMomentum().t());
292 
293  if ( mom.T() > 0. ) {
294  // The particle to be propagated
295  myPart = BaseParticlePropagator(RawParticle(mom,pos),0.,0.,4.);
296  myPart.setCharge(myTrack.charge());
297 
298  // Propagate to Preshower layer 1
299  myPart.propagateToPreshowerLayer1(false);
300  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
301  myTrack.setLayer1(myPart,myPart.getSuccess());
302 
303  // Propagate to Preshower Layer 2
304  myPart.propagateToPreshowerLayer2(false);
305  if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
306  myTrack.setLayer2(myPart,myPart.getSuccess());
307 
308  // Propagate to Ecal Endcap
309  myPart.propagateToEcalEntrance(false);
310  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
311  myTrack.setEcal(myPart,myPart.getSuccess());
312 
313  // Propagate to HCAL entrance
314  myPart.propagateToHcalEntrance(false);
315  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
316  myTrack.setHcal(myPart,myPart.getSuccess());
317 
318  // Attempt propagation to HF for low pt and high eta
319  if ( myPart.cos2ThetaV()>0.8 || mom.T() < 3. ) {
320  // Propagate to VFCAL entrance
321  myPart.propagateToVFcalEntrance(false);
322  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
323  myTrack.setVFcal(myPart,myPart.getSuccess());
324 
325  // Otherwise propagate to the HCAL exit and HO.
326  } else {
327  // Propagate to HCAL exit
328  myPart.propagateToHcalExit(false);
329  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
330  myTrack.setHcalExit(myPart,myPart.getSuccess());
331  // Propagate to HOLayer entrance
332  myPart.setMagneticField(0);
333  myPart.propagateToHOLayer(false);
334  if ( myTrack.notYetToEndVertex(myPart.vertex()) )
335  myTrack.setHO(myPart,myPart.getSuccess());
336  }
337  }
338  }
339 }
340 
341 void
342 FBaseSimEvent::addParticles(const HepMC::GenEvent& myGenEvent) {
343 
345  int genEventSize = myGenEvent.particles_size();
346  std::vector<int> myGenVertices(genEventSize, static_cast<int>(0));
347 
348  // If no particles, no work to be done !
349  if ( myGenEvent.particles_empty() ) return;
350 
351  // Are there particles in the FSimEvent already ?
352  int offset = nGenParts();
353 
354  // Primary vertex
355  HepMC::GenVertex* primaryVertex = *(myGenEvent.vertices_begin());
356 
357  // unit transformation (needs review)
358  XYZTLorentzVector primaryVertexPosition(primaryVertex->position().x()/10.,
359  primaryVertex->position().y()/10.,
360  primaryVertex->position().z()/10.,
361  primaryVertex->position().t()/10.);
362 
363  // Set the main vertex
364  myFilter->setMainVertex(primaryVertexPosition);
365 
366  // This is the main vertex index
367  int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
368 
369  HepMC::GenEvent::particle_const_iterator piter;
370  HepMC::GenEvent::particle_const_iterator pbegin = myGenEvent.particles_begin();
371  HepMC::GenEvent::particle_const_iterator pend = myGenEvent.particles_end();
372 
373  int initialBarcode = 0;
374  if ( pbegin != pend ) initialBarcode = (*pbegin)->barcode();
375  // Loop on the particles of the generated event
376  for ( piter = pbegin; piter != pend; ++piter ) {
377 
378  // This is the generated particle pointer - for the signal event only
379  HepMC::GenParticle* p = *piter;
380 
381  if ( !offset ) {
382  (*theGenParticles)[nGenParticles++] = p;
384  theGenSize *= 2;
385  theGenParticles->resize(theGenSize);
386  }
387 
388  }
389 
390  // Reject particles with late origin vertex (i.e., coming from late decays)
391  // This should not happen, but one never knows what users may be up to!
392  // For example exotic particles might decay late - keep the decay products in the case.
393  XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
394  HepMC::GenVertex* productionVertex = p->production_vertex();
395  if ( productionVertex ) {
396  unsigned productionMother = productionVertex->particles_in_size();
397  if ( productionMother ) {
398  unsigned motherId = (*(productionVertex->particles_in_const_begin()))->pdg_id();
399  if ( abs(motherId) < 1000000 )
400  productionVertexPosition =
401  XYZTLorentzVector(productionVertex->position().x()/10.,
402  productionVertex->position().y()/10.,
403  productionVertex->position().z()/10.,
404  productionVertex->position().t()/10.);
405  }
406  }
407  if ( !myFilter->accept(productionVertexPosition) ) continue;
408 
409  int abspdgId = abs(p->pdg_id());
410  HepMC::GenVertex* endVertex = p->end_vertex();
411 
412  // Keep only:
413  // 1) Stable particles (watch out! New status code = 1001!)
414  bool testStable = p->status()%1000==1;
415  // Declare stable standard particles that decay after a macroscopic path length
416  // (except if exotic)
417  if ( p->status() == 2 && abspdgId < 1000000) {
418  if ( endVertex ) {
419  XYZTLorentzVector decayPosition =
420  XYZTLorentzVector(endVertex->position().x()/10.,
421  endVertex->position().y()/10.,
422  endVertex->position().z()/10.,
423  endVertex->position().t()/10.);
424  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
425  if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
426  }
427  }
428 
429  // 2) or particles with stable daughters (watch out! New status code = 1001!)
430  bool testDaugh = false;
431  if ( !testStable &&
432  p->status() == 2 &&
433  endVertex &&
434  endVertex->particles_out_size() ) {
435  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
436  endVertex->particles_out_const_begin();
437  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
438  endVertex->particles_out_const_end();
439  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
440  HepMC::GenParticle* daugh = *firstDaughterIt;
441  if ( daugh->status()%1000==1 ) {
442  // Check that it is not a "prompt electron or muon brem":
443  if (abspdgId == 11 || abspdgId == 13) {
444  if ( endVertex ) {
445  XYZTLorentzVector endVertexPosition = XYZTLorentzVector(endVertex->position().x()/10.,
446  endVertex->position().y()/10.,
447  endVertex->position().z()/10.,
448  endVertex->position().t()/10.);
449  // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
450  if ( endVertexPosition.Perp2() < lateVertexPosition ) {
451  break;
452  }
453  }
454  }
455  testDaugh=true;
456  break;
457  }
458  }
459  }
460 
461  // 3) or particles that fly more than one micron.
462  double dist = 0.;
463  if ( !testStable && !testDaugh && p->production_vertex() ) {
465  productionVertexPosition(p->production_vertex()->position().x()/10.,
466  p->production_vertex()->position().y()/10.,
467  p->production_vertex()->position().z()/10.,
468  p->production_vertex()->position().t()/10.);
469  dist = (primaryVertexPosition-productionVertexPosition).Vect().Mag2();
470  }
471  bool testDecay = ( dist > 1e-8 ) ? true : false;
472 
473  // Save the corresponding particle and vertices
474  if ( testStable || testDaugh || testDecay ) {
475 
476  /*
477  const HepMC::GenParticle* mother = p->production_vertex() ?
478  *(p->production_vertex()->particles_in_const_begin()) : 0;
479  */
480 
481  int motherBarcode = p->production_vertex() &&
482  p->production_vertex()->particles_in_const_begin() !=
483  p->production_vertex()->particles_in_const_end() ?
484  (*(p->production_vertex()->particles_in_const_begin()))->barcode() : 0;
485 
486  int originVertex =
487  motherBarcode && myGenVertices[motherBarcode-initialBarcode] ?
488  myGenVertices[motherBarcode-initialBarcode] : mainVertex;
489 
490  XYZTLorentzVector momentum(p->momentum().px(),
491  p->momentum().py(),
492  p->momentum().pz(),
493  p->momentum().e());
494  RawParticle part(momentum, vertex(originVertex).position());
495  part.setID(p->pdg_id());
496 
497  // Add the particle to the event and to the various lists
498 
499  int theTrack = testStable && p->end_vertex() ?
500  // The particle is scheduled to decay
501  addSimTrack(&part,originVertex, nGenParts()-offset,p->end_vertex()) :
502  // The particle is not scheduled to decay
503  addSimTrack(&part,originVertex, nGenParts()-offset);
504 
505  if (
506  // This one deals with particles with no end vertex
507  !p->end_vertex() ||
508  // This one deals with particles that have a pre-defined
509  // decay proper time, but have not decayed yet
510  ( testStable && p->end_vertex() && !p->end_vertex()->particles_out_size() )
511  // In both case, just don't add a end vertex in the FSimEvent
512  ) continue;
513 
514  // Add the vertex to the event and to the various lists
515  XYZTLorentzVector decayVertex =
516  XYZTLorentzVector(p->end_vertex()->position().x()/10.,
517  p->end_vertex()->position().y()/10.,
518  p->end_vertex()->position().z()/10.,
519  p->end_vertex()->position().t()/10.);
520  // vertex(mainVertex).position();
521  int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
522 
523  if ( theVertex != -1 ) myGenVertices[p->barcode()-initialBarcode] = theVertex;
524 
525  // There we are !
526  }
527  }
528 
529 }
530 
531 
532 int
533 FBaseSimEvent::addSimTrack(const RawParticle* p, int iv, int ig,
534  const HepMC::GenVertex* ev) {
535 
536  // Check that the particle is in the Famos "acceptance"
537  // Keep all primaries of pile-up events, though
538  if ( !myFilter->accept(p) && ig >= -1 ) return -1;
539 
540  // The new track index
541  int trackId = nSimTracks++;
543  theTrackSize *= 2;
544  theSimTracks->resize(theTrackSize);
545  }
546 
547  // Attach the particle to the origin vertex, and to the mother
548  vertex(iv).addDaughter(trackId);
549  if ( !vertex(iv).noParent() ) {
550  track(vertex(iv).parent().id()).addDaughter(trackId);
551 
552  if ( ig == -1 ) {
553  int motherId = track(vertex(iv).parent().id()).genpartIndex();
554  if ( motherId < -1 ) ig = motherId;
555  }
556  }
557 
558  // Some transient information for FAMOS internal use
559  (*theSimTracks)[trackId] = ev ?
560  // A proper decay time is scheduled
561  FSimTrack(p,iv,ig,trackId,this,
562  ev->position().t()/10.
563  * p->PDGmass()
564  / std::sqrt(p->momentum().Vect().Mag2())) :
565  // No proper decay time is scheduled
566  FSimTrack(p,iv,ig,trackId,this);
567 
568  return trackId;
569 
570 }
571 
572 int
574 
575  // Check that the vertex is in the Famos "acceptance"
576  if ( !myFilter->accept(v) ) return -1;
577 
578  // The number of vertices
579  int vertexId = nSimVertices++;
581  theVertexSize *= 2;
582  theSimVertices->resize(theVertexSize);
584  }
585 
586  // Attach the end vertex to the particle (if accepted)
587  if ( im !=-1 ) track(im).setEndVertex(vertexId);
588 
589  // Some transient information for FAMOS internal use
590  (*theSimVertices)[vertexId] = FSimVertex(v,im,vertexId,this);
591 
592  (*theFSimVerticesType)[vertexId] = FSimVertexType(type);
593 
594  return vertexId;
595 
596 }
597 
598 void
599 FBaseSimEvent::printMCTruth(const HepMC::GenEvent& myGenEvent) {
600 
601  std::cout << "Id Gen Name eta phi pT E Vtx1 "
602  << " x y z "
603  << "Moth Vtx2 eta phi R Z Da1 Da2 Ecal?" << std::endl;
604 
605  for ( HepMC::GenEvent::particle_const_iterator
606  piter = myGenEvent.particles_begin();
607  piter != myGenEvent.particles_end();
608  ++piter ) {
609 
610  HepMC::GenParticle* p = *piter;
611  /* */
612  int partId = p->pdg_id();
614 
615  if ( pdt->particle(ParticleID(partId)) !=0 ) {
616  name = (pdt->particle(ParticleID(partId)))->name();
617  } else {
618  name = "none";
619  }
620 
621  XYZTLorentzVector momentum1(p->momentum().px(),
622  p->momentum().py(),
623  p->momentum().pz(),
624  p->momentum().e());
625 
626  int vertexId1 = 0;
627 
628  if ( !p->production_vertex() ) continue;
629 
630  XYZVector vertex1 (p->production_vertex()->position().x()/10.,
631  p->production_vertex()->position().y()/10.,
632  p->production_vertex()->position().z()/10.);
633  vertexId1 = p->production_vertex()->barcode();
634 
635  std::cout.setf(std::ios::fixed, std::ios::floatfield);
636  std::cout.setf(std::ios::right, std::ios::adjustfield);
637 
638  std::cout << std::setw(4) << p->barcode() << " "
639  << name;
640 
641  for(unsigned int k=0;k<11-name.length() && k<12; k++) std::cout << " ";
642 
643  double eta = momentum1.eta();
644  if ( eta > +10. ) eta = +10.;
645  if ( eta < -10. ) eta = -10.;
646  std::cout << std::setw(6) << std::setprecision(2) << eta << " "
647  << std::setw(6) << std::setprecision(2) << momentum1.phi() << " "
648  << std::setw(7) << std::setprecision(2) << momentum1.pt() << " "
649  << std::setw(7) << std::setprecision(2) << momentum1.e() << " "
650  << std::setw(4) << vertexId1 << " "
651  << std::setw(6) << std::setprecision(1) << vertex1.x() << " "
652  << std::setw(6) << std::setprecision(1) << vertex1.y() << " "
653  << std::setw(6) << std::setprecision(1) << vertex1.z() << " ";
654 
655  const HepMC::GenParticle* mother =
656  *(p->production_vertex()->particles_in_const_begin());
657 
658  if ( mother )
659  std::cout << std::setw(4) << mother->barcode() << " ";
660  else
661  std::cout << " " ;
662 
663  if ( p->end_vertex() ) {
664  XYZTLorentzVector vertex2(p->end_vertex()->position().x()/10.,
665  p->end_vertex()->position().y()/10.,
666  p->end_vertex()->position().z()/10.,
667  p->end_vertex()->position().t()/10.);
668  int vertexId2 = p->end_vertex()->barcode();
669 
670  std::vector<const HepMC::GenParticle*> children;
671  HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
672  p->end_vertex()->particles_out_const_begin();
673  HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
674  p->end_vertex()->particles_out_const_end();
675  for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
676  children.push_back(*firstDaughterIt);
677  }
678 
679  std::cout << std::setw(4) << vertexId2 << " "
680  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
681  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
682  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
683  << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
684  for ( unsigned id=0; id<children.size(); ++id )
685  std::cout << std::setw(4) << children[id]->barcode() << " ";
686  }
687  std::cout << std::endl;
688 
689  }
690 
691 }
692 
693 void
695 
696  std::cout << " Id Gen Name eta phi pT E Vtx1 "
697  << " x y z "
698  << "Moth Vtx2 eta phi R Z Daughters Ecal?" << std::endl;
699 
700  for( int i=0; i<(int)nTracks(); i++ )
701  std::cout << track(i) << std::endl;
702 
703  for( int i=0; i<(int)nVertices(); i++ )
704  std::cout << "i = " << i << " " << vertexType(i) << std::endl;
705 
706 
707 
708 }
709 
710 void
712 
713  nSimTracks = 0;
714  nSimVertices = 0;
715  nGenParticles = 0;
717 
718 }
719 
720 void
722  (*theChargedTracks)[nChargedParticleTracks++] = id;
725  theChargedSize *= 2;
727  }
728 }
729 
730 int
732  if (id>=0 && id<(int)nChargedParticleTracks)
733  return (*theChargedTracks)[id];
734  else
735  return -1;
736 }
737 
738 /*
739 const SimTrack &
740 FBaseSimEvent::embdTrack(int i) const {
741  return (*theSimTracks)[i].simTrack();
742 }
743 
744 const SimVertex &
745 FBaseSimEvent::embdVertex(int i) const {
746  return (*theSimVertices)[i].simVertex();
747 }
748 */
749 
750 const HepMC::GenParticle*
752  return (*theGenParticles)[i];
753 }
754 
755 /*
756 FSimTrack&
757 FBaseSimEvent::track(int id) const {
758  return (*theSimTracks)[id];
759 }
760 
761 
762 FSimVertex&
763 FBaseSimEvent::vertex(int id) const {
764  return (*theSimVertices)[id];
765 }
766 */
const ParticleDataTable * pdt
double lateVertexPosition
type
Definition: HCALResponse.h:21
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
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
HepPDT::ParticleDataTable ParticleDataTable
std::vector< FSimTrack > * theSimTracks
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
KineParticleFilter * myFilter
The particle filter.
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
unsigned int theGenSize
bool notYetToEndVertex(const XYZTLorentzVector &pos) const
Compare the end vertex position with another position.
Definition: FSimTrack.cc:33
FSimVertex & vertex(int id) const
Return vertex with given Id.
unsigned int nVertices() const
Number of vertices.
Definition: FBaseSimEvent.h:89
FSimVertexType & vertexType(int id) const
Return vertex with given Id.
bool propagateToVFcalEntrance(bool first=true)
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
~FBaseSimEvent()
usual virtual destructor
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:18
unsigned int nTracks() const
Number of tracks.
Definition: FBaseSimEvent.h:84
bool propagateToHcalExit(bool first=true)
unsigned int nSimVertices
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
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
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
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
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
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:94
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.
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:83
tuple cout
Definition: gather_cfg.py:145
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 propagateToHOLayer(bool first=true)
bool propagateToPreshowerLayer2(bool first=true)
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