File indexing completed on 2025-11-04 09:29:51
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
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039 #include "SteppingAction.hh"
0040 #include "SteppingMessenger.hh"
0041
0042 #include "DetectorConstruction.hh"
0043 #include "PrimaryGeneratorAction.hh"
0044 #include "RunAction.hh"
0045
0046 #include "G4Alpha.hh"
0047 #include "G4AnalysisManager.hh"
0048 #include "G4DNAGenericIonsManager.hh"
0049 #include "G4Electron.hh"
0050 #include "G4Event.hh"
0051 #include "G4EventManager.hh"
0052 #include "G4Gamma.hh"
0053 #include "G4Proton.hh"
0054 #include "G4SteppingManager.hh"
0055 #include "G4SystemOfUnits.hh"
0056
0057
0058
0059 SteppingAction::SteppingAction() : G4UserSteppingAction()
0060 {
0061 fSteppingMessenger = new SteppingMessenger(this);
0062 }
0063
0064
0065
0066 SteppingAction::~SteppingAction()
0067 {
0068 delete fSteppingMessenger;
0069 }
0070
0071
0072
0073 void SteppingAction::UserSteppingAction(const G4Step* step)
0074 {
0075
0076 if (!step->GetPostStepPoint()) return;
0077 if (!step->GetPostStepPoint()->GetProcessDefinedStep()) return;
0078
0079
0080 if (step->GetTrack()->GetCurrentStepNumber()>=1 && fKill) step->GetTrack()->SetTrackStatus(fKillTrackAndSecondaries);
0081
0082
0083
0084 G4double flagParticle = -1.;
0085 G4double flagProcess = -1.;
0086 G4double x, y, z, xp, yp, zp;
0087
0088
0089
0090
0091
0092 G4ParticleDefinition* partDef = step->GetTrack()->GetDynamicParticle()->GetDefinition();
0093
0094 if (partDef == G4Gamma::GammaDefinition()) flagParticle = 0;
0095
0096 if (partDef == G4Electron::ElectronDefinition()) flagParticle = 1;
0097
0098 if (partDef == G4Proton::ProtonDefinition()) flagParticle = 2;
0099
0100 if (partDef == G4Alpha::AlphaDefinition()) flagParticle = 4;
0101
0102 G4DNAGenericIonsManager* instance;
0103 instance = G4DNAGenericIonsManager::Instance();
0104
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114 if (partDef == instance->GetIon("hydrogen")) flagParticle = 3;
0115
0116 if (partDef == instance->GetIon("alpha+")) flagParticle = 5;
0117
0118 if (partDef == instance->GetIon("helium")) flagParticle = 6;
0119
0120
0121 if (partDef->GetPDGCharge()>4) flagParticle = 7;
0122
0123
0124
0125
0126
0127
0128
0129
0130
0131
0132
0133
0134
0135
0136
0137
0138
0139
0140
0141
0142
0143 G4StepPoint* preStep = step->GetPreStepPoint();
0144 G4StepPoint* postStep = step->GetPostStepPoint();
0145 G4int procID = postStep->GetProcessDefinedStep()->GetProcessSubType();
0146
0147 const G4String& processName = postStep->GetProcessDefinedStep()->GetProcessName();
0148
0149
0150 if (flagParticle == 0) {
0151 if (procID == 12)
0152 flagProcess = 81;
0153 else if (procID == 13)
0154 flagProcess = 82;
0155 else if (procID == 14)
0156 flagProcess = 83;
0157 else if (procID == 11)
0158 flagProcess = 84;
0159 }
0160
0161
0162 else if (flagParticle == 1) {
0163 if (procID == 58)
0164 flagProcess = 10;
0165 else if (procID == 51)
0166 flagProcess = 11;
0167 else if (procID == 52)
0168 flagProcess = 12;
0169 else if (procID == 53)
0170 flagProcess = 13;
0171 else if (procID == 55)
0172 flagProcess = 14;
0173 else if (procID == 54)
0174 flagProcess = 15;
0175 else if (procID == 10)
0176 flagProcess = 110;
0177 else if (procID == 1)
0178 flagProcess = 120;
0179 else if (procID == 2)
0180 flagProcess = 130;
0181 }
0182
0183
0184 else if (flagParticle == 2) {
0185 if (procID == 51)
0186 flagProcess = 21;
0187 else if (procID == 52)
0188 flagProcess = 22;
0189 else if (procID == 53)
0190 flagProcess = 23;
0191 else if (procID == 56)
0192 flagProcess = 24;
0193 else if (procID == 10)
0194 flagProcess = 210;
0195 else if (procID == 1)
0196 flagProcess = 220;
0197 else if (procID == 2)
0198 flagProcess = 230;
0199 else if (procID == 8)
0200 flagProcess = 240;
0201 }
0202
0203
0204 else if (flagParticle == 3) {
0205 if (procID == 51)
0206 flagProcess = 31;
0207 else if (procID == 52)
0208 flagProcess = 32;
0209 else if (procID == 53)
0210 flagProcess = 33;
0211 else if (procID == 57)
0212 flagProcess = 35;
0213 }
0214
0215
0216 else if (flagParticle == 4) {
0217 if (procID == 51)
0218 flagProcess = 41;
0219 else if (procID == 52)
0220 flagProcess = 42;
0221 else if (procID == 53)
0222 flagProcess = 43;
0223 else if (procID == 56)
0224 flagProcess = 44;
0225 else if (procID == 10)
0226 flagProcess = 410;
0227 else if (procID == 1)
0228 flagProcess = 420;
0229 else if (procID == 2)
0230 flagProcess = 430;
0231 else if (procID == 8)
0232 flagProcess = 440;
0233 }
0234
0235
0236 else if (flagParticle == 5) {
0237 if (procID == 51)
0238 flagProcess = 51;
0239 else if (procID == 52)
0240 flagProcess = 52;
0241 else if (procID == 53)
0242 flagProcess = 53;
0243 else if (procID == 56)
0244 flagProcess = 54;
0245 else if (procID == 57)
0246 flagProcess = 55;
0247 else if (procID == 10)
0248 flagProcess = 510;
0249 else if (procID == 1)
0250 flagProcess = 520;
0251 else if (procID == 2)
0252 flagProcess = 530;
0253 else if (procID == 8)
0254 flagProcess = 540;
0255 }
0256
0257
0258 else if (flagParticle == 6) {
0259 if (procID == 51)
0260 flagProcess = 61;
0261 else if (procID == 52)
0262 flagProcess = 62;
0263 else if (procID == 53)
0264 flagProcess = 63;
0265 else if (procID == 57)
0266 flagProcess = 65;
0267 }
0268
0269 else if (flagParticle == 7) {
0270
0271 if (processName == "Capture")
0272 flagProcess = 1;
0273
0274
0275
0276 else if (procID == 210)
0277 flagProcess = 2;
0278
0279 else if (procID == 53)
0280 flagProcess = 73;
0281 else if (procID == 10)
0282 flagProcess = 710;
0283 else if (procID == 1)
0284 flagProcess = 720;
0285 else if (procID == 2)
0286 flagProcess = 730;
0287 else if (procID == 8)
0288 flagProcess = 740;
0289 }
0290
0291
0292
0293
0294
0295
0296
0297
0298
0299
0300
0301
0302
0303
0304
0305
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328
0329
0330
0331 if (processName != "Transportation") {
0332 x = preStep->GetPosition().x() / nanometer;
0333 y = preStep->GetPosition().y() / nanometer;
0334 z = preStep->GetPosition().z() / nanometer;
0335
0336 xp = postStep->GetPosition().x() / nanometer;
0337 yp = postStep->GetPosition().y() / nanometer;
0338 zp = postStep->GetPosition().z() / nanometer;
0339
0340 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0341
0342 analysisManager->FillNtupleDColumn(0, flagParticle);
0343 analysisManager->FillNtupleDColumn(1, flagProcess);
0344 analysisManager->FillNtupleDColumn(2, xp);
0345 analysisManager->FillNtupleDColumn(3, yp);
0346 analysisManager->FillNtupleDColumn(4, zp);
0347 analysisManager->FillNtupleDColumn(5, step->GetTotalEnergyDeposit() / eV);
0348
0349 analysisManager->FillNtupleDColumn(
0350 6, std::sqrt((x - xp) * (x - xp) + (y - yp) * (y - yp) + (z - zp) * (z - zp)));
0351
0352 analysisManager->FillNtupleDColumn(
0353 7, (preStep->GetKineticEnergy() - postStep->GetKineticEnergy()) / eV);
0354
0355 analysisManager->FillNtupleDColumn(8, preStep->GetKineticEnergy() / eV);
0356
0357 analysisManager->FillNtupleDColumn(9, preStep->GetMomentumDirection()
0358 * postStep->GetMomentumDirection());
0359
0360 analysisManager->FillNtupleIColumn(
0361 10, G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID());
0362
0363 analysisManager->FillNtupleIColumn(11, step->GetTrack()->GetTrackID());
0364
0365 analysisManager->FillNtupleIColumn(12, step->GetTrack()->GetParentID());
0366
0367 analysisManager->FillNtupleIColumn(13, step->GetTrack()->GetCurrentStepNumber());
0368
0369 analysisManager->AddNtupleRow();
0370 }
0371 }