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 fTrackerMaterialName(""), 49 fEmCaloMaterialName(""), 50 fHadCaloMaterialName(""), 51 fCubicVolumeScoringTrackerShell(1.0), 52 fCubicVolumeScoringEmCaloShell(1.0), 53 fCubicVolumeScoringHadCaloShell(1.0) 54 { 55 fSteppingArray.fill(0.0); 56 fTrackingArray1.fill(0); 57 fTrackingArray2.fill(0.0); 58 } 59 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 61 62 void Run::RecordEvent(const G4Event* anEvent) 63 { 64 // This method is called automatically by the Geant4 kernel (not by the user!) at the end 65 // of each event : in MT-mode, it is called only for the working thread that handled the event. 66 G4int nEvt = anEvent->GetEventID(); 67 if (nEvt % 10 == 0) G4cout << " Event#=" << nEvt << G4endl; 68 G4Run::RecordEvent(anEvent); 69 } 70 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 72 73 void Run::Merge(const G4Run* aRun) 74 { 75 // This method is called automatically by the Geant4 kernel (not by the user!) only in the case 76 // of multithreaded mode and only for working threads. 77 const Run* localRun = static_cast<const Run*>(aRun); 78 fPrimaryParticleId = localRun->GetPrimaryParticleId(); 79 fPrimaryParticleEnergy = localRun->GetPrimaryParticleEnergy(); 80 fPrimaryParticleDirection = localRun->GetPrimaryParticleDirection(); 81 fTrackerMaterialName = localRun->GetTrackerMaterialName(); 82 fEmCaloMaterialName = localRun->GetEmCaloMaterialName(); 83 fHadCaloMaterialName = localRun->GetHadCaloMaterialName(); 84 fCubicVolumeScoringTrackerShell = localRun->GetCubicVolumeScoringTrackerShell(); 85 fCubicVolumeScoringEmCaloShell = localRun->GetCubicVolumeScoringEmCaloShell(); 86 fCubicVolumeScoringHadCaloShell = localRun->GetCubicVolumeScoringHadCaloShell(); 87 fNumEvents += localRun->GetNumberOfEvent(); 88 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) { 89 fSteppingArray[i] += localRun->GetSteppingArray()[i]; 90 } 91 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 92 fTrackingArray1[i] += localRun->GetTrackingArray1()[i]; 93 fTrackingArray2[i] += localRun->GetTrackingArray2()[i]; 94 } 95 G4Run::Merge(aRun); 96 } 97 98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 99 100 void Run::PrintInfo() const 101 { 102 // This method is called by RunAction::EndOfRunAction. 103 // In MT-mode, only the master thread calls it. 104 const G4double floatingNumberOfEvents = 105 std::max(1.0, fNumEvents > 0 ? fNumEvents * 1.0 : GetNumberOfEvent() * 1.0); 106 // The fluence in a scoring shell is defined as sum of step lengths in that shell 107 // divided by the cubic-volume of the scoring shell. 108 const G4double conversionFactor = CLHEP::cm * CLHEP::cm; // From mm^-2 to cm^-2 109 const G4double factorTracker = 110 conversionFactor / (0.5 * fCubicVolumeScoringTrackerShell * floatingNumberOfEvents); 111 const G4double factorEmCalo = 112 conversionFactor / (0.5 * fCubicVolumeScoringEmCaloShell * floatingNumberOfEvents); 113 const G4double factorHadCalo = 114 conversionFactor / (0.5 * fCubicVolumeScoringHadCaloShell * floatingNumberOfEvents); 115 G4cout << std::setprecision(6) << G4endl << G4endl 116 << " =============== Run::PrintInfo() =============== \t RunID = " << GetRunID() 117 << G4endl << " Primary particle PDG code = " << fPrimaryParticleId << G4endl 118 << " Primary particle kinetic energy = " << fPrimaryParticleEnergy / CLHEP::GeV << " GeV" 119 << G4endl << " Primary particle direction = " << fPrimaryParticleDirection << G4endl 120 << " Tracker material = " << fTrackerMaterialName << G4endl 121 << " EmCalo material = " << fEmCaloMaterialName << G4endl 122 << " HadCalo material = " << fHadCaloMaterialName << G4endl 123 << " Cubic-volume scoring tracker shell = " << fCubicVolumeScoringTrackerShell << " mm^3" 124 << G4endl << " Cubic-volume scoring emCalo shell = " << fCubicVolumeScoringEmCaloShell 125 << " mm^3" << G4endl 126 << " Cubic-volume scoring hadCalo shell = " << fCubicVolumeScoringHadCaloShell << " mm^3" 127 << G4endl << " Number of events = " << floatingNumberOfEvents << G4endl 128 << " Conversion factor: fluence from mm^-2 to cm^-2 = " << conversionFactor << G4endl 129 << " Particle fluence in unit of cm^-2 :" << G4endl; 130 for (G4int i = 0; i < SteppingAction::fkNumberScoringShells; ++i) { 131 G4double factor = factorTracker; 132 if (i == 1) 133 factor = factorEmCalo; 134 else if (i == 2) 135 factor = factorHadCalo; 136 for (G4int j = 0; j < SteppingAction::fkNumberKinematicRegions; ++j) { 137 for (G4int k = 0; k < SteppingAction::fkNumberScoringPositions; ++k) { 138 for (G4int ll = 0; ll < SteppingAction::fkNumberParticleTypes; ++ll) { 139 G4int index = SteppingAction::GetIndex(i, j, k, ll); 140 // G4cout << "(i, j, k, ll)=(" << i << ", " << j << ", " << k << ", " 141 // << ll << ") ->" << index; 142 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12) 143 << SteppingAction::fkArrayScoringShellNames[i] << " " << std::setw(12) 144 << SteppingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12) 145 << SteppingAction::fkArrayScoringPositionNames[k] << " " << std::setw(12) 146 << SteppingAction::fkArrayParticleTypeNames[ll] << " " << std::setw(8) 147 << factor * fSteppingArray[index] << G4endl; 148 } 149 } 150 } 151 } 152 G4cout << " ------------------------------------------------------------- " << G4endl 153 << " Extra information: particle production \t \t <N> <E_kin> <Sum_Ekin> [MeV]" 154 << G4endl; 155 const G4double normalization = 1.0 / floatingNumberOfEvents; 156 for (G4int i = 0; i < TrackingAction::fkNumberScoringVolumes; ++i) { 157 for (G4int j = 0; j < TrackingAction::fkNumberKinematicRegions; ++j) { 158 for (G4int k = 0; k < TrackingAction::fkNumberParticleTypes; ++k) { 159 G4int index = TrackingAction::GetIndex(i, j, k); 160 // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ") ->" << index; 161 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12) 162 << TrackingAction::fkArrayScoringVolumeNames[i] << " " << std::setw(12) 163 << TrackingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12) 164 << TrackingAction::fkArrayParticleTypeNames[k] << " " << std::setw(8) 165 << normalization * fTrackingArray1[index] << " " << std::setw(8) 166 << (fTrackingArray1[index] > 0 ? fTrackingArray2[index] / fTrackingArray1[index] 167 : 0.0) 168 << " " << std::setw(8) << normalization * fTrackingArray2[index] << G4endl; 169 } 170 } 171 } 172 G4cout << " ============================================================= " << G4endl << G4endl; 173 } 174 175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 176 177 void Run::SetSteppingArray( 178 const std::array<G4double, SteppingAction::fkNumberCombinations>& inputArray) 179 { 180 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) { 181 fSteppingArray[i] = inputArray[i]; 182 } 183 } 184 185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 186 187 void Run::SetTrackingArray1( 188 const std::array<G4long, TrackingAction::fkNumberCombinations>& inputArray) 189 { 190 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 191 fTrackingArray1[i] = inputArray[i]; 192 } 193 } 194 195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 196 197 void Run::SetTrackingArray2( 198 const std::array<G4double, TrackingAction::fkNumberCombinations>& inputArray) 199 { 200 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) { 201 fTrackingArray2[i] = inputArray[i]; 202 } 203 } 204 205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 206