Geant4 Cross Reference

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