CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
G4SimEvent.cc
Go to the documentation of this file.
3 
4 #include "G4SystemOfUnits.hh"
5 
6 class IdSort{
7 public:
8  bool operator()(const SimTrack& a, const SimTrack& b) {
9  return a.trackId() < b.trackId();
10  }
11 };
12 
13 
14 G4SimEvent::G4SimEvent() : hepMCEvent(0),
15  weight_(0),
16  collisionPoint_(math::XYZTLorentzVectorD(0.,0.,0.,0.)),
17  nparam_(0),param_(0) {}
18 
20 {
21 
22 /*
23  while ( !g4tracks.empty() )
24  {
25  delete g4tracks.back() ;
26  g4tracks.pop_back() ;
27  }
28  while ( !g4vertices.empty() )
29  {
30  delete g4vertices.back() ;
31  g4vertices.pop_back() ;
32  }
33 */
34 
35  // per suggestion by Chris Jones, it's faster
36  // that delete back() and pop_back()
37  //
38  unsigned int i = 0 ;
39 
40  for ( i=0; i<g4tracks.size(); i++ )
41  {
42  delete g4tracks[i] ;
43  g4tracks[i] = 0 ;
44  }
45  g4tracks.clear() ;
46 
47  for ( i=0; i<g4vertices.size(); i++ )
48  {
49  delete g4vertices[i] ;
50  g4vertices[i] = 0 ;
51  }
52  g4vertices.clear();
53 }
54 
56 {
57  for (unsigned int i=0; i<g4tracks.size(); i++)
58  {
59  G4SimTrack * trk = g4tracks[i];
60  int ip = trk->part();
61  math::XYZTLorentzVectorD p( trk->momentum().x()/GeV,
62  trk->momentum().y()/GeV,
63  trk->momentum().z()/GeV,
64  trk->energy()/GeV ) ;
65  int iv = trk->ivert();
66  int ig = trk->igenpart();
67  int id = trk->id();
68  math::XYZVectorD tkpos( trk->trackerSurfacePosition().x()/cm,
69  trk->trackerSurfacePosition().y()/cm,
70  trk->trackerSurfacePosition().z()/cm ) ;
72  trk->trackerSurfaceMomentum().y()/GeV,
73  trk->trackerSurfaceMomentum().z()/GeV,
74  trk->trackerSurfaceMomentum().e()/GeV ) ;
75  // ip = particle ID as PDG
76  // pp = 4-momentum
77  // iv = corresponding G4SimVertex index
78  // ig = corresponding GenParticle index
79  SimTrack t = SimTrack(ip,p,iv,ig,tkpos,tkmom);
80  t.setTrackId(id);
81  t.setEventId(EncodedEventId(0));
82  c.push_back(t);
83  }
84  std::stable_sort(c.begin(),c.end(),IdSort());
85 
86 }
87 
89 {
90  for (unsigned int i=0; i<g4vertices.size(); i++)
91  {
92  G4SimVertex * vtx = g4vertices[i];
93  //
94  // starting 1_1_0_pre3, SimVertex stores in cm !!!
95  //
96  math::XYZVectorD v3( vtx->vertexPosition().x()/cm,
97  vtx->vertexPosition().y()/cm,
98  vtx->vertexPosition().z()/cm ) ;
99  float t = vtx->vertexGlobalTime()/second;
100  int iv = vtx->parentIndex();
101  // vv = position
102  // t = global time
103  // iv = index of the parent in the SimEvent SimTrack container (-1 if no parent)
104  SimVertex v = SimVertex(v3,t,iv,i);
106  c.push_back(v);
107  }
108 }
109 
const double energy() const
Definition: G4SimTrack.h:35
int i
Definition: DBlmapReader.cc:9
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
const math::XYZVectorD & trackerSurfacePosition() const
Definition: G4SimTrack.h:41
bool operator()(const SimTrack &a, const SimTrack &b)
Definition: G4SimEvent.cc:8
const int part() const
Definition: G4SimTrack.h:33
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:45
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
U second(std::pair< T, U > const &p)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:9
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:46
const math::XYZVectorD & vertexPosition() const
index of the parent (-1 if no parent)
Definition: G4SimVertex.h:16
void load(edm::SimTrackContainer &c) const
Definition: G4SimEvent.cc:55
struct @254 param_
unsigned int trackId() const
Definition: CoreSimTrack.h:49
std::vector< SimVertex > SimVertexContainer
int const ivert() const
Definition: G4SimTrack.h:36
double b
Definition: hdecay.h:120
int const igenpart() const
Definition: G4SimTrack.h:37
const int parentIndex() const
Definition: G4SimVertex.h:18
virtual ~G4SimEvent()
Definition: G4SimEvent.cc:19
double a
Definition: hdecay.h:121
const double vertexGlobalTime() const
Definition: G4SimVertex.h:17
const int id() const
Definition: G4SimTrack.h:32
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: G4SimTrack.h:42
std::vector< SimTrack > SimTrackContainer
mathSSE::Vec4< T > v
void setEventId(EncodedEventId e)
Definition: CoreSimVertex.h:28
const math::XYZVectorD & momentum() const
Definition: G4SimTrack.h:34