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