Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/electromagnetic/TestEm16/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 electromagnetic/TestEm16/src/SteppingAction.cc
 27 /// \brief Implementation of the SteppingAction class
 28 //
 29 //
 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 
 33 #include "SteppingAction.hh"
 34 
 35 #include "HistoManager.hh"
 36 #include "Run.hh"
 37 
 38 #include "G4EmProcessSubType.hh"
 39 #include "G4ParticleTypes.hh"
 40 #include "G4RunManager.hh"
 41 #include "G4SteppingManager.hh"
 42 #include "G4SystemOfUnits.hh"
 43 #include "G4UnitsTable.hh"
 44 #include "G4VProcess.hh"
 45 
 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 47 
 48 SteppingAction::SteppingAction() : G4UserSteppingAction() {}
 49 
 50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 51 
 52 void SteppingAction::UserSteppingAction(const G4Step* aStep)
 53 {
 54   Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun());
 55 
 56   const G4VProcess* process = aStep->GetPostStepPoint()->GetProcessDefinedStep();
 57   if (process == nullptr) return;
 58 
 59   G4int eventID = G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID();
 60   G4Track* trk = aStep->GetTrack();
 61 
 62   thread_local G4int iCalled = 0;
 63   const G4int nprint = 10;  // set to 10 to get debug print for first 10 calls
 64   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
 65 
 66   if (fSynchrotronRadiation == process->GetProcessSubType()) {
 67     const G4StepPoint* PrePoint = aStep->GetPreStepPoint();
 68     const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
 69     size_t lp = (*secondary).size();
 70     if (lp) {
 71       G4double Egamma = (*secondary)[lp - 1]->GetTotalEnergy();
 72       run->f_n_gam_sync++;
 73       run->f_e_gam_sync += Egamma;
 74       run->f_e_gam_sync2 += Egamma * Egamma;
 75       if (Egamma > run->f_e_gam_sync_max) run->f_e_gam_sync_max = Egamma;
 76       run->f_lam_gam_sync += aStep->GetStepLength();
 77       if (iCalled < nprint) {
 78         G4double Eelec = PrePoint->GetTotalEnergy();
 79         G4ThreeVector Pelec = PrePoint->GetMomentum();
 80         G4ThreeVector PGamma = (*secondary)[lp - 1]->GetMomentum();
 81         G4bool IsGamma = ((*secondary)[lp - 1]->GetDefinition() == G4Gamma::Gamma());
 82         const G4int wid = 8;
 83         if (analysisManager->GetVerboseLevel() > 1 && iCalled < nprint)
 84           G4cout << __FILE__ << " line " << __LINE__ << " " << __FUNCTION__
 85                  << " processName=" << process->GetProcessName()
 86                  << " Step Length=" << std::setw(wid)
 87                  << G4BestUnit(aStep->GetStepLength(), "Length") << " Eelec=" << std::setw(wid)
 88                  << G4BestUnit(Eelec, "Energy") << " Pelec=" << G4BestUnit(Pelec, "Energy")
 89                  << " IsGamma=" << IsGamma << " Egamma=" << std::setw(8)
 90                  << G4BestUnit(Egamma, "Energy") << " PGamma=" << std::setw(8)
 91                  << G4BestUnit(PGamma, "Energy") << " #secondaries lp=" << lp << '\n';
 92       }
 93       analysisManager->FillH1(1, Egamma);
 94       analysisManager->FillH1(2, Egamma,
 95                               Egamma / keV);  // power spectrum : gamma weighted with its energy
 96       analysisManager->FillH1(3, aStep->GetStepLength());
 97     }
 98   }
 99 
100   thread_local G4ThreeVector IncomingPhotonDirection;
101   if (fGammaReflection == process->GetProcessSubType()) {
102     const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
103     size_t lp = (*secondary).size();
104     if (lp) {
105       G4double Egamma = (*secondary)[lp - 1]->GetTotalEnergy();
106       ++run->f_n_Xray_Refl;
107       analysisManager->FillH1(1, Egamma, 1. / run->GetNumberOfEventToBeProcessed());
108       if (analysisManager->GetVerboseLevel() > 1 && run->f_n_Xray_Refl < nprint)
109         G4cout << __FILE__ << " line " << __LINE__ << " " << __FUNCTION__
110                << " iCalled=" << std::setw(3) << iCalled << " eventID=" << std::setw(3) << eventID
111                << " f_n_Xray_Refl=" << run->f_n_Xray_Refl
112                << " ProcessName=" << process->GetProcessName()
113                << " direction=" << trk->GetMomentumDirection() << " Egamma=" << Egamma / keV
114                << " keV" << G4endl;
115       IncomingPhotonDirection = trk->GetMomentumDirection();
116     }
117   }
118 
119   const G4VProcess* creator = trk->GetCreatorProcess();
120   G4int CreatorProcessSubType = 0;
121   if (creator != nullptr) CreatorProcessSubType = creator->GetProcessSubType();
122   if (CreatorProcessSubType == fGammaReflection) {  // transportation after XrayReflection
123     if (IncomingPhotonDirection != G4ThreeVector(0, 0, 0)) {
124       G4double cos_angle =
125         IncomingPhotonDirection * trk->GetMomentumDirection();  // incoming * outgoing direction
126       G4double IncidentAngle = std::acos(cos_angle) / 2.;
127       analysisManager->FillH1(4, IncidentAngle, 1. / run->GetNumberOfEventToBeProcessed());
128       if (analysisManager->GetVerboseLevel() > 1 && iCalled < nprint)
129         G4cout << __FILE__ << " line " << __LINE__ << " " << __FUNCTION__
130                << " iCalled=" << std::setw(3) << iCalled << " eventID=" << std::setw(3) << eventID
131                << " CreatorProcname=" << creator->GetProcessName()
132                << " ProcessName=" << process->GetProcessName()
133                << " IncomingPhotonDirection=" << IncomingPhotonDirection
134                << "  outgoing PhotonDirection=" << trk->GetMomentumDirection()
135                << " IncidentAngle=" << IncidentAngle
136                << " NumberOfEventToBeProcessed=" << run->GetNumberOfEventToBeProcessed()
137                << " RunID=" << run->GetRunID() << G4endl;
138     }
139   }
140 
141   ++iCalled;
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145