Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/ParticleFluence/Sphere/src/Run.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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