G4OCCT 0.1.0
Geant4 interface to Open CASCADE Technology (OCCT) geometry definitions
Loading...
Searching...
No Matches
G4OCCTSteppingAction.cc
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-2.1-or-later
2// Copyright (C) 2026 G4OCCT Contributors
3
5
8#include "G4OCCTRunAction.hh"
10
11#include <G4AnalysisManager.hh>
12#include <G4Step.hh>
13#include <G4SystemOfUnits.hh>
14#include <G4Track.hh>
15#include <G4VPhysicalVolume.hh>
16
18 G4OCCTTrackingAction* trackingAction,
19 const G4OCCTRunAction* runAction)
20 : fEventAction(eventAction), fTrackingAction(trackingAction), fRunAction(runAction) {}
21
23 const G4double edep = step->GetTotalEnergyDeposit();
24 const G4double stepLength = step->GetStepLength();
25
26 fEventAction->AddStep(edep, stepLength);
27 fTrackingAction->AddEdepToCurrentTrack(edep);
28
29 const G4int ntId = fRunAction->GetStepsNtupleId();
30 if (ntId < 0)
31 return;
32
33 const G4Track* track = step->GetTrack();
34 const G4StepPoint* pre = step->GetPreStepPoint();
35 const G4ThreeVector& pos = pre->GetPosition();
36 const G4VPhysicalVolume* vol = pre->GetPhysicalVolume();
37 const G4String volName = vol ? vol->GetLogicalVolume()->GetName() : "OutOfWorld";
38 const G4int eventId = fEventAction->GetEventId();
39
40 auto* am = G4AnalysisManager::Instance();
41 G4int col = 0;
42 am->FillNtupleIColumn(ntId, col++, eventId);
43 am->FillNtupleIColumn(ntId, col++, track->GetTrackID());
44 am->FillNtupleIColumn(ntId, col++, track->GetCurrentStepNumber());
45 am->FillNtupleSColumn(ntId, col++, volName);
46 am->FillNtupleSColumn(ntId, col++, track->GetDefinition()->GetParticleName());
47 am->FillNtupleDColumn(ntId, col++, edep / MeV);
48 am->FillNtupleDColumn(ntId, col++, stepLength / mm);
49 am->FillNtupleDColumn(ntId, col++, pos.x() / mm);
50 am->FillNtupleDColumn(ntId, col++, pos.y() / mm);
51 am->FillNtupleDColumn(ntId, col++, pos.z() / mm);
52 am->AddNtupleRow(ntId);
53}
Event action that accumulates per-event totals and fills the events ntuple.
Run action that manages CSV output via G4AnalysisManager.
Stepping action that fills the per-step ntuple.
Tracking action that fills the per-track ntuple.
Event action that accumulates per-event quantities.
G4int GetEventId() const
Returns the current event ID; cached in BeginOfEventAction.
void AddStep(G4double edep, G4double stepLength)
Called from G4OCCTSteppingAction for each step.
Run action that opens and closes a G4AnalysisManager CSV output file.
G4int GetStepsNtupleId() const
G4OCCTSteppingAction(G4OCCTEventAction *eventAction, G4OCCTTrackingAction *trackingAction, const G4OCCTRunAction *runAction)
void UserSteppingAction(const G4Step *step) override
Tracking action that records one row per G4Track into the "tracks" ntuple.
void AddEdepToCurrentTrack(G4double edep)
Called from G4OCCTSteppingAction to accumulate energy deposit for the current track.