Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/ParticleFluence/Calo/include/SteppingAction.hh

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.hh
 27 /// \brief Definition of the SteppingAction class
 28 //
 29 //
 30 
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 
 34 #ifndef SteppingAction_H
 35 #define SteppingAction_H 1
 36 
 37 #include "G4ThreeVector.hh"
 38 #include "G4UserSteppingAction.hh"
 39 #include "globals.hh"
 40 
 41 #include <array>
 42 
 43 class Run;
 44 
 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 46 
 47 class SteppingAction : public G4UserSteppingAction
 48 {
 49   public:
 50     SteppingAction();
 51     ~SteppingAction() override = default;
 52 
 53     void UserSteppingAction(const G4Step*) override;
 54     // This is the main method where the step lengths of particles inside
 55     // the scoring volumes are collected, and then the corresponding fluences
 56     // are filled up in the Run object where they are stored (and then
 57     // printed out at the end of the Run).
 58     // (For simplicity and brevity, we avoid histograms and compute instead
 59     //  some statistics ourself, which will be print-out at the end of the run.)
 60 
 61     void Initialize();
 62     // This method is called by RunAction::BeginOfRunAction for the
 63     // initialization of the stepping-action at the beginning of each Run.
 64     // This is necessary because different runs can have different primary particle
 65     // types, kinetic energies, and detector configurations.
 66 
 67     void SetRunPointer(Run* inputValue = nullptr) { fRunPtr = inputValue; }
 68     // This method is called by RunAction::BeginOfRunAction for providing to the
 69     // stepping-action the pointer to the run object at the beginning of each Run.
 70     // This pointer is then used to pass the information collected by the stepping-action
 71     // to the run object.
 72 
 73     G4double GetCubicVolumeScoringUpDown() const { return fCubicVolumeScoringUpDown; }
 74     G4double GetCubicVolumeScoringSide() const { return fCubicVolumeScoringSide; }
 75     // The cubic-volumes of the scoring volumes are needed to get the fluence from
 76     // the sum of step lengths in those scoring volumes.
 77     // Notice that two of the three scoring volumes - upstream and downstream -
 78     // have the same cubic-volume, that we call "fCubicVolumeScoringUpDown".
 79 
 80     static const G4int fkNumberScoringVolumes = 3;  // downstream, side, upstream
 81     static const G4int fkNumberKinematicRegions = 3;  // all, below 20 MeV, above 20 MeV
 82     static const G4int fkNumberParticleTypes = 11;  // all, e, gamma, mu, nu, pi, n, p, ions,
 83                                                     // other-mesons, other-baryons
 84     static const G4int fkNumberCombinations =
 85       fkNumberScoringVolumes * fkNumberKinematicRegions * fkNumberParticleTypes;
 86     static const std::array<G4String, fkNumberScoringVolumes> fkArrayScoringVolumeNames;
 87     static const std::array<G4String, fkNumberKinematicRegions> fkArrayKinematicRegionNames;
 88     static const std::array<G4String, fkNumberParticleTypes> fkArrayParticleTypeNames;
 89     static G4int GetIndex(const G4int iScoringVolume, const G4int iKinematicRegion,
 90                           const G4int iParticleType);
 91 
 92   private:
 93     Run* fRunPtr;  // Pointer to the Run object
 94     G4int fPrimaryParticleId;
 95     G4double fPrimaryParticleEnergy;
 96     G4ThreeVector fPrimaryParticleDirection;
 97     G4String fAbsorberMaterialName;
 98     G4String fActiveMaterialName;
 99     G4bool fIsFirstStepOfTheEvent;
100     G4bool fIsFirstStepInAbsorberLayer;
101     G4bool fIsFirstStepInActiveLayer;
102     G4bool fIsFirstStepInScoringUpDown;
103     G4bool fIsFirstStepInScoringSide;
104     G4double fCubicVolumeScoringUpDown;
105     G4double fCubicVolumeScoringSide;
106 
107     std::array<G4double, fkNumberCombinations> fArraySumStepLengths;
108     // Array to collect the sum of step lengths in the scoring volumes for the whole run,
109     // according to the various cases (kinematical region and particle type).
110     // Note that the fluence in a scoring volume is defined as sum of step lengths
111     // in that scoring volume divided by the cubic-volume of that scoring volume.
112 };
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
116 #endif
117