CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
G4SimEvent Class Reference

#include <G4SimEvent.h>

Public Member Functions

void add (G4SimTrack *t)
 
void add (G4SimVertex *v)
 
void collisionPoint (const math::XYZTLorentzVectorD &v)
 
const math::XYZTLorentzVectorDcollisionPoint () const
 
 G4SimEvent ()
 
const G4SimTrackg4track (int i) const
 
const G4SimVertexg4vertex (int i) const
 
void hepEvent (const HepMC::GenEvent *r)
 
const HepMC::GenEvent * hepEvent () const
 
void load (edm::SimTrackContainer &c) const
 
void load (edm::SimVertexContainer &c) const
 
unsigned int nGenParts () const
 
void nparam (int n)
 
const int nparam () const
 
unsigned int nTracks () const
 
unsigned int nVertices () const
 
void param (const std::vector< float > &p)
 
const std::vector< float > & param () const
 
void weight (float w)
 
float weight () const
 
virtual ~G4SimEvent ()
 

Protected Attributes

math::XYZTLorentzVectorD collisionPoint_
 
std::vector< G4SimTrack * > g4tracks
 
std::vector< G4SimVertex * > g4vertices
 
const HepMC::GenEvent * hepMCEvent
 
int nparam_
 
std::vector< float > param_
 
float weight_
 

Detailed Description

Definition at line 14 of file G4SimEvent.h.

Constructor & Destructor Documentation

G4SimEvent::G4SimEvent ( )

Definition at line 13 of file G4SimEvent.cc.

13  : hepMCEvent(nullptr),
14  weight_(0),
16  nparam_(0),param_(0) {}
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
int nparam_
Definition: G4SimEvent.h:42
float weight_
Definition: G4SimEvent.h:40
const HepMC::GenEvent * hepMCEvent
Definition: G4SimEvent.h:39
std::vector< float > param_
Definition: G4SimEvent.h:43
math::XYZTLorentzVectorD collisionPoint_
Definition: G4SimEvent.h:41
G4SimEvent::~G4SimEvent ( )
virtual

Definition at line 18 of file G4SimEvent.cc.

References g4tracks, g4vertices, mps_fire::i, and nullptr.

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 }
#define nullptr
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:45

Member Function Documentation

void G4SimEvent::add ( G4SimTrack t)
inline

Definition at line 34 of file G4SimEvent.h.

References g4tracks.

Referenced by SimTrackManager::getOrCreateVertex(), SimTrackManager::reallyStoreTracks(), and counter.Counter::register().

34 { g4tracks.push_back(t); }
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
void G4SimEvent::add ( G4SimVertex v)
inline

Definition at line 35 of file G4SimEvent.h.

References g4vertices.

Referenced by counter.Counter::register().

35 { g4vertices.push_back(v); }
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:45
void G4SimEvent::collisionPoint ( const math::XYZTLorentzVectorD v)
inline

Definition at line 28 of file G4SimEvent.h.

References collisionPoint_, and findQualityFiles::v.

Referenced by RunManager::produce().

28 { collisionPoint_ = v ; }
math::XYZTLorentzVectorD collisionPoint_
Definition: G4SimEvent.h:41
const math::XYZTLorentzVectorD& G4SimEvent::collisionPoint ( ) const
inline

Definition at line 29 of file G4SimEvent.h.

References collisionPoint_.

29 { return collisionPoint_; }
math::XYZTLorentzVectorD collisionPoint_
Definition: G4SimEvent.h:41
const G4SimTrack& G4SimEvent::g4track ( int  i) const
inline

Definition at line 36 of file G4SimEvent.h.

References g4tracks.

36 { return *g4tracks[i-1]; }
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
const G4SimVertex& G4SimEvent::g4vertex ( int  i) const
inline

Definition at line 37 of file G4SimEvent.h.

References g4vertices.

37 { return *g4vertices[i-1]; }
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:45
void G4SimEvent::hepEvent ( const HepMC::GenEvent *  r)
inline

Definition at line 24 of file G4SimEvent.h.

References hepMCEvent, and alignCSCRings::r.

Referenced by RunManager::produce().

24 { hepMCEvent = r; }
const HepMC::GenEvent * hepMCEvent
Definition: G4SimEvent.h:39
const HepMC::GenEvent* G4SimEvent::hepEvent ( ) const
inline

Definition at line 25 of file G4SimEvent.h.

References hepMCEvent.

25 { return hepMCEvent; }
const HepMC::GenEvent * hepMCEvent
Definition: G4SimEvent.h:39
void G4SimEvent::load ( edm::SimTrackContainer c) const

Definition at line 40 of file G4SimEvent.cc.

References G4SimTrack::energy(), g4tracks, GeV, mps_fire::i, G4SimTrack::id(), G4SimTrack::igenpart(), G4SimTrack::ivert(), G4SimTrack::momentum(), AlCaHLTBitMon_ParallelJobs::p, G4SimTrack::part(), mcMuonSeeds_cfi::SimTrack, lumiQTWidget::t, G4SimTrack::trackerSurfaceMomentum(), and G4SimTrack::trackerSurfacePosition().

