Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/ParticleFluence/Calo/src/SteppingAction.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 SteppingAction.cc
 27 /// \brief Implementation of the SteppingAction class
 28 //
 29 //
 30 
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 
 34 #include "SteppingAction.hh"
 35 
 36 #include "Run.hh"
 37 
 38 #include "G4IonTable.hh"
 39 #include "G4LossTableManager.hh"
 40 #include "G4ParticleDefinition.hh"
 41 #include "G4ParticleTypes.hh"
 42 #include "G4Step.hh"
 43 #include "G4StepPoint.hh"
 44 #include "G4SystemOfUnits.hh"
 45 #include "G4TouchableHistory.hh"
 46 #include "G4Track.hh"
 47 #include "G4VPhysicalVolume.hh"
 48 #include "G4VSolid.hh"
 49 #include "G4VTouchable.hh"
 50 
 51 const std::array<G4String, SteppingAction::fkNumberScoringVolumes>
 52   SteppingAction::fkArrayScoringVolumeNames = {"downstream", "side", "upstream"};
 53 
 54 const std::array<G4String, SteppingAction::fkNumberKinematicRegions>
 55   SteppingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
 56 
 57 const std::array<G4String, SteppingAction::fkNumberParticleTypes>
 58   SteppingAction::fkArrayParticleTypeNames = {"all",      "electron",   "gamma",      "muon",
 59                                               "neutrino", "pion",       "neutron",    "proton",
 60                                               "ion",      "otherMeson", "otherBaryon"};
 61 
 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 63 
 64 G4int SteppingAction::GetIndex(const G4int iScoringVolume, const G4int iKinematicRegion,
 65                                const G4int iParticleType)
 66 {
 67   G4int index = -1;
 68   if (iScoringVolume >= 0 && iScoringVolume < fkNumberScoringVolumes && iKinematicRegion >= 0
 69       && iKinematicRegion < fkNumberKinematicRegions && iParticleType >= 0
 70       && iParticleType < fkNumberParticleTypes)
 71   {
 72     index = iScoringVolume * fkNumberKinematicRegions * fkNumberParticleTypes
 73             + iKinematicRegion * fkNumberParticleTypes + iParticleType;
 74   }
 75   if (index < 0 || index >= fkNumberCombinations) {
 76     G4cerr << "SteppingAction::GetIndex : WRONG index=" << index << "  set it to 0 !" << G4endl;
 77     index = 0;
 78   }
 79   return index;
 80 }
 81 
 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 83 
 84 SteppingAction::SteppingAction() : G4UserSteppingAction()
 85 {
 86   Initialize();
 87 }
 88 
 89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 90 
 91 void SteppingAction::Initialize()
 92 {
 93   // Initialization needed at the beginning of each Run
 94   fPrimaryParticleId = 0;
 95   fPrimaryParticleEnergy = 0.0;
 96   fPrimaryParticleDirection = G4ThreeVector(0.0, 0.0, 1.0);
 97   fAbsorberMaterialName = "";
 98   fActiveMaterialName = "";
 99   fIsFirstStepOfTheEvent = true;
100   fIsFirstStepInAbsorberLayer = true;
101   fIsFirstStepInActiveLayer = true;
102   fIsFirstStepInScoringUpDown = true;
103   fIsFirstStepInScoringSide = true;
104   fCubicVolumeScoringUpDown = 1.0;
105   fCubicVolumeScoringSide = 1.0;
106   for (G4int i = 0; i < fkNumberCombinations; ++i) {
107     fArraySumStepLengths[i] = 0.0;
108   }
109   /*
110   for ( G4int i = 0; i < fkNumberCombinations; ++i ) fArraySumStepLengths[i] = 999.9;
111   G4cout << " fkNumberCombinations=" << fkNumberCombinations << G4endl;
112   for ( G4int i = 0; i < fkNumberScoringVolumes; ++i ) {
113     for ( G4int j = 0; j < fkNumberKinematicRegions; ++j ) {
114       for ( G4int k = 0; k < fkNumberParticleTypes; ++k ) {
115         G4int index = GetIndex( i, j, k );
116         G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
117         if ( fArraySumStepLengths[ index ] < 1.0 ) G4cout << " <=== REPEATED!";
118         else                                       fArraySumStepLengths[ index ] = 0.0;
119         G4cout << G4endl;
120       }
121     }
122   }
123   for ( G4int i = 0; i < fkNumberCombinations; ++i ) {
124     if ( fArraySumStepLengths[i] > 999.0 ) G4cout << " i=" << i << " NOT COVERED !" << G4endl;
125   }
126   */
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
131 void SteppingAction::UserSteppingAction(const G4Step* theStep)
132 {
133   // Get information on the primary particle
134   if (fIsFirstStepOfTheEvent) {
135     if (theStep->GetTrack()->GetParentID() == 0) {
136       fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
137       fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
138       fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
139       if (fRunPtr) {
140         fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
141         fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
142         fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
143       }
144       fIsFirstStepOfTheEvent = false;
145     }
146   }
147   // Get information on the materials of the calorimeter
148   if (fIsFirstStepInAbsorberLayer
149       && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiAbsorber")
150   {
151     fAbsorberMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
152     if (fRunPtr) fRunPtr->SetAbsorberMaterialName(fAbsorberMaterialName);
153     fIsFirstStepInAbsorberLayer = false;
154   }
155   if (fIsFirstStepInActiveLayer
156       && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiActive")
157   {
158     fActiveMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
159     if (fRunPtr) fRunPtr->SetActiveMaterialName(fActiveMaterialName);
160     fIsFirstStepInActiveLayer = false;
161   }
162   // Get information on step lengths in the scoring volumes
163   G4int iScoringVolume = -1;
164   if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringDownstream") {
165     iScoringVolume = 0;
166     if (fIsFirstStepInScoringUpDown) {
167       fCubicVolumeScoringUpDown =
168         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
169       if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
170       fIsFirstStepInScoringUpDown = false;
171     }
172   }
173   else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringSide") {
174     iScoringVolume = 1;
175     if (fIsFirstStepInScoringSide) {
176       fCubicVolumeScoringSide =
177         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
178       if (fRunPtr) fRunPtr->SetCubicVolumeScoringSide(fCubicVolumeScoringSide);
179       fIsFirstStepInScoringSide = false;
180     }
181   }
182   else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringUpstream") {
183     iScoringVolume = 2;
184     if (fIsFirstStepInScoringUpDown) {
185       fCubicVolumeScoringUpDown =
186         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
187       if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
188       fIsFirstStepInScoringUpDown = false;
189     }
190   }
191   if (iScoringVolume >= 0) {
192     // In the case of the upstream scoring volume, consider only particles whose direction
193     // is opposite with respect to the primary particle (this is needed, in particular,
194     // for avoiding to account the incoming, primary beam particle in the "upstream" fluence).
195     if (iScoringVolume == 2
196         && fPrimaryParticleDirection.dot(theStep->GetPreStepPoint()->GetMomentumDirection()) > 0.0)
197       return;
198     G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
199     G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
200                      ? 0
201                      : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
202     /*
203     G4cout << std::setprecision(6)
204            << theStep->GetTrack()->GetDefinition()->GetParticleName() << "  absPdg=" << absPdg
205            << "  Ekin[MeV]=" << theStep->GetPreStepPoint()->GetKineticEnergy()
206            << "  (rho,z)[mm]=(" << theStep->GetTrack()->GetPosition().perp()
207            << "," << theStep->GetTrack()->GetPosition().z() << ")"
208            << "  " << theStep->GetTrack()->GetVolume()->GetName()
209            << "  " << theStep->GetTrack()->GetMaterial()->GetName()
210            << "  L[mm]=" << stepLength << "  "
211            << ( fPrimaryParticleDirection.dot(
212                   theStep->GetPreStepPoint()->GetMomentumDirection() ) > 0.0
213               ? "forward" : "backward" )
214            << G4endl;
215     */
216     // Three kinematical regions:  [0] : any value ;  [1] : below 20 MeV ;  [2] : above 20 MeV
217     G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
218     G4int iParticleType = -1;
219     if (absPdg == 11)
220       iParticleType = 1;  // electron (and positron)
221     else if (absPdg == 22)
222       iParticleType = 2;  // gamma
223     else if (absPdg == 13)
224       iParticleType = 3;  // muons (mu- and mu+)
225     else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
226       iParticleType = 4;  // neutrinos
227     // (and anti-neutrinos), all flavors
228     else if (absPdg == 111 || absPdg == 211)
229       iParticleType = 5;  // (charged) pions
230     else if (absPdg == 2112)
231       iParticleType = 6;  // neutron (and anti-neutron)
232     else if (absPdg == 2212)
233       iParticleType = 7;  // proton  (and anti-proton)
234     else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||  // ions (and anti-ions)
235              G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
236       iParticleType = 8;
237     else if (absPdg < 1000)
238       iParticleType = 9;  // other mesons (e.g. kaons) (Note: this works
239                           // in most cases, but not always!)
240     else if (absPdg > 1000)
241       iParticleType = 10;  // other baryons (e.g. hyperons, anti-hyperons,
242                            // etc.)
243     // Consider the specific case : scoring volume, kinematic region and particle type
244     G4int index = GetIndex(iScoringVolume, iKinematicRegion, iParticleType);
245     fArraySumStepLengths[index] += stepLength;
246     // Consider the "all" particle case, with the same scoring volume and kinematic region
247     index = GetIndex(iScoringVolume, iKinematicRegion, 0);
248     fArraySumStepLengths[index] += stepLength;
249     // Consider the "any" kinematic region case, with the same scoring volume and particle type
250     index = GetIndex(iScoringVolume, 0, iParticleType);
251     fArraySumStepLengths[index] += stepLength;
252     // Consider the "any" kinematic region and "all" particle, with the same scoring volume
253     index = GetIndex(iScoringVolume, 0, 0);
254     fArraySumStepLengths[index] += stepLength;
255     if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
256   }
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260