Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-06-01 07:54:45

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file Run.cc
0027 /// \brief Implementation of the Run class
0028 
0029 #include "Run.hh"
0030 
0031 #include "G4Run.hh"
0032 #include "G4RunManager.hh"
0033 #include "G4SystemOfUnits.hh"
0034 
0035 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0036 
0037 Run::Run()
0038   : G4Run(),
0039     fNumEvents(0),
0040     fPrimaryParticleId(0),
0041     fPrimaryParticleEnergy(0.0),
0042     fPrimaryParticleDirection(G4ThreeVector(0.0, 0.0, 0.0)),
0043     fTargetMaterialName(""),
0044     fCubicVolumeScoringUpDown(1.0),
0045     fCubicVolumeScoringSide(1.0)
0046 {
0047   fSteppingArray.fill(0.0);
0048   fTrackingArray1.fill(0);
0049   fTrackingArray2.fill(0.0);
0050 }
0051 
0052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0053 
0054 void Run::RecordEvent(const G4Event* anEvent)
0055 {
0056   // This method is called automatically by the Geant4 kernel (not by the user!) at the end
0057   // of each event : in MT-mode, it is called only for the working thread that handled the event.
0058   G4int nEvt = anEvent->GetEventID();
0059   if (nEvt % 10 == 0) G4cout << " Event#=" << nEvt << G4endl;
0060   G4Run::RecordEvent(anEvent);
0061 }
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 void Run::Merge(const G4Run* aRun)
0066 {
0067   // This method is called automatically by the Geant4 kernel (not by the user!) only in the case
0068   // of multithreaded mode and only for working threads.
0069   const Run* localRun = static_cast<const Run*>(aRun);
0070   fPrimaryParticleId = localRun->GetPrimaryParticleId();
0071   fPrimaryParticleEnergy = localRun->GetPrimaryParticleEnergy();
0072   fPrimaryParticleDirection = localRun->GetPrimaryParticleDirection();
0073   fTargetMaterialName = localRun->GetTargetMaterialName();
0074   fCubicVolumeScoringUpDown = localRun->GetCubicVolumeScoringUpDown();
0075   fCubicVolumeScoringSide = localRun->GetCubicVolumeScoringSide();
0076   fNumEvents += localRun->GetNumberOfEvent();
0077   for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0078     fSteppingArray[i] += localRun->GetSteppingArray()[i];
0079   }
0080   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0081     fTrackingArray1[i] += localRun->GetTrackingArray1()[i];
0082     fTrackingArray2[i] += localRun->GetTrackingArray2()[i];
0083   }
0084   G4Run::Merge(aRun);
0085 }
0086 
0087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0088 
0089 void Run::PrintInfo() const
0090 {
0091   // This method is called by RunAction::EndOfRunAction. In MT-mode, only the master thread
0092   // calls it.
0093   const G4double floatingNumberOfEvents =
0094     std::max(1.0, fNumEvents > 0 ? fNumEvents * 1.0 : GetNumberOfEvent() * 1.0);
0095   // The fluence in the scoring volume is defined as sum of step lengths in that volume
0096   // divided by the volume of that scoring volume.
0097   const G4double conversionFactor = CLHEP::cm * CLHEP::cm;  // From  mm^-2  to  cm^-2
0098   const G4double factorUpDown =
0099     conversionFactor / (fCubicVolumeScoringUpDown * floatingNumberOfEvents);
0100   const G4double factorSide = conversionFactor / (fCubicVolumeScoringSide * floatingNumberOfEvents);
0101   G4cout << std::setprecision(6) << G4endl << G4endl
0102          << " ===============  Run::PrintInfo()  =============== \t RunID = " << GetRunID()
0103          << G4endl << " Primary particle PDG code = " << fPrimaryParticleId << G4endl
0104          << " Primary particle kinetic energy = " << fPrimaryParticleEnergy / CLHEP::GeV << " GeV"
0105          << G4endl << " Primary particle direction = " << fPrimaryParticleDirection << G4endl
0106          << " Target material = " << fTargetMaterialName << G4endl
0107          << " Cubic-volume scoring up-down = " << fCubicVolumeScoringUpDown << " mm^3" << G4endl
0108          << " Cubic-volume scoring side    = " << fCubicVolumeScoringSide << " mm^3" << G4endl
0109          << " Number of events = " << floatingNumberOfEvents << G4endl
0110          << " Conversion factor: fluence from mm^-2 to cm^-2 = " << conversionFactor << G4endl
0111          << " Particle fluence in unit of cm^-2 :" << G4endl;
0112   for (G4int i = 0; i < SteppingAction::fkNumberScoringVolumes; ++i) {
0113     G4double factor = (i == 1 ? factorSide : factorUpDown);
0114     for (G4int j = 0; j < SteppingAction::fkNumberKinematicRegions; ++j) {
0115       for (G4int k = 0; k < SteppingAction::fkNumberParticleTypes; ++k) {
0116         G4int index = SteppingAction::GetIndex(i, j, k);
0117         // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
0118         G4cout << "   case=" << std::setw(3) << index << "   " << std::setw(12)
0119                << SteppingAction::fkArrayScoringVolumeNames[i] << "   " << std::setw(12)
0120                << SteppingAction::fkArrayKinematicRegionNames[j] << "   " << std::setw(12)
0121                << SteppingAction::fkArrayParticleTypeNames[k] << "   " << std::setw(8)
0122                << factor * fSteppingArray[index] << G4endl;
0123       }
0124     }
0125   }
0126   G4cout << " ------------------------------------------------------------- " << G4endl
0127          << " Extra information: particle production \t \t      <N>      <E_kin> <Sum_Ekin> [MeV]"
0128          << G4endl;
0129   const G4double normalization = 1.0 / floatingNumberOfEvents;
0130   for (G4int i = 0; i < TrackingAction::fkNumberScoringVolumes; ++i) {
0131     for (G4int j = 0; j < TrackingAction::fkNumberKinematicRegions; ++j) {
0132       for (G4int k = 0; k < TrackingAction::fkNumberParticleTypes; ++k) {
0133         G4int index = TrackingAction::GetIndex(i, j, k);
0134         // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
0135         G4cout << "   case=" << std::setw(3) << index << "   " << std::setw(12)
0136                << TrackingAction::fkArrayScoringVolumeNames[i] << "   " << std::setw(12)
0137                << TrackingAction::fkArrayKinematicRegionNames[j] << "   " << std::setw(12)
0138                << TrackingAction::fkArrayParticleTypeNames[k] << "   " << std::setw(8)
0139                << normalization * fTrackingArray1[index] << "   " << std::setw(8)
0140                << (fTrackingArray1[index] > 0 ? fTrackingArray2[index] / fTrackingArray1[index]
0141                                               : 0.0)
0142                << "   " << std::setw(8) << normalization * fTrackingArray2[index] << G4endl;
0143       }
0144     }
0145   }
0146   G4cout << " ============================================================= " << G4endl << G4endl;
0147 }
0148 
0149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0150 
0151 void Run::SetSteppingArray(
0152   const std::array<G4double, SteppingAction::fkNumberCombinations>& inputArray)
0153 {
0154   for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0155     fSteppingArray[i] = inputArray[i];
0156   }
0157 }
0158 
0159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0160 
0161 void Run::SetTrackingArray1(
0162   const std::array<G4long, TrackingAction::fkNumberCombinations>& inputArray)
0163 {
0164   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0165     fTrackingArray1[i] = inputArray[i];
0166   }
0167 }
0168 
0169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0170 
0171 void Run::SetTrackingArray2(
0172   const std::array<G4double, TrackingAction::fkNumberCombinations>& inputArray)
0173 {
0174   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0175     fTrackingArray2[i] = inputArray[i];
0176   }
0177 }
0178 
0179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......