Referenced by OscarProducer::produce(), and OscarMTProducer::produce().

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 }
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 part() const
Definition: G4SimTrack.h:33
int ivert() const
Definition: G4SimTrack.h:36
const math::XYZVectorD & trackerSurfacePosition() const
Definition: G4SimTrack.h:41
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
int igenpart() const
Definition: G4SimTrack.h:37
int id() const
Definition: G4SimTrack.h:32
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: G4SimTrack.h:42
const math::XYZVectorD & momentum() const
Definition: G4SimTrack.h:34
void G4SimEvent::load ( edm::SimVertexContainer c) const

Definition at line 73 of file G4SimEvent.cc.

References g4vertices, mps_fire::i, G4SimVertex::parentIndex(), G4SimVertex::processType(), edm::second(), CoreSimVertex::setEventId(), SimVertex::setProcessType(), mcMuonSeeds_cfi::SimVertex, lumiQTWidget::t, findQualityFiles::v, G4SimVertex::vertexGlobalTime(), G4SimVertex::vertexPosition(), and badGlobalMuonTaggersAOD_cff::vtx.

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 }
int parentIndex() const
Definition: G4SimVertex.h:18
void setProcessType(unsigned int ty)
Definition: SimVertex.h:39
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
double vertexGlobalTime() const
Definition: G4SimVertex.h:17
void setEventId(EncodedEventId e)
Definition: CoreSimVertex.h:28
unsigned int G4SimEvent::nGenParts ( ) const
inline

Definition at line 23 of file G4SimEvent.h.

References hepMCEvent.

Referenced by RunManager::produce().

23 { return hepMCEvent->particles_size(); }
const HepMC::GenEvent * hepMCEvent
Definition: G4SimEvent.h:39
void G4SimEvent::nparam ( int  n)
inline

Definition at line 30 of file G4SimEvent.h.

References gen::n, and nparam_.

30 { nparam_ = n; }
int nparam_
Definition: G4SimEvent.h:42
const int G4SimEvent::nparam ( ) const
inline

Definition at line 31 of file G4SimEvent.h.

References nparam_.

31 { return nparam_; }
int nparam_
Definition: G4SimEvent.h:42
unsigned int G4SimEvent::nTracks ( ) const
inline

Definition at line 21 of file G4SimEvent.h.

References g4tracks.

Referenced by RunManager::produce().

21 { return g4tracks.size(); }
std::vector< G4SimTrack * > g4tracks
Definition: G4SimEvent.h:44
unsigned int G4SimEvent::nVertices ( ) const
inline

Definition at line 22 of file G4SimEvent.h.

References g4vertices.

Referenced by RunManager::produce().

22 { return g4vertices.size(); }
std::vector< G4SimVertex * > g4vertices
Definition: G4SimEvent.h:45
void G4SimEvent::param ( const std::vector< float > &  p)
inline

Definition at line 32 of file G4SimEvent.h.

References AlCaHLTBitMon_ParallelJobs::p, and param_.

const std::vector<float>& G4SimEvent::param ( ) const
inline

Definition at line 33 of file G4SimEvent.h.

References param_.

33 { return param_; }
std::vector< float > param_
Definition: G4SimEvent.h:43
void G4SimEvent::weight ( float  w)
inline

Definition at line 26 of file G4SimEvent.h.

References w, and weight_.

Referenced by RunManager::produce().

26 { weight_ = w; }
const double w
Definition: UKUtility.cc:23
float weight_
Definition: G4SimEvent.h:40
float G4SimEvent::weight ( ) const
inline

Definition at line 27 of file G4SimEvent.h.

References weight_.

27 { return weight_; }
float weight_
Definition: G4SimEvent.h:40

Member Data Documentation

math::XYZTLorentzVectorD G4SimEvent::collisionPoint_
protected

Definition at line 41 of file G4SimEvent.h.

Referenced by collisionPoint().

std::vector<G4SimTrack *> G4SimEvent::g4tracks
protected

Definition at line 44 of file G4SimEvent.h.

Referenced by add(), g4track(), load(), nTracks(), and ~G4SimEvent().

std::vector<G4SimVertex *> G4SimEvent::g4vertices
protected

Definition at line 45 of file G4SimEvent.h.

Referenced by add(), g4vertex(), load(), nVertices(), and ~G4SimEvent().

const HepMC::GenEvent* G4SimEvent::hepMCEvent
protected

Definition at line 39 of file G4SimEvent.h.

Referenced by hepEvent(), and nGenParts().

int G4SimEvent::nparam_
protected

Definition at line 42 of file G4SimEvent.h.

Referenced by nparam().

std::vector<float> G4SimEvent::param_
protected

Definition at line 43 of file G4SimEvent.h.

Referenced by param().

float G4SimEvent::weight_
protected

Definition at line 40 of file G4SimEvent.h.

Referenced by weight().