Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/ParticleFluence/Sphere/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::fkNumberKinematicRegions>
 52   SteppingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
 53 
 54 const std::array<G4String, SteppingAction::fkNumberScoringPositions>
 55   SteppingAction::fkArrayScoringPositionNames = {"forward", "backward"};
 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 iKinematicRegion, const G4int iScoringPosition,
 65                                const G4int iParticleType)
 66 {
 67   G4int index = -1;
 68   if (iKinematicRegion >= 0 && iKinematicRegion < fkNumberKinematicRegions && iScoringPosition >= 0
 69       && iScoringPosition < fkNumberScoringPositions && iParticleType >= 0
 70       && iParticleType < fkNumberParticleTypes)
 71   {
 72     index = iKinematicRegion * fkNumberScoringPositions * fkNumberParticleTypes
 73             + iScoringPosition * 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   fTargetMaterialName = "";
 98   fIsFirstStepOfTheEvent = true;
 99   fIsFirstStepInTarget = true;
100   fIsFirstStepInScoringShell = true;
101   fCubicVolumeScoringShell = 1.0;
102   for (G4int i = 0; i < fkNumberCombinations; ++i) {
103     fArraySumStepLengths[i] = 0.0;
104   }
105   /*
106   for ( G4int i = 0; i < fkNumberCombinations; ++i ) fArraySumStepLengths[i] = 999.9;
107   G4cout << " fkNumberCombinations=" << fkNumberCombinations << G4endl;
108   for ( G4int i = 0; i < fkNumberKinematicRegions; ++i ) {
109     for ( G4int j = 0; j < fkNumberScoringPositions; ++j ) {
110       for ( G4int k = 0; k < fkNumberParticleTypes; ++k ) {
111         G4int index = GetIndex( i, j, k );
112         G4cout << "(i, j, k)=(" << i << ", " << j << ", " << k << ")  ->" << index;
113         if ( fArraySumStepLengths[ index ] < 1.0 ) G4cout << " <=== REPEATED!";
114         else                                       fArraySumStepLengths[ index ] = 0.0;
115         G4cout << G4endl;
116       }
117     }
118   }
119   for ( G4int i = 0; i < fkNumberCombinations; ++i ) {
120     if ( fArraySumStepLengths[i] > 999.0 ) G4cout << " i=" << i << " NOT COVERED !" << G4endl;
121   }
122   */
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
127 void SteppingAction::UserSteppingAction(const G4Step* theStep)
128 {
129   // Get information on the primary particle
130   if (fIsFirstStepOfTheEvent) {
131     if (theStep->GetTrack()->GetParentID() == 0) {
132       fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
133       fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
134       fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
135       if (fRunPtr) {
136         fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
137         fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
138         fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
139       }
140       fIsFirstStepOfTheEvent = false;
141     }
142   }
143   // Get information on the target material
144   if (fIsFirstStepInTarget
145       && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiSphere")
146   {
147     fTargetMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
148     if (fRunPtr) fRunPtr->SetTargetMaterialName(fTargetMaterialName);
149     fIsFirstStepInTarget = false;
150   }
151   // Get information on step lengths in the scoring shell
152   if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringShell") {
153     if (fIsFirstStepInScoringShell) {
154       fCubicVolumeScoringShell =
155         theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
156       if (fRunPtr) fRunPtr->SetCubicVolumeScoringShell(fCubicVolumeScoringShell);
157       fIsFirstStepInScoringShell = false;
158     }
159     G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
160     G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
161                      ? 0
162                      : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
163     /*
164     G4cout << theStep->GetTrack()->GetDefinition()->GetParticleName() << "  absPdg=" << absPdg
165            << "  Ekin[MeV]=" << theStep->GetPreStepPoint()->GetKineticEnergy()
166            << "  r[mm]=" << theStep->GetTrack()->GetPosition().mag()
167            << "  z[mm]=" << theStep->GetTrack()->GetPosition().z()
168            << "  " << theStep->GetTrack()->GetVolume()->GetName()
169            << "  " << theStep->GetTrack()->GetMaterial()->GetName()
170            << "  L[mm]=" << stepLength << "  "
171            << ( fPrimaryParticleDirection.dot( theStep->GetTrack()->GetPosition().unit() ) > 0.0
172                 ? "forward" : "backward" ) << G4endl;
173     */
174     // Three kinematical regions:  [0] : any value ;  [1] : below 20 MeV ;  [2] : above 20 MeV
175     G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
176     // Two scoring positions:  [0] : forward hemisphere ;  [1] : backward hemisphere
177     // (with respect to the primary particle initial direction)
178     G4int iScoringPosition =
179       fPrimaryParticleDirection.dot(theStep->GetTrack()->GetPosition().unit()) > 0.0 ? 0 : 1;
180     G4int iParticleType = -1;
181     if (absPdg == 11)
182       iParticleType = 1;  // electron (and positron)
183     else if (absPdg == 22)
184       iParticleType = 2;  // gamma
185     else if (absPdg == 13)
186       iParticleType = 3;  // muons (mu- and mu+)
187     else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
188       iParticleType = 4;  // neutrinos (and
189                           // anti-neutrinos), all flavors
190     else if (absPdg == 111 || absPdg == 211)
191       iParticleType = 5;  // (charged) pions
192     else if (absPdg == 2112)
193       iParticleType = 6;  // neutron (and anti-neutron)
194     else if (absPdg == 2212)
195       iParticleType = 7;  // proton  (and anti-proton)
196     else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||  // ions (and anti-ions)
197              G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
198       iParticleType = 8;
199     else if (absPdg < 1000)
200       iParticleType = 9;  // other mesons (e.g. kaons) (Note: this works
201                           // in most cases, but not always!)
202     else if (absPdg > 1000)
203       iParticleType = 10;  // other baryons (e.g. hyperons, anti-hyperons,
204                            // etc.)
205     // Consider the specific case : kinematic region, scoring position, and particle type
206     G4int index = GetIndex(iKinematicRegion, iScoringPosition, iParticleType);
207     fArraySumStepLengths[index] += stepLength;
208     // Consider the "all" particle case, with the same kinematic region and scoring position
209     index = GetIndex(iKinematicRegion, iScoringPosition, 0);
210     fArraySumStepLengths[index] += stepLength;
211     // Consider the "any" kinematic region case, with the same scoring position and particle type
212     index = GetIndex(0, iScoringPosition, iParticleType);
213     fArraySumStepLengths[index] += stepLength;
214     // Consider the "any" kinematic region and "all" particle, with the same scoring position
215     index = GetIndex(0, iScoringPosition, 0);
216     fArraySumStepLengths[index] += stepLength;
217     if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
218   }
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222