File indexing completed on 2025-01-18 09:17: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
0030 #include "RunAction.hh"
0031
0032 #include "DetectorConstruction.hh"
0033 #include "PrimaryGeneratorAction.hh"
0034
0035 #include "G4AccumulableManager.hh"
0036 #include "G4LogicalVolume.hh"
0037 #include "G4ParticleDefinition.hh"
0038 #include "G4ParticleGun.hh"
0039 #include "G4Run.hh"
0040 #include "G4RunManager.hh"
0041 #include "G4SystemOfUnits.hh"
0042 #include "G4UnitsTable.hh"
0043
0044 namespace B1
0045 {
0046
0047
0048
0049 RunAction::RunAction()
0050 {
0051
0052
0053 const G4double milligray = 1.e-3 * gray;
0054 const G4double microgray = 1.e-6 * gray;
0055 const G4double nanogray = 1.e-9 * gray;
0056 const G4double picogray = 1.e-12 * gray;
0057
0058 new G4UnitDefinition("milligray", "milliGy", "Dose", milligray);
0059 new G4UnitDefinition("microgray", "microGy", "Dose", microgray);
0060 new G4UnitDefinition("nanogray", "nanoGy", "Dose", nanogray);
0061 new G4UnitDefinition("picogray", "picoGy", "Dose", picogray);
0062
0063
0064 G4AccumulableManager* accumulableManager = G4AccumulableManager::Instance();
0065 accumulableManager->RegisterAccumulable(fEdep);
0066 accumulableManager->RegisterAccumulable(fEdep2);
0067 }
0068
0069
0070
0071 void RunAction::BeginOfRunAction(const G4Run*)
0072 {
0073
0074 G4RunManager::GetRunManager()->SetRandomNumberStore(false);
0075
0076
0077 G4AccumulableManager* accumulableManager = G4AccumulableManager::Instance();
0078 accumulableManager->Reset();
0079 }
0080
0081
0082
0083 void RunAction::EndOfRunAction(const G4Run* run)
0084 {
0085 G4int nofEvents = run->GetNumberOfEvent();
0086 if (nofEvents == 0) return;
0087
0088
0089 G4AccumulableManager* accumulableManager = G4AccumulableManager::Instance();
0090 accumulableManager->Merge();
0091
0092
0093
0094 G4double edep = fEdep.GetValue();
0095 G4double edep2 = fEdep2.GetValue();
0096
0097 G4double rms = edep2 - edep * edep / nofEvents;
0098 if (rms > 0.)
0099 rms = std::sqrt(rms);
0100 else
0101 rms = 0.;
0102
0103 const auto detConstruction = static_cast<const DetectorConstruction*>(
0104 G4RunManager::GetRunManager()->GetUserDetectorConstruction());
0105 G4double mass = detConstruction->GetScoringVolume()->GetMass();
0106 G4double dose = edep / mass;
0107 G4double rmsDose = rms / mass;
0108
0109
0110
0111
0112 const auto generatorAction = static_cast<const PrimaryGeneratorAction*>(
0113 G4RunManager::GetRunManager()->GetUserPrimaryGeneratorAction());
0114 G4String runCondition;
0115 if (generatorAction) {
0116 const G4ParticleGun* particleGun = generatorAction->GetParticleGun();
0117 runCondition += particleGun->GetParticleDefinition()->GetParticleName();
0118 runCondition += " of ";
0119 G4double particleEnergy = particleGun->GetParticleEnergy();
0120 runCondition += G4BestUnit(particleEnergy, "Energy");
0121 }
0122
0123
0124
0125 if (IsMaster()) {
0126 G4cout << G4endl << "--------------------End of Global Run-----------------------";
0127 }
0128 else {
0129 G4cout << G4endl << "--------------------End of Local Run------------------------";
0130 }
0131
0132 G4cout << G4endl << " The run consists of " << nofEvents << " " << runCondition << G4endl
0133 << " Cumulated dose per run, in scoring volume : " << G4BestUnit(dose, "Dose")
0134 << " rms = " << G4BestUnit(rmsDose, "Dose") << G4endl
0135 << "------------------------------------------------------------" << G4endl << G4endl;
0136 }
0137
0138
0139
0140 void RunAction::AddEdep(G4double edep)
0141 {
0142 fEdep += edep;
0143 fEdep2 += edep * edep;
0144 }
0145
0146
0147
0148 }