File indexing completed on 2025-02-23 09:21:52
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027 #include "SteppingAction.hh"
0028
0029 #include "DetectorConstruction.hh"
0030
0031 #include "G4AnalysisManager.hh"
0032 #include "G4DNAElastic.hh"
0033 #include "G4DNAElectronSolvation.hh"
0034 #include "G4Event.hh"
0035 #include "G4EventManager.hh"
0036 #include "G4ITTrackHolder.hh"
0037 #include "G4Molecule.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "G4Track.hh"
0040 #include "globals.hh"
0041
0042 #include <map>
0043
0044 using MapOfDelayedLists = std::map<double, std::map<int, G4TrackList*>>;
0045
0046
0047
0048 SteppingAction::SteppingAction(DetectorConstruction* fpDet)
0049 : G4UserSteppingAction(), fpDetector(fpDet)
0050 {}
0051
0052
0053
0054 SteppingAction::~SteppingAction() {}
0055
0056
0057
0058 void SteppingAction::UserSteppingAction(const G4Step* step)
0059 {
0060 SetupFlags(step);
0061
0062 if (fVolumeType == DNAVolumeType::physWorld) {
0063 step->GetTrack()->SetTrackStatus(fStopAndKill);
0064 }
0065
0066 if (fVolumeType == DNAVolumeType::VoxelStraight) {
0067 return;
0068 }
0069
0070 G4double dE = step->GetTotalEnergyDeposit() / eV;
0071 const G4VProcess* pProcess = step->GetPostStepPoint()->GetProcessDefinedStep();
0072
0073
0074 if ((0 > dE) || (nullptr != dynamic_cast<const G4DNAElectronSolvation*>(pProcess))) {
0075 return;
0076 }
0077
0078 G4double x = step->GetPreStepPoint()->GetPosition().x() / nanometer;
0079 G4double y = step->GetPreStepPoint()->GetPosition().y() / nanometer;
0080 G4double z = step->GetPreStepPoint()->GetPosition().z() / nanometer;
0081
0082 G4double copyNo = G4double(step->GetPreStepPoint()->GetTouchable()->GetVolume()->GetCopyNo());
0083
0084 G4double eventId =
0085 G4double(G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID());
0086
0087 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0088 analysisManager->FillNtupleDColumn(1, 0, x);
0089 analysisManager->FillNtupleDColumn(1, 1, y);
0090 analysisManager->FillNtupleDColumn(1, 2, z);
0091 analysisManager->FillNtupleDColumn(1, 3, dE);
0092 analysisManager->FillNtupleDColumn(
0093 1, 4,
0094 (step->GetPreStepPoint()->GetKineticEnergy() - step->GetPostStepPoint()->GetKineticEnergy())
0095 / eV);
0096 analysisManager->FillNtupleIColumn(1, 5, G4int(fVolumeType));
0097 analysisManager->FillNtupleDColumn(1, 6, copyNo);
0098 analysisManager->FillNtupleIColumn(1, 7, eventId);
0099 analysisManager->AddNtupleRow(1);
0100
0101 MapOfDelayedLists delayList = G4ITTrackHolder::Instance()->GetDelayedLists();
0102 for (auto& delayedmap_it : delayList) {
0103 for (auto& trackList : delayedmap_it.second) {
0104 if (nullptr == trackList.second) {
0105 return;
0106 }
0107 G4TrackList::iterator itt = trackList.second->begin();
0108 G4TrackList::iterator endd = trackList.second->end();
0109 for (; itt != endd; ++itt) {
0110 G4Track* track = *itt;
0111 if (track->GetParentID() != step->GetTrack()->GetTrackID()
0112 || track->GetPosition() != step->GetPostStepPoint()->GetPosition())
0113 {
0114 return;
0115 }
0116
0117 track->SetTrackStatus(fStopAndKill);
0118 }
0119 }
0120 }
0121 }
0122
0123
0124 void SteppingAction::SetupFlags(const G4Step* step)
0125 {
0126 fVolumeType = SetupVolumeType(step->GetPreStepPoint()->GetPhysicalVolume());
0127 }
0128
0129
0130
0131 DNAVolumeType SteppingAction::SetupVolumeType(const G4VPhysicalVolume* pPhyVolume)
0132 {
0133 auto it = (fpDetector->GetGeoDataMap()).find(pPhyVolume);
0134
0135 if (it == (fpDetector->GetGeoDataMap()).end()) {
0136 G4ExceptionDescription exceptionDescription;
0137 exceptionDescription << pPhyVolume->GetName() << " is wrong physical volume";
0138 G4Exception("SteppingAction::SetupVolumeFlag()", "SteppingAction01", FatalException,
0139 exceptionDescription);
0140 }
0141
0142 return it->second;
0143 }
0144