File indexing completed on 2026-03-30 07:51:13
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "RE01EventAction.hh"
0030
0031 #include "RE01CalorimeterHit.hh"
0032 #include "RE01TrackerHit.hh"
0033 #include "RE01Trajectory.hh"
0034
0035 #include "G4Event.hh"
0036 #include "G4EventManager.hh"
0037 #include "G4HCofThisEvent.hh"
0038 #include "G4PrimaryParticle.hh"
0039 #include "G4PrimaryVertex.hh"
0040 #include "G4SDManager.hh"
0041 #include "G4SystemOfUnits.hh"
0042 #include "G4TrajectoryContainer.hh"
0043 #include "G4UImanager.hh"
0044 #include "G4VHitsCollection.hh"
0045 #include "G4VVisManager.hh"
0046 #include "G4ios.hh"
0047
0048
0049 RE01EventAction::RE01EventAction() : G4UserEventAction(), fTrackerCollID(-1), fCalorimeterCollID(-1)
0050 {
0051 ;
0052 }
0053
0054
0055 RE01EventAction::~RE01EventAction()
0056 {
0057 ;
0058 }
0059
0060
0061 void RE01EventAction::BeginOfEventAction(const G4Event*)
0062 {
0063 G4SDManager* SDman = G4SDManager::GetSDMpointer();
0064 if (fTrackerCollID < 0 || fCalorimeterCollID < 0) {
0065 G4String colNam;
0066 fTrackerCollID = SDman->GetCollectionID(colNam = "trackerCollection");
0067 fCalorimeterCollID = SDman->GetCollectionID(colNam = "calCollection");
0068 }
0069 }
0070
0071
0072 void RE01EventAction::EndOfEventAction(const G4Event* evt)
0073 {
0074 G4cout << ">>> Summary of Event " << evt->GetEventID() << G4endl;
0075 if (evt->GetNumberOfPrimaryVertex() == 0) {
0076 G4cout << "Event is empty." << G4endl;
0077 return;
0078 }
0079
0080 if (fTrackerCollID < 0 || fCalorimeterCollID < 0) return;
0081
0082 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
0083 RE01TrackerHitsCollection* THC = 0;
0084 RE01CalorimeterHitsCollection* CHC = 0;
0085 if (HCE) {
0086 THC = (RE01TrackerHitsCollection*)(HCE->GetHC(fTrackerCollID));
0087 CHC = (RE01CalorimeterHitsCollection*)(HCE->GetHC(fCalorimeterCollID));
0088 }
0089
0090 if (THC) {
0091 int n_hit = THC->entries();
0092 G4cout << G4endl;
0093 G4cout << "Tracker hits "
0094 << "--------------------------------------------------------------" << G4endl;
0095 G4cout << n_hit << " hits are stored in RE01TrackerHitsCollection." << G4endl;
0096 if (fpEventManager->GetVerboseLevel() > 0) {
0097 G4cout << "List of hits in tracker" << G4endl;
0098 for (int i = 0; i < n_hit; i++) {
0099 (*THC)[i]->Print();
0100 }
0101 }
0102 }
0103 if (CHC) {
0104 int n_hit = CHC->entries();
0105 G4cout << G4endl;
0106 G4cout << "Calorimeter hits "
0107 << "--------------------------------------------------------------" << G4endl;
0108 G4cout << n_hit << " hits are stored in RE01CalorimeterHitsCollection." << G4endl;
0109 G4double totE = 0;
0110 for (int i = 0; i < n_hit; i++) {
0111 totE += (*CHC)[i]->GetEdep();
0112 }
0113 G4cout << " Total energy deposition in calorimeter : " << totE / GeV << " (GeV)" << G4endl;
0114 }
0115
0116
0117 G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
0118 G4int n_trajectories = 0;
0119 if (trajectoryContainer) n_trajectories = trajectoryContainer->entries();
0120
0121 G4cout << G4endl;
0122 G4cout << "Trajectories in tracker "
0123 << "--------------------------------------------------------------" << G4endl;
0124 if (fpEventManager->GetVerboseLevel() > 0) {
0125 for (G4int i = 0; i < n_trajectories; i++) {
0126 RE01Trajectory* trj = (RE01Trajectory*)((*(evt->GetTrajectoryContainer()))[i]);
0127 trj->ShowTrajectory();
0128 }
0129 }
0130
0131 G4cout << G4endl;
0132 G4cout << "Primary particles "
0133 << "--------------------------------------------------------------" << G4endl;
0134 G4int n_vertex = evt->GetNumberOfPrimaryVertex();
0135 for (G4int iv = 0; iv < n_vertex; iv++) {
0136 G4PrimaryVertex* pv = evt->GetPrimaryVertex(iv);
0137 G4cout << G4endl;
0138 G4cout << "Primary vertex " << G4ThreeVector(pv->GetX0(), pv->GetY0(), pv->GetZ0())
0139 << " at t = " << (pv->GetT0()) / ns << " [ns]" << G4endl;
0140 if (fpEventManager->GetVerboseLevel() > 0) {
0141 G4PrimaryParticle* pp = pv->GetPrimary();
0142 while (pp) {
0143 PrintPrimary(pp, 0);
0144 pp = pp->GetNext();
0145 }
0146 }
0147 }
0148 }
0149
0150
0151 void RE01EventAction::PrintPrimary(G4PrimaryParticle* pp, G4int ind)
0152 {
0153 for (G4int ii = 0; ii <= ind; ii++) {
0154 G4cout << " ";
0155 }
0156 G4cout << "==PDGcode " << pp->GetPDGcode() << " ";
0157 if (pp->GetG4code() != 0) {
0158 G4cout << "(" << pp->GetG4code()->GetParticleName() << ")";
0159 }
0160 else {
0161 G4cout << "is not defined in G4";
0162 }
0163 G4cout << " " << pp->GetMomentum() / GeV << " [GeV] ";
0164 if (pp->GetTrackID() < 0) {
0165 G4cout << G4endl;
0166 }
0167 else {
0168 G4cout << ">>> G4Track ID " << pp->GetTrackID() << G4endl;
0169 }
0170
0171 G4PrimaryParticle* daughter = pp->GetDaughter();
0172 while (daughter) {
0173 PrintPrimary(daughter, ind + 1);
0174 daughter = daughter->GetNext();
0175 }
0176 }