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 SteppingAction.cc
0027 /// \brief Implementation of the SteppingAction class
0028 //
0029 //
0030 
0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0032 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0033 
0034 #include "SteppingAction.hh"
0035 
0036 #include "Run.hh"
0037 
0038 #include "G4IonTable.hh"
0039 #include "G4LossTableManager.hh"
0040 #include "G4ParticleDefinition.hh"
0041 #include "G4ParticleTypes.hh"
0042 #include "G4Step.hh"
0043 #include "G4StepPoint.hh"
0044 #include "G4SystemOfUnits.hh"
0045 #include "G4TouchableHistory.hh"
0046 #include "G4Track.hh"
0047 #include "G4VPhysicalVolume.hh"
0048 #include "G4VSolid.hh"
0049 #include "G4VTouchable.hh"
0050 
0051 const std::array<G4String, SteppingAction::fkNumberScoringVolumes>
0052   SteppingAction::fkArrayScoringVolumeNames = {"downstream", "side", "upstream"};
0053 
0054 const std::array<G4String, SteppingAction::fkNumberKinematicRegions>
0055   SteppingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
0056 
0057 const std::array<G4String, SteppingAction::fkNumberParticleTypes>
0058   SteppingAction::fkArrayParticleTypeNames = {"all",      "electron",   "gamma",      "muon",
0059                                               "neutrino", "pion",       "neutron",    "proton",
0060                                               "ion",      "otherMeson", "otherBaryon"};
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0063 
0064 G4int SteppingAction::GetIndex(const G4int iScoringVolume, const G4int iKinematicRegion,
0065                                const G4int iParticleType)
0066 {
0067   G4int index = -1;
0068   if (iScoringVolume >= 0 && iScoringVolume < fkNumberScoringVolumes && iKinematicRegion >= 0
0069       && iKinematicRegion < fkNumberKinematicRegions && iParticleType >= 0
0070       && iParticleType < fkNumberParticleTypes)
0071   {
0072     index = iScoringVolume * fkNumberKinematicRegions * fkNumberParticleTypes
0073             + iKinematicRegion * fkNumberParticleTypes + iParticleType;
0074   }
0075   if (index < 0 || index >= fkNumberCombinations) {
0076     G4cerr << "SteppingAction::GetIndex : WRONG index=" << index << "  set it to 0 !" << G4endl;
0077     index = 0;
0078   }
0079   return index;
0080 }
0081 
0082 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0083 
0084 SteppingAction::SteppingAction() : G4UserSteppingAction()
0085 {
0086   Initialize();
0087 }
0088 
0089 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0090 
0091 void SteppingAction::Initialize()
0092 {
0093   // Initialization needed at the beginning of each Run
0094   fPrimaryParticleId = 0;
0095   fPrimaryParticleEnergy = 0.0;
0096   fPrimaryParticleDirection = G4ThreeVector(0.0, 0.0, 1.0);
0097   fAbsorberMaterialName = "";
0098   fActiveMaterialName = "";
0099   fIsFirstStepOfTheEvent = true;
0100   fIsFirstStepInAbsorberLayer = true;
0101   fIsFirstStepInActiveLayer = true;
0102   fIsFirstStepInScoringUpDown = true;
0103   fIsFirstStepInScoringSide = true;
0104   fCubicVolumeScoringUpDown = 1.0;
0105   fCubicVolumeScoringSide = 1.0;
0106   for (G4int i = 0; i < fkNumberCombinations; ++i) {
0107     fArraySumStepLengths[i] = 0.0;
0108   }
0109   /*
0110   for ( G4int i = 0; i < fkNumberCombinations; ++i ) fArraySumStepLengths[i] = 999.9;
0111   G4cout << " fkNumberCombinations=" << fkNumberCombinations << G4endl;
0112   for ( G4int i = 0; i < fkNumberScoringVolumes; ++i ) {
0113     for ( G4int j = 0; j < fkNumberKinematicRegions; ++j ) {
0114       for ( G4int k = 0; k < fkNumberParticleTypes; ++k ) {
0115         G4int index = GetIndex( i, j, k );
0116         G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
0117         if ( fArraySumStepLengths[ index ] < 1.0 ) G4cout << " <=== REPEATED!";
0118         else                                       fArraySumStepLengths[ index ] = 0.0;
0119         G4cout << G4endl;
0120       }
0121     }
0122   }
0123   for ( G4int i = 0; i < fkNumberCombinations; ++i ) {
0124     if ( fArraySumStepLengths[i] > 999.0 ) G4cout << " i=" << i << " NOT COVERED !" << G4endl;
0125   }
0126   */
0127 }
0128 
0129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0130 
0131 void SteppingAction::UserSteppingAction(const G4Step* theStep)
0132 {
0133   // Get information on the primary particle
0134   if (fIsFirstStepOfTheEvent) {
0135     if (theStep->GetTrack()->GetParentID() == 0) {
0136       fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0137       fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
0138       fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
0139       if (fRunPtr) {
0140         fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
0141         fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
0142         fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
0143       }
0144       fIsFirstStepOfTheEvent = false;
0145     }
0146   }
0147   // Get information on the materials of the calorimeter
0148   if (fIsFirstStepInAbsorberLayer
0149       && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiAbsorber")
0150   {
0151     fAbsorberMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0152     if (fRunPtr) fRunPtr->SetAbsorberMaterialName(fAbsorberMaterialName);
0153     fIsFirstStepInAbsorberLayer = false;
0154   }
0155   if (fIsFirstStepInActiveLayer
0156       && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiActive")
0157   {
0158     fActiveMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0159     if (fRunPtr) fRunPtr->SetActiveMaterialName(fActiveMaterialName);
0160     fIsFirstStepInActiveLayer = false;
0161   }
0162   // Get information on step lengths in the scoring volumes
0163   G4int iScoringVolume = -1;
0164   if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringDownstream") {
0165     iScoringVolume = 0;
0166     if (fIsFirstStepInScoringUpDown) {
0167       fCubicVolumeScoringUpDown =
0168         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0169       if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0170       fIsFirstStepInScoringUpDown = false;
0171     }
0172   }
0173   else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringSide") {
0174     iScoringVolume = 1;
0175     if (fIsFirstStepInScoringSide) {
0176       fCubicVolumeScoringSide =
0177         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0178       if (fRunPtr) fRunPtr->SetCubicVolumeScoringSide(fCubicVolumeScoringSide);
0179       fIsFirstStepInScoringSide = false;
0180     }
0181   }
0182   else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringUpstream") {
0183     iScoringVolume = 2;
0184     if (fIsFirstStepInScoringUpDown) {
0185       fCubicVolumeScoringUpDown =
0186         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0187       if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0188       fIsFirstStepInScoringUpDown = false;
0189     }
0190   }
0191   if (iScoringVolume >= 0) {
0192     // In the case of the upstream scoring volume, consider only particles whose direction
0193     // is opposite with respect to the primary particle (this is needed, in particular,
0194     // for avoiding to account the incoming, primary beam particle in the "upstream" fluence).
0195     if (iScoringVolume == 2
0196         && fPrimaryParticleDirection.dot(theStep->GetPreStepPoint()->GetMomentumDirection()) > 0.0)
0197       return;
0198     G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
0199     G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
0200                      ? 0
0201                      : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
0202     /*
0203     G4cout << std::setprecision(6)
0204            << theStep->GetTrack()->GetDefinition()->GetParticleName() << "  absPdg=" << absPdg
0205            << "  Ekin[MeV]=" << theStep->GetPreStepPoint()->GetKineticEnergy()
0206            << "  (rho,z)[mm]=(" << theStep->GetTrack()->GetPosition().perp()
0207            << "," << theStep->GetTrack()->GetPosition().z() << ")"
0208            << "  " << theStep->GetTrack()->GetVolume()->GetName()
0209            << "  " << theStep->GetTrack()->GetMaterial()->GetName()
0210            << "  L[mm]=" << stepLength << "  "
0211            << ( fPrimaryParticleDirection.dot(
0212                   theStep->GetPreStepPoint()->GetMomentumDirection() ) > 0.0
0213               ? "forward" : "backward" )
0214            << G4endl;
0215     */
0216     // Three kinematical regions:  [0] : any value ;  [1] : below 20 MeV ;  [2] : above 20 MeV
0217     G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
0218     G4int iParticleType = -1;
0219     if (absPdg == 11)
0220       iParticleType = 1;  // electron (and positron)
0221     else if (absPdg == 22)
0222       iParticleType = 2;  // gamma
0223     else if (absPdg == 13)
0224       iParticleType = 3;  // muons (mu- and mu+)
0225     else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
0226       iParticleType = 4;  // neutrinos
0227     // (and anti-neutrinos), all flavors
0228     else if (absPdg == 111 || absPdg == 211)
0229       iParticleType = 5;  // (charged) pions
0230     else if (absPdg == 2112)
0231       iParticleType = 6;  // neutron (and anti-neutron)
0232     else if (absPdg == 2212)
0233       iParticleType = 7;  // proton  (and anti-proton)
0234     else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||  // ions (and anti-ions)
0235              G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
0236       iParticleType = 8;
0237     else if (absPdg < 1000)
0238       iParticleType = 9;  // other mesons (e.g. kaons) (Note: this works
0239                           // in most cases, but not always!)
0240     else if (absPdg > 1000)
0241       iParticleType = 10;  // other baryons (e.g. hyperons, anti-hyperons,
0242                            // etc.)
0243     // Consider the specific case : scoring volume, kinematic region and particle type
0244     G4int index = GetIndex(iScoringVolume, iKinematicRegion, iParticleType);
0245     fArraySumStepLengths[index] += stepLength;
0246     // Consider the "all" particle case, with the same scoring volume and kinematic region
0247     index = GetIndex(iScoringVolume, iKinematicRegion, 0);
0248     fArraySumStepLengths[index] += stepLength;
0249     // Consider the "any" kinematic region case, with the same scoring volume and particle type
0250     index = GetIndex(iScoringVolume, 0, iParticleType);
0251     fArraySumStepLengths[index] += stepLength;
0252     // Consider the "any" kinematic region and "all" particle, with the same scoring volume
0253     index = GetIndex(iScoringVolume, 0, 0);
0254     fArraySumStepLengths[index] += stepLength;
0255     if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
0256   }
0257 }
0258 
0259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......