|
|
|||
File indexing completed on 2026-05-19 07:40:39
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 /// \file SteppingAction.cc 0027 /// \brief Implementation of the SteppingAction class 0028 0029 // This example is provided by the Geant4-DNA collaboration 0030 // Any report or published results obtained using the Geant4-DNA software 0031 // shall cite the following Geant4-DNA collaboration publications: 0032 // Med. Phys. 45 (2018) e722-e739 0033 // Phys. Med. 31 (2015) 861-874 0034 // Med. Phys. 37 (2010) 4692-4708 0035 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 0036 // 0037 // The Geant4-DNA web site is available at http://geant4-dna.org 0038 // 0039 0040 #include "SteppingAction.hh" 0041 0042 #include "HistoManager.hh" 0043 #include "Run.hh" 0044 0045 #include "G4RunManager.hh" 0046 0047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0048 0049 SteppingAction::SteppingAction() : G4UserSteppingAction() {} 0050 0051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0052 0053 SteppingAction::~SteppingAction() {} 0054 0055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0056 0057 void SteppingAction::UserSteppingAction(const G4Step* aStep) 0058 { 0059 const G4StepPoint* endPoint = aStep->GetPostStepPoint(); 0060 G4String procName = endPoint->GetProcessDefinedStep()->GetProcessName(); 0061 0062 Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun()); 0063 0064 G4bool transmit = (endPoint->GetStepStatus() <= fGeomBoundary); 0065 if (transmit) { 0066 run->CountProcesses(procName); 0067 } 0068 else { 0069 // Count processes and sum track length 0070 G4double stepLength = aStep->GetStepLength(); 0071 run->CountProcesses(procName); 0072 run->SumTrack(stepLength); 0073 } 0074 0075 // Call analysis manager 0076 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 0077 0078 // Scattered primary particle 0079 G4int id = 1; 0080 if (aStep->GetTrack()->GetTrackStatus() == fAlive) { 0081 G4double energy = endPoint->GetKineticEnergy(); 0082 analysisManager->FillH1(id, energy); 0083 id = 2; 0084 G4ThreeVector direction = endPoint->GetMomentumDirection(); 0085 // Assuming incident particle shot along x 0086 G4double costeta = direction.x(); 0087 analysisManager->FillH1(id, costeta); 0088 } 0089 0090 // Secondaries 0091 const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep(); 0092 for (size_t lp = 0; lp < (*secondary).size(); lp++) { 0093 G4double charge = (*secondary)[lp]->GetDefinition()->GetPDGCharge(); 0094 if (charge != 0.) { 0095 id = 3; 0096 } 0097 else { 0098 id = 5; 0099 } 0100 G4double energy = (*secondary)[lp]->GetKineticEnergy(); 0101 analysisManager->FillH1(id, energy); 0102 0103 ++id; 0104 G4ThreeVector direction = (*secondary)[lp]->GetMomentumDirection(); 0105 G4double costeta = direction.x(); 0106 analysisManager->FillH1(id, costeta); 0107 0108 // Energy tranferred to charged secondaries 0109 if (charge != 0.) { 0110 run->SumeTransf(energy); 0111 } 0112 } 0113 0114 // Kill event after first interaction 0115 G4RunManager::GetRunManager()->AbortEvent(); 0116 }
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|