CMS 3D CMS Logo

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 
14  weight_(0),
15  collisionPoint_(math::XYZTLorentzVectorD(0.,0.,0.,0.)),
16  nparam_(0),param_(0) {}
17 
19 {
20  // per suggestion by Chris Jones, it's faster
21  // that delete back() and pop_back()
22  //
23  unsigned int i = 0 ;
24 
25  for ( i=0; i<g4tracks.size(); i++ )
26  {
27  delete g4tracks[i] ;
28  g4tracks[i] = nullptr ;
29  }
30  g4tracks.clear() ;
31 
32  for ( i=0; i<g4vertices.size(); i++ )
33  {
34  delete g4vertices[i] ;
35  g4vertices[i] = nullptr ;
36  }
37  g4vertices.clear();
38 }
39 
41 {
42  for (unsigned int i=0; i<g4tracks.size(); i++)
43  {
44  G4SimTrack * trk = g4tracks[i];
45  int ip = trk->part();
47  trk->momentum().y()/GeV,
48  trk->momentum().z()/GeV,
49  trk->energy()/GeV ) ;
50  int iv = trk->ivert();
51  int ig = trk->igenpart();
52  int id = trk->id();
53  math::XYZVectorD tkpos( trk->trackerSurfacePosition().x()/cm,
54  trk->trackerSurfacePosition().y()/cm,
55  trk->trackerSurfacePosition().z()/cm ) ;
57  trk->trackerSurfaceMomentum().y()/GeV,
58  trk->trackerSurfaceMomentum().z()/GeV,
59  trk->trackerSurfaceMomentum().e()/GeV ) ;
60  // ip = particle ID as PDG
61  // pp = 4-momentum
62  // iv = corresponding G4SimVertex index
63  // ig = corresponding GenParticle index
64  SimTrack t = SimTrack(ip,p,iv,ig,tkpos,tkmom);
65  t.setTrackId(id);
66  t.setEventId(EncodedEventId(0));
67  c.push_back(t);
68  }
69  std::stable_sort(c.begin(),c.end(),IdSort());
70 
71 }
72 
74 {
75  for (unsigned int i=0; i<g4vertices.size(); i++)
76  {
78  //
79  // starting 1_1_0_pre3, SimVertex stores in cm !!!
80  //
81  math::XYZVectorD v3( vtx->vertexPosition().x()/cm,
82  vtx->vertexPosition().y()/cm,
83  vtx->vertexPosition().z()/cm ) ;
84  float t = vtx->vertexGlobalTime()/second;
85  int iv = vtx->parentIndex();
86  // vv = position
87  // t = global time
88  // iv = index of the parent in the SimEvent SimTrack container (-1 if no parent)
89  SimVertex v = SimVertex(v3,t,iv,i);
90  v.setProcessType(vtx->processType());
92  c.push_back(v);
93  }
94 }
95 
double energy() const
Definition: G4SimTrack.h:35
const double GeV
Definition: MathUtil.h:16
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
int parentIndex() const
Definition: G4SimVertex.h:18
int part() const
Definition: G4SimTrack.h:33
int ivert() const
Definition: G4SimTrack.h:36
const math::XYZVectorD & trackerSurfacePosition() const
Definition: G4SimTrack.h:41
void setProcessType(unsigned int ty)
Definition: SimVertex.h:39
bool operator()(const SimTrack &a, const SimTrack &b)
Definition: G4SimEvent.cc:8
#define nullptr
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
U second(std::pair< T, U > const &p)
unsigned int processType() const
Definition: G4SimVertex.h:19
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:45
const math::XYZVectorD & vertexPosition() const
index of the parent (-1 if no parent)
Definition: G4SimVertex.h:16
Definition: Error.h:16
void load(edm::SimTrackContainer &c) const
Definition: G4SimEvent.cc:40
struct @656 param_
unsigned int trackId() const
Definition: CoreSimTrack.h:34
std::vector< SimVertex > SimVertexContainer
double b
Definition: hdecay.h:120
int igenpart() const
Definition: G4SimTrack.h:37
virtual ~G4SimEvent()
Definition: G4SimEvent.cc:18
int id() const
Definition: G4SimTrack.h:32
double a
Definition: hdecay.h:121
double vertexGlobalTime() const
Definition: G4SimVertex.h:17
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: G4SimTrack.h:42
std::vector< SimTrack > SimTrackContainer
void setEventId(EncodedEventId e)
Definition: CoreSimVertex.h:28
const math::XYZVectorD & momentum() const
Definition: G4SimTrack.h:34