Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/ParticleFluence/Sphere/src/TrackingAction.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 TrackingAction.cc
 27 /// \brief Implementation of the TrackingAction class
 28 //
 29 //
 30 
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 
 34 #include "TrackingAction.hh"
 35 
 36 #include "Run.hh"
 37 
 38 #include "G4IonTable.hh"
 39 #include "G4ParticleDefinition.hh"
 40 #include "G4ParticleTypes.hh"
 41 #include "G4Step.hh"
 42 #include "G4StepPoint.hh"
 43 #include "G4SystemOfUnits.hh"
 44 #include "G4Track.hh"
 45 
 46 const std::array<G4String, TrackingAction::fkNumberScoringVolumes>
 47   TrackingAction::fkArrayScoringVolumeNames = {"sphere"};
 48 
 49 const std::array<G4String, TrackingAction::fkNumberKinematicRegions>
 50   TrackingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
 51 
 52 const std::array<G4String, TrackingAction::fkNumberParticleTypes>
 53   TrackingAction::fkArrayParticleTypeNames = {"all",      "electron",   "gamma",      "muon",
 54                                               "neutrino", "pion",       "neutron",    "proton",
 55                                               "ion",      "otherMeson", "otherBaryon"};
 56 
 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 58 
 59 G4int TrackingAction::GetIndex(const G4int iScoringVolume, const G4int iKinematicRegion,
 60                                const G4int iParticleType)
 61 {
 62   G4int index = -1;
 63   if (iScoringVolume >= 0 && iScoringVolume < fkNumberScoringVolumes && iKinematicRegion >= 0
 64       && iKinematicRegion < fkNumberKinematicRegions && iParticleType >= 0
 65       && iParticleType < fkNumberParticleTypes)
 66   {
 67     index = iScoringVolume * fkNumberKinematicRegions * fkNumberParticleTypes
 68             + iKinematicRegion * fkNumberParticleTypes + iParticleType;
 69   }
 70   return index;
 71 }
 72 
 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 74 
 75 TrackingAction::TrackingAction() : G4UserTrackingAction()
 76 {
 77   Initialize();
 78 }
 79 
 80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 81 
 82 void TrackingAction::Initialize()
 83 {
 84   // Initialization needed at the beginning of each Run
 85   fArrayMultiplicities.fill(0);
 86   fArraySumKineticEnergies.fill(0.0);
 87 }
 88 
 89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 90 
 91 void TrackingAction::PreUserTrackingAction(const G4Track* aTrack)
 92 {
 93   // This method is called not only once when a particle is created,
 94   // but also each time it is resumed, in the case the track gets suspended,
 95   // as it happens in the case of neutrons with _HP Physics Lists.
 96   // To be sure that we collect information about a track one and only once,
 97   // we require that the current step be the first one.
 98   if (aTrack == nullptr || aTrack->GetCurrentStepNumber() != 0 || aTrack->GetDefinition() == nullptr
 99       || aTrack->GetLogicalVolumeAtVertex() == nullptr
100       || aTrack->GetLogicalVolumeAtVertex()->GetName() != "logicSphere")
101     return;
102   G4int iScoringVolume = 0;
103   // Three kinematical regions:  [0] : any value ;  [1] : below 20 MeV ;  [2] : above 20 MeV
104   G4int iKinematicRegion = aTrack->GetKineticEnergy() < 20.0 ? 1 : 2;
105   G4int absPdg = std::abs(aTrack->GetDefinition()->GetPDGEncoding());
106   G4int iParticleType = -1;
107   if (absPdg == 11)
108     iParticleType = 1;  // electron (and positron)
109   else if (absPdg == 22)
110     iParticleType = 2;  // gamma
111   else if (absPdg == 13)
112     iParticleType = 3;  // muons (mu- and mu+)
113   else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
114     iParticleType = 4;
115   // neutrinos (and anti-neutrinos), all flavors
116   else if (absPdg == 111 || absPdg == 211)
117     iParticleType = 5;  // (charged) pions
118   else if (absPdg == 2112)
119     iParticleType = 6;  // neutron (and anti-neutron)
120   else if (absPdg == 2212)
121     iParticleType = 7;  // proton  (and anti-proton)
122   else if (G4IonTable::IsIon(aTrack->GetDefinition())
123            || G4IonTable::IsAntiIon(aTrack->GetDefinition()))
124     iParticleType = 8;
125   // ions (and anti-ions)
126   else if (absPdg < 1000)
127     iParticleType = 9;  // other mesons (e.g. kaons)
128   // (Note: this works in most cases, but not always!)
129   else if (absPdg > 1000)
130     iParticleType = 10;  // other baryons (e.g. hyperons,
131                          // anti-hyperons, etc.)
132   // Consider the specific case : scoring volume, kinematic region and particle type
133   G4int index = GetIndex(iScoringVolume, iKinematicRegion, iParticleType);
134   ++fArrayMultiplicities[index];
135   fArraySumKineticEnergies[index] += aTrack->GetKineticEnergy();
136   // Consider the "all" particle case, with the same scoring volume and kinematic region
137   index = GetIndex(iScoringVolume, iKinematicRegion, 0);
138   ++fArrayMultiplicities[index];
139   fArraySumKineticEnergies[index] += aTrack->GetKineticEnergy();
140   // Consider the "any" kinematic region case, with the same scoring volume and particle type
141   index = GetIndex(iScoringVolume, 0, iParticleType);
142   ++fArrayMultiplicities[index];
143   fArraySumKineticEnergies[index] += aTrack->GetKineticEnergy();
144   // Consider the "any" kinematic region and "all" particle, with the same scoring volume
145   index = GetIndex(iScoringVolume, 0, 0);
146   ++fArrayMultiplicities[index];
147   fArraySumKineticEnergies[index] += aTrack->GetKineticEnergy();
148   if (fRunPtr) {
149     fRunPtr->SetTrackingArray1(fArrayMultiplicities);
150     fRunPtr->SetTrackingArray2(fArraySumKineticEnergies);
151   }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
156 void TrackingAction::PostUserTrackingAction(const G4Track* /* aTrack  */) {}
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159