Geant4  10.02.p03
RE01EventAction Class Reference

#include <RE01EventAction.hh>

Inheritance diagram for RE01EventAction:
Collaboration diagram for RE01EventAction:

Public Member Functions

 RE01EventAction ()
 
virtual ~RE01EventAction ()
 
virtual void BeginOfEventAction (const G4Event *)
 
virtual void EndOfEventAction (const G4Event *)
 
- Public Member Functions inherited from G4UserEventAction
 G4UserEventAction ()
 
virtual ~G4UserEventAction ()
 
virtual void SetEventManager (G4EventManager *value)
 

Private Member Functions

void PrintPrimary (G4PrimaryParticle *pp, G4int ind)
 

Private Attributes

G4int fTrackerCollID
 
G4int fCalorimeterCollID
 

Additional Inherited Members

- Protected Attributes inherited from G4UserEventAction
G4EventManagerfpEventManager
 

Detailed Description

Definition at line 40 of file RE01EventAction.hh.

Constructor & Destructor Documentation

◆ RE01EventAction()

RE01EventAction::RE01EventAction ( )

Definition at line 51 of file RE01EventAction.cc.

◆ ~RE01EventAction()

RE01EventAction::~RE01EventAction ( )
virtual

Definition at line 57 of file RE01EventAction.cc.

58 {;}

Member Function Documentation

◆ BeginOfEventAction()

void RE01EventAction::BeginOfEventAction ( const G4Event )
virtual

Reimplemented from G4UserEventAction.

Definition at line 61 of file RE01EventAction.cc.

62 {
65  {
66  G4String colNam;
67  fTrackerCollID = SDman->GetCollectionID(colNam="trackerCollection");
68  fCalorimeterCollID = SDman->GetCollectionID(colNam="calCollection");
69  }
70 }
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
Here is the call graph for this function:

◆ EndOfEventAction()

void RE01EventAction::EndOfEventAction ( const G4Event evt)
virtual

Reimplemented from G4UserEventAction.

Definition at line 73 of file RE01EventAction.cc.

74 {
75  G4cout << ">>> Summary of Event " << evt->GetEventID() << G4endl;
76 
77  if(fTrackerCollID<0||fCalorimeterCollID<0) return;
78 
79  G4HCofThisEvent * HCE = evt->GetHCofThisEvent();
82  if(HCE)
83  {
86  }
87 
88  if(THC)
89  {
90  int n_hit = THC->entries();
91  G4cout << G4endl;
92  G4cout << "Tracker hits " <<
93  "--------------------------------------------------------------"
94  << G4endl;
95  G4cout << n_hit << " hits are stored in RE01TrackerHitsCollection."
96  << G4endl;
98  {
99  G4cout << "List of hits in tracker" << G4endl;
100  for(int i=0;i<n_hit;i++)
101  { (*THC)[i]->Print(); }
102  }
103  }
104  if(CHC)
105  {
106  int n_hit = CHC->entries();
107  G4cout << G4endl;
108  G4cout << "Calorimeter hits "<<
109  "--------------------------------------------------------------"
110  << G4endl;
111  G4cout << n_hit << " hits are stored in RE01CalorimeterHitsCollection."
112  << G4endl;
113  G4double totE = 0;
114  for(int i=0;i<n_hit;i++)
115  { totE += (*CHC)[i]->GetEdep(); }
116  G4cout << " Total energy deposition in calorimeter : "
117  << totE / GeV << " (GeV)" << G4endl;
118  }
119 
120  // get number of stored trajectories
121  G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
122  G4int n_trajectories = 0;
123  if (trajectoryContainer) n_trajectories = trajectoryContainer->entries();
124  // extract the trajectories and print them out
125  G4cout << G4endl;
126  G4cout << "Trajectories in tracker "<<
127  "--------------------------------------------------------------"
128  << G4endl;
130  {
131  for(G4int i=0; i<n_trajectories; i++)
132  {
133  RE01Trajectory* trj =
134  (RE01Trajectory*)((*(evt->GetTrajectoryContainer()))[i]);
135  trj->ShowTrajectory();
136  }
137  }
138 
139  G4cout << G4endl;
140  G4cout << "Primary particles "<<
141  "--------------------------------------------------------------"
142  << G4endl;
143  G4int n_vertex = evt->GetNumberOfPrimaryVertex();
144  for(G4int iv=0;iv<n_vertex;iv++)
145  {
146  G4PrimaryVertex* pv = evt->GetPrimaryVertex(iv);
147  G4cout << G4endl;
148  G4cout << "Primary vertex "
149  << G4ThreeVector(pv->GetX0(),pv->GetY0(),pv->GetZ0())
150  << " at t = " << (pv->GetT0())/ns << " [ns]" << G4endl;
152  {
154  while(pp)
155  {
156  PrintPrimary(pp,0);
157  pp = pp->GetNext();
158  }
159  }
160  }
161 }
G4double GetT0() const
void PrintPrimary(G4PrimaryParticle *pp, G4int ind)
CLHEP::Hep3Vector G4ThreeVector
G4PrimaryParticle * GetPrimary(G4int i=0) const
G4double GetZ0() const
G4EventManager * fpEventManager
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfPrimaryVertex() const
Definition: G4Event.hh:164
G4GLOB_DLL std::ostream G4cout
G4PrimaryParticle * GetNext() const
G4double GetX0() const
static const double GeV
Definition: G4SIunits.hh:214
G4TrajectoryContainer * GetTrajectoryContainer() const
Definition: G4Event.hh:189
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:167
G4double GetY0() const
G4int GetEventID() const
Definition: G4Event.hh:151
G4int GetVerboseLevel()
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define ns
Definition: xmlparse.cc:614
virtual void ShowTrajectory(std::ostream &os=G4cout) const
Here is the call graph for this function:

◆ PrintPrimary()

void RE01EventAction::PrintPrimary ( G4PrimaryParticle pp,
G4int  ind 
)
private

Definition at line 164 of file RE01EventAction.cc.

165 {
166  for(G4int ii=0;ii<=ind;ii++)
167  { G4cout << " "; }
168  G4cout << "==PDGcode " << pp->GetPDGcode() << " ";
169  if(pp->GetG4code()!=0)
170  { G4cout << "(" << pp->GetG4code()->GetParticleName() << ")"; }
171  else
172  { G4cout << "is not defined in G4"; }
173  G4cout << " " << pp->GetMomentum()/GeV << " [GeV] ";
174  if(pp->GetTrackID()<0)
175  { G4cout << G4endl; }
176  else
177  { G4cout << ">>> G4Track ID " << pp->GetTrackID() << G4endl; }
178 
179  G4PrimaryParticle* daughter = pp->GetDaughter();
180  while(daughter)
181  {
182  PrintPrimary(daughter,ind+1);
183  daughter = daughter->GetNext();
184  }
185 }
G4ParticleDefinition * GetG4code() const
void PrintPrimary(G4PrimaryParticle *pp, G4int ind)
G4int GetPDGcode() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4PrimaryParticle * GetNext() const
G4PrimaryParticle * GetDaughter() const
static const double GeV
Definition: G4SIunits.hh:214
G4int GetTrackID() const
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fCalorimeterCollID

G4int RE01EventAction::fCalorimeterCollID
private

Definition at line 52 of file RE01EventAction.hh.

◆ fTrackerCollID

G4int RE01EventAction::fTrackerCollID
private

Definition at line 51 of file RE01EventAction.hh.


The documentation for this class was generated from the following files: