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