Loading...
Searching...
No Matches
F04SteppingAction.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27/// \file field/field04/src/F04SteppingAction.cc
28/// \brief Implementation of the F04SteppingAction class
29//
30
31#include "G4Track.hh"
32
33#include "F04SteppingAction.hh"
34#include "G4SteppingManager.hh"
35
36#include "G4ParticleTypes.hh"
37#include "G4LogicalVolumeStore.hh"
38
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
44{
45 G4Track* theTrack = theStep->GetTrack();
46
47 // Get pointers to test volumes (only once)
48 if ( ! fTargetVolume ) {
50 = G4LogicalVolumeStore::GetInstance()->GetVolume("Target");
52 = G4LogicalVolumeStore::GetInstance()->GetVolume("TestPlane");
53 }
54
55 if (theTrack->GetParentID()==0) {
56 //This is a primary track
57 G4LogicalVolume* theVolume
58 = theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
59 if (theVolume != fTargetVolume) {
60 theTrack->SetTrackStatus(fStopAndKill);
61 return;
62 }
63 }
64
65// G4ParticleDefinition* particleType = theTrack->GetDefinition();
66
67 // check if it is entering the test volume
68 G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
69 G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
70 G4LogicalVolume* thePreLV = thePrePV->GetLogicalVolume();
71
72 G4LogicalVolume* thePostLV = nullptr;
73 G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
74
75 if (thePostPoint) {
76 G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
77 if (thePostPV) thePostLV = thePostPV->GetLogicalVolume();
78 }
79
80 if (thePostLV == fTestPlaneVolume &&
81 thePreLV != fTestPlaneVolume) {
82
83// G4double x = theTrack->GetPosition().x();
84// G4double y = theTrack->GetPosition().y();
85
86 G4ThreeVector theMomentumDirection = theTrack->
87 GetDynamicParticle()->GetMomentumDirection();
88// G4double theTotalMomentum = theTrack->GetDynamicParticle()->
89// GetTotalMomentum();
90
91 // then kill the track
92 theTrack->SetTrackStatus(fStopAndKill);
93 return;
94
95 }
96
97// G4double z = theTrack->GetPosition().z();
98
99 auto trackInformation =
100 (F04UserTrackInformation*)theTrack->GetUserInformation();
101
102 if (trackInformation->GetTrackStatusFlag() != reverse) {
103 if (thePreLV != fTargetVolume ) {
104 if ( theTrack-> GetMomentumDirection().z()>0.0 &&
105 theTrack->GetVertexMomentumDirection().z()<0.0 )
106 {
107 trackInformation->SetTrackStatusFlag(reverse);
108 }
109 }
110 }
111
112 // check if it is alive
113 if (theTrack->GetTrackStatus() == fAlive) { return; }
114
115 if (thePostPoint->GetProcessDefinedStep() != nullptr) {
116 if (thePostPoint->GetProcessDefinedStep()->GetProcessName() != "Decay")
117 return;
118 }
119
120}
Definition of the F04SteppingAction class.
Definition of the F04UserTrackInformation class.
G4LogicalVolume * fTargetVolume
void UserSteppingAction(const G4Step *) override
G4LogicalVolume * fTestPlaneVolume

Applications | User Support | Publications | Collaboration