File indexing completed on 2026-04-05 07:50:35
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 #include "SteppingAction.hh"
0030
0031 #include "DetectorConstruction.hh"
0032 #include "HistoManager.hh"
0033 #include "Run.hh"
0034 #include "TrackingAction.hh"
0035
0036 #include "G4Gamma.hh"
0037 #include "G4RunManager.hh"
0038 #include "G4SteppingManager.hh"
0039 #include "G4UnitsTable.hh"
0040
0041
0042
0043 SteppingAction::SteppingAction(DetectorConstruction* det, TrackingAction* TrAct)
0044 : fDetector(det), fTrackAction(TrAct), fWall(0), fCavity(0)
0045 {
0046 first = true;
0047 fTrackSegm = 0.;
0048 fDirectionIn = G4ThreeVector(0., 0., 0.);
0049 }
0050
0051
0052
0053 SteppingAction::~SteppingAction() {}
0054
0055
0056
0057 void SteppingAction::UserSteppingAction(const G4Step* step)
0058 {
0059
0060 if (first) {
0061 fWall = fDetector->GetWall();
0062 fCavity = fDetector->GetCavity();
0063 first = false;
0064 }
0065
0066 Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun());
0067
0068
0069 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0070
0071
0072
0073 G4StepPoint* point1 = step->GetPreStepPoint();
0074 G4VPhysicalVolume* volume = point1->GetTouchableHandle()->GetVolume();
0075
0076
0077
0078 G4StepPoint* point2 = step->GetPostStepPoint();
0079 const G4VProcess* process = point2->GetProcessDefinedStep();
0080 if (process) run->CountProcesses(process->GetProcessName());
0081
0082
0083
0084 if (volume == fCavity) {
0085 G4double edep = step->GetTotalEnergyDeposit();
0086 if (edep > 0.) fTrackAction->AddEdepCavity(edep);
0087 }
0088
0089
0090
0091 if (step->GetTrack()->GetDefinition() == G4Gamma::Gamma()) return;
0092
0093
0094
0095 G4int id;
0096 G4double steplen = step->GetStepLength();
0097 if (volume == fWall) {
0098 run->StepInWall(steplen);
0099 id = 9;
0100 }
0101 else {
0102 run->StepInCavity(steplen);
0103 id = 10;
0104 }
0105 analysisManager->FillH1(id, steplen);
0106
0107
0108
0109 if ((volume == fWall) && (point2->GetStepStatus() == fGeomBoundary)) {
0110 fDirectionIn = point1->GetMomentumDirection();
0111 }
0112
0113
0114
0115 if (volume == fWall) return;
0116
0117 G4double ekin1 = point1->GetKineticEnergy();
0118 G4double ekin2 = point2->GetKineticEnergy();
0119
0120
0121
0122 if (point1->GetStepStatus() == fGeomBoundary) {
0123 fTrackSegm = 0.;
0124 G4ThreeVector vertex = step->GetTrack()->GetVertexPosition();
0125 analysisManager->FillH1(4, vertex.z());
0126 run->FlowInCavity(0, ekin1);
0127 analysisManager->FillH1(5, ekin1);
0128 if (steplen > 0.) {
0129 G4ThreeVector directionOut = (point2->GetPosition() - point1->GetPosition()).unit();
0130 G4ThreeVector normal =
0131 point1->GetTouchableHandle()->GetSolid()->SurfaceNormal(point1->GetPosition());
0132 analysisManager->FillH1(6, std::acos(-fDirectionIn * normal));
0133 analysisManager->FillH1(7, std::acos(-directionOut * normal));
0134 }
0135 }
0136
0137
0138
0139 if (step->GetTrack()->GetCurrentStepNumber() == 1) fTrackSegm = 0.;
0140 fTrackSegm += steplen;
0141 if (ekin2 <= 0.) {
0142 run->AddTrakCavity(fTrackSegm);
0143 analysisManager->FillH1(8, fTrackSegm);
0144 }
0145
0146
0147
0148 if (point2->GetStepStatus() == fGeomBoundary) {
0149 run->FlowInCavity(1, ekin2);
0150 run->AddTrakCavity(fTrackSegm);
0151 analysisManager->FillH1(8, fTrackSegm);
0152 }
0153 }
0154
0155