Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 /// \file Run.cc 27 /// \brief Implementation of the Run class 28 // 29 // 30 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 33 34 #include "Run.hh" 35 36 #include "G4Run.hh" 37 #include "G4RunManager.hh" 38 #include "G4SystemOfUnits.hh" 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 Run::Run() 43 : G4Run(), 44 fNumEvents(0), 45 fPrimaryParticleId(0), 46 fPrimaryParticleEnergy(0.0), 47 fPrimaryParticleDirection(G4ThreeVector(0.0, 0.0, 0.0)), 48 fTargetMaterialName(""), 49 fCubicVolumeScoringUpDown(1.0), 50 fCubicVolumeScoringSide(1.0) 51 { 52 fSteppingArray.fill(0.0); 53 fTrackingArray1.fill(0); 54 fTrackingArray2.fill(0.0); 55 } 56 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 58 59 void Run::RecordEvent(const G4Event* anEvent) 60 { 61 // This method is called automatically by the Geant4 kernel (not by the user!) at the end 62 // of each event : in MT-mode, it is called only for the working thread that handled the event. 63 G4int nEvt = anEvent->GetEventID(); 64 if (nEvt % 10 == 0) G4cout << " Event#=" << nEvt << G4endl; 65 G4Run::RecordEvent(anEvent); 66 } 67 68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 69 70 void Run::Merge(const G4Run* aRun) 71 { 72 // This method is called automatically by the Geant4 kernel (not by the user!) only in the case 73 // of multithreaded mode and only for working threads. 74 const Run* localRun = static_cast<const Run*>(aRun); 75 fPrimaryParticleId = localRun->GetPrimaryParticleId(); 76 fPrimaryParticleEnergy = localRun->GetPrimaryParticleEnergy(); 77 fPrimaryParticleDirection = localRun->GetPrimaryParticleDirection(); 78 fTargetMaterialName = localRun->GetTargetMaterialName(); 79 fCubicVolumeScoringUpDown = localRun->GetCubicVolumeScoringUpDown(); 80 fCubicVolumeScoringSide = localRun->GetCubicVolumeScoringSide(); 81 fNumEvents += localRun->GetNumberOfEvent(); 82 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) { 83 fSteppingArray[i] += localRun->GetSteppingArray()[i]; 84 } 85 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 86 fTrackingArray1[i] += localRun->GetTrackingArray1()[i]; 87 fTrackingArray2[i] += localRun->GetTrackingArray2()[i]; 88 } 89 G4Run::Merge(aRun); 90 } 91 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 93 94 void Run::PrintInfo() const 95 { 96 // This method is called by RunAction::EndOfRunAction. In MT-mode, only the master thread 97 // calls it. 98 const G4double floatingNumberOfEvents = 99 std::max(1.0, fNumEvents > 0 ? fNumEvents * 1.0 : GetNumberOfEvent() * 1.0); 100 // The fluence in the scoring volume is defined as sum of step lengths in that volume 101 // divided by the volume of that scoring volume. 102 const G4double conversionFactor = CLHEP::cm * CLHEP::cm; // From mm^-2 to cm^-2 103 const G4double factorUpDown = 104 conversionFactor / (fCubicVolumeScoringUpDown * floatingNumberOfEvents); 105 const G4double factorSide = conversionFactor / (fCubicVolumeScoringSide * floatingNumberOfEvents); 106 G4cout << std::setprecision(6) << G4endl << G4endl 107 << " =============== Run::PrintInfo() =============== \t RunID = " << GetRunID() 108 << G4endl << " Primary particle PDG code = " << fPrimaryParticleId << G4endl 109 << " Primary particle kinetic energy = " << fPrimaryParticleEnergy / CLHEP::GeV << " GeV" 110 << G4endl << " Primary particle direction = " << fPrimaryParticleDirection << G4endl 111 << " Target material = " << fTargetMaterialName << G4endl 112 << " Cubic-volume scoring up-down = " << fCubicVolumeScoringUpDown << " mm^3" << G4endl 113 << " Cubic-volume scoring side = " << fCubicVolumeScoringSide << " mm^3" << G4endl 114 << " Number of events = " << floatingNumberOfEvents << G4endl 115 << " Conversion factor: fluence from mm^-2 to cm^-2 = " << conversionFactor << G4endl 116 << " Particle fluence in unit of cm^-2 :" << G4endl; 117 for (G4int i = 0; i < SteppingAction::fkNumberScoringVolumes; ++i) { 118 G4double factor = (i == 1 ? factorSide : factorUpDown); 119 for (G4int j = 0; j < SteppingAction::fkNumberKinematicRegions; ++j) { 120 for (G4int k = 0; k < SteppingAction::fkNumberParticleTypes; ++k) { 121 G4int index = SteppingAction::GetIndex(i, j, k); 122 // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ") ->" << index; 123 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12) 124 << SteppingAction::fkArrayScoringVolumeNames[i] << " " << std::setw(12) 125 << SteppingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12) 126 << SteppingAction::fkArrayParticleTypeNames[k] << " " << std::setw(8) 127 << factor * fSteppingArray[index] << G4endl; 128 } 129 } 130 } 131 G4cout << " ------------------------------------------------------------- " << G4endl 132 << " Extra information: particle production \t \t <N> <E_kin> <Sum_Ekin> [MeV]" 133 << G4endl; 134 const G4double normalization = 1.0 / floatingNumberOfEvents; 135 for (G4int i = 0; i < TrackingAction::fkNumberScoringVolumes; ++i) { 136 for (G4int j = 0; j < TrackingAction::fkNumberKinematicRegions; ++j) { 137 for (G4int k = 0; k < TrackingAction::fkNumberParticleTypes; ++k) { 138 G4int index = TrackingAction::GetIndex(i, j, k); 139 // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ") ->" << index; 140 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12) 141 << TrackingAction::fkArrayScoringVolumeNames[i] << " " << std::setw(12) 142 << TrackingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12) 143 << TrackingAction::fkArrayParticleTypeNames[k] << " " << std::setw(8) 144 << normalization * fTrackingArray1[index] << " " << std::setw(8) 145 << (fTrackingArray1[index] > 0 ? fTrackingArray2[index] / fTrackingArray1[index] 146 : 0.0) 147 << " " << std::setw(8) << normalization * fTrackingArray2[index] << G4endl; 148 } 149 } 150 } 151 G4cout << " ============================================================= " << G4endl << G4endl; 152 } 153 154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 155 156 void Run::SetSteppingArray( 157 const std::array<G4double, SteppingAction::fkNumberCombinations>& inputArray) 158 { 159 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) { 160 fSteppingArray[i] = inputArray[i]; 161 } 162 } 163 164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 165 166 void Run::SetTrackingArray1( 167 const std::array<G4long, TrackingAction::fkNumberCombinations>& inputArray) 168 { 169 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 170 fTrackingArray1[i] = inputArray[i]; 171 } 172 } 173 174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 175 176 void Run::SetTrackingArray2( 177 const std::array<G4double, TrackingAction::fkNumberCombinations>& inputArray) 178 { 179 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 180 fTrackingArray2[i] = inputArray[i]; 181 } 182 } 183 184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 185