Back to home page

EIC code displayed by LXR

 
 

    


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

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