Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-10 08:06:16

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 //
0030 
0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0032 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0033 
0034 #include "Run.hh"
0035 
0036 #include "G4Run.hh"
0037 #include "G4RunManager.hh"
0038 #include "G4SystemOfUnits.hh"
0039 
0040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0041 
0042 Run::Run()
0043   : G4Run(),
0044     fNumEvents(0),
0045     fPrimaryParticleId(0),
0046     fPrimaryParticleEnergy(0.0),
0047     fPrimaryParticleDirection(G4ThreeVector(0.0, 0.0, 0.0)),
0048     fAbsorberMaterialName(""),
0049     fActiveMaterialName(""),
0050     fCubicVolumeScoringUpDown(1.0),
0051     fCubicVolumeScoringSide(1.0)
0052 {
0053   fSteppingArray.fill(0.0);
0054   fTrackingArray1.fill(0);
0055   fTrackingArray2.fill(0.0);
0056 }
0057 
0058 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0059 
0060 void Run::RecordEvent(const G4Event* anEvent)
0061 {
0062   // This method is called automatically by the Geant4 kernel (not by the user!) at the end
0063   // of each event : in MT-mode, it is called only for the working thread that handled the event.
0064   G4int nEvt = anEvent->GetEventID();
0065   if (nEvt % 10 == 0) G4cout << " Event#=" << nEvt << G4endl;
0066   G4Run::RecordEvent(anEvent);
0067 }
0068 
0069 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0070 
0071 void Run::Merge(const G4Run* aRun)
0072 {
0073   // This method is called automatically by the Geant4 kernel (not by the user!) only in the case
0074   // of multithreaded mode and only for working threads.
0075   const Run* localRun = static_cast<const Run*>(aRun);
0076   fPrimaryParticleId = localRun->GetPrimaryParticleId();
0077   fPrimaryParticleEnergy = localRun->GetPrimaryParticleEnergy();
0078   fPrimaryParticleDirection = localRun->GetPrimaryParticleDirection();
0079   fAbsorberMaterialName = localRun->GetAbsorberMaterialName();
0080   fActiveMaterialName = localRun->GetActiveMaterialName();
0081   fCubicVolumeScoringUpDown = localRun->GetCubicVolumeScoringUpDown();
0082   fCubicVolumeScoringSide = localRun->GetCubicVolumeScoringSide();
0083   fNumEvents += localRun->GetNumberOfEvent();
0084   for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0085     fSteppingArray[i] += localRun->GetSteppingArray()[i];
0086   }
0087   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0088     fTrackingArray1[i] += localRun->GetTrackingArray1()[i];
0089     fTrackingArray2[i] += localRun->GetTrackingArray2()[i];
0090   }
0091   G4Run::Merge(aRun);
0092 }
0093 
0094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0095 
0096 void Run::PrintInfo() const
0097 {
0098   // This method is called by RunAction::EndOfRunAction. In MT-mode, only the master thread
0099   // calls it.
0100   const G4double floatingNumberOfEvents =
0101     std::max(1.0, fNumEvents > 0 ? fNumEvents * 1.0 : GetNumberOfEvent() * 1.0);
0102   // The fluence in the scoring volume is defined as sum of step lengths in that volume
0103   // divided by the volume of that scoring volume.
0104   const G4double conversionFactor = CLHEP::cm * CLHEP::cm;  // From  mm^-2  to  cm^-2
0105   const G4double factorUpDown =
0106     conversionFactor / (fCubicVolumeScoringUpDown * floatingNumberOfEvents);
0107   const G4double factorSide = conversionFactor / (fCubicVolumeScoringSide * floatingNumberOfEvents);
0108   G4cout << std::setprecision(6) << G4endl << G4endl
0109          << " ===============  Run::PrintInfo()  =============== \t RunID = " << GetRunID()
0110          << G4endl << " Primary particle PDG code = " << fPrimaryParticleId << G4endl
0111          << " Primary particle kinetic energy = " << fPrimaryParticleEnergy / CLHEP::GeV << " GeV"
0112          << G4endl << " Primary particle direction = " << fPrimaryParticleDirection << G4endl
0113          << " Absorber material = " << fAbsorberMaterialName << G4endl
0114          << " Active material   = " << fActiveMaterialName << G4endl
0115          << " Cubic-volume scoring up-down = " << fCubicVolumeScoringUpDown << " mm^3" << G4endl
0116          << " Cubic-volume scoring side    = " << fCubicVolumeScoringSide << " mm^3" << G4endl
0117          << " Number of events = " << floatingNumberOfEvents << G4endl
0118          << " Conversion factor: fluence from mm^-2 to cm^-2 = " << conversionFactor << G4endl
0119          << " Particle fluence in unit of cm^-2 :" << G4endl;
0120   for (G4int i = 0; i < SteppingAction::fkNumberScoringVolumes; ++i) {
0121     G4double factor = (i == 1 ? factorSide : factorUpDown);
0122     for (G4int j = 0; j < SteppingAction::fkNumberKinematicRegions; ++j) {
0123       for (G4int k = 0; k < SteppingAction::fkNumberParticleTypes; ++k) {
0124         G4int index = SteppingAction::GetIndex(i, j, k);
0125         // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
0126         G4cout << "   case=" << std::setw(3) << index << "   " << std::setw(12)
0127                << SteppingAction::fkArrayScoringVolumeNames[i] << "   " << std::setw(12)
0128                << SteppingAction::fkArrayKinematicRegionNames[j] << "   " << std::setw(12)
0129                << SteppingAction::fkArrayParticleTypeNames[k] << "   " << std::setw(8)
0130                << factor * fSteppingArray[index] << G4endl;
0131       }
0132     }
0133   }
0134   G4cout << " ------------------------------------------------------------- " << G4endl
0135          << " Extra information: particle production \t \t      <N>      <E_kin> <Sum_Ekin> [MeV]"
0136          << G4endl;
0137   const G4double normalization = 1.0 / floatingNumberOfEvents;
0138   for (G4int i = 0; i < TrackingAction::fkNumberScoringVolumes; ++i) {
0139     for (G4int j = 0; j < TrackingAction::fkNumberKinematicRegions; ++j) {
0140       for (G4int k = 0; k < TrackingAction::fkNumberParticleTypes; ++k) {
0141         G4int index = TrackingAction::GetIndex(i, j, k);
0142         // G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
0143         G4cout << "   case=" << std::setw(3) << index << "   " << std::setw(12)
0144                << TrackingAction::fkArrayScoringVolumeNames[i] << "   " << std::setw(12)
0145                << TrackingAction::fkArrayKinematicRegionNames[j] << "   " << std::setw(12)
0146                << TrackingAction::fkArrayParticleTypeNames[k] << "   " << std::setw(8)
0147                << normalization * fTrackingArray1[index] << "   " << std::setw(8)
0148                << (fTrackingArray1[index] > 0 ? fTrackingArray2[index] / fTrackingArray1[index]
0149                                               : 0.0)
0150                << "   " << std::setw(8) << normalization * fTrackingArray2[index] << G4endl;
0151       }
0152     }
0153   }
0154   G4cout << " ============================================================= " << G4endl << G4endl;
0155 }
0156 
0157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0158 
0159 void Run::SetSteppingArray(
0160   const std::array<G4double, SteppingAction::fkNumberCombinations>& inputArray)
0161 {
0162   for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0163     fSteppingArray[i] = inputArray[i];
0164   }
0165 }
0166 
0167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0168 
0169 void Run::SetTrackingArray1(
0170   const std::array<G4long, TrackingAction::fkNumberCombinations>& inputArray)
0171 {
0172   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0173     fTrackingArray1[i] = inputArray[i];
0174   }
0175 }
0176 
0177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0178 
0179 void Run::SetTrackingArray2(
0180   const std::array<G4double, TrackingAction::fkNumberCombinations>& inputArray)
0181 {
0182   for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0183     fTrackingArray2[i] = inputArray[i];
0184   }
0185 }
0186 
0187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......