Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointSteppingAction Class Reference

#include <G4AdjointSteppingAction.hh>

+ Inheritance diagram for G4AdjointSteppingAction:

Public Member Functions

 G4AdjointSteppingAction ()
 
 ~G4AdjointSteppingAction ()
 
void UserSteppingAction (const G4Step *)
 
void SetExtSourceEMax (G4double Emax)
 
void SetStartEvent (G4bool aBool)
 
G4bool GetDidAdjParticleReachTheExtSource ()
 
G4ThreeVector GetLastMomentum ()
 
G4ThreeVector GetLastPosition ()
 
G4double GetLastEkin ()
 
G4double GetLastWeight ()
 
void SetPrimWeight (G4double weight)
 
G4ParticleDefinitionGetLastPartDef ()
 
void SetUserAdjointSteppingAction (G4UserSteppingAction *anAction)
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
void SetSteppingManagerPointer (G4SteppingManager *pValue)
 
virtual void UserSteppingAction (const G4Step *)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager
 

Detailed Description

Definition at line 69 of file G4AdjointSteppingAction.hh.

Constructor & Destructor Documentation

◆ G4AdjointSteppingAction()

G4AdjointSteppingAction::G4AdjointSteppingAction ( )

Definition at line 43 of file G4AdjointSteppingAction.cc.

44 : ext_sourceEMax(0.), start_event(false),
45 did_adj_part_reach_ext_source(false), last_ekin(0.), last_weight(0.),
46 prim_weight(0.), last_part_def(0), theUserAdjointSteppingAction(0)
47{
48 theG4AdjointCrossSurfChecker = G4AdjointCrossSurfChecker::GetInstance();
49}
static G4AdjointCrossSurfChecker * GetInstance()

◆ ~G4AdjointSteppingAction()

G4AdjointSteppingAction::~G4AdjointSteppingAction ( )

Definition at line 52 of file G4AdjointSteppingAction.cc.

53{;}

Member Function Documentation

◆ GetDidAdjParticleReachTheExtSource()

G4bool G4AdjointSteppingAction::GetDidAdjParticleReachTheExtSource ( )
inline

Definition at line 79 of file G4AdjointSteppingAction.hh.

79{return did_adj_part_reach_ext_source;}

Referenced by G4AdjointSimManager::GetDidAdjParticleReachTheExtSource().

◆ GetLastEkin()

G4double G4AdjointSteppingAction::GetLastEkin ( )
inline

Definition at line 82 of file G4AdjointSteppingAction.hh.

82{return last_ekin;}

Referenced by G4AdjointSimManager::RegisterAtEndOfAdjointTrack().

◆ GetLastMomentum()

G4ThreeVector G4AdjointSteppingAction::GetLastMomentum ( )
inline

Definition at line 80 of file G4AdjointSteppingAction.hh.

80{return last_momentum;}

Referenced by G4AdjointSimManager::RegisterAtEndOfAdjointTrack().

◆ GetLastPartDef()

G4ParticleDefinition * G4AdjointSteppingAction::GetLastPartDef ( )
inline

Definition at line 85 of file G4AdjointSteppingAction.hh.

85{return last_part_def;}

Referenced by G4AdjointSimManager::RegisterAtEndOfAdjointTrack().

◆ GetLastPosition()

G4ThreeVector G4AdjointSteppingAction::GetLastPosition ( )
inline

Definition at line 81 of file G4AdjointSteppingAction.hh.

81{return last_pos;}

Referenced by G4AdjointSimManager::RegisterAtEndOfAdjointTrack().

◆ GetLastWeight()

G4double G4AdjointSteppingAction::GetLastWeight ( )
inline

Definition at line 83 of file G4AdjointSteppingAction.hh.

83{return last_weight;}

Referenced by G4AdjointSimManager::RegisterAtEndOfAdjointTrack().

◆ SetExtSourceEMax()

void G4AdjointSteppingAction::SetExtSourceEMax ( G4double  Emax)
inline

Definition at line 77 of file G4AdjointSteppingAction.hh.

77{ext_sourceEMax=Emax;}

Referenced by G4AdjointSimManager::SetExtSourceEmax().

◆ SetPrimWeight()

void G4AdjointSteppingAction::SetPrimWeight ( G4double  weight)
inline

Definition at line 84 of file G4AdjointSteppingAction.hh.

84{prim_weight=weight;}

Referenced by G4AdjointSimManager::RegisterAdjointPrimaryWeight().

◆ SetStartEvent()

void G4AdjointSteppingAction::SetStartEvent ( G4bool  aBool)
inline

Definition at line 78 of file G4AdjointSteppingAction.hh.

78{start_event =aBool;}

◆ SetUserAdjointSteppingAction()

void G4AdjointSteppingAction::SetUserAdjointSteppingAction ( G4UserSteppingAction anAction)
inline

Definition at line 86 of file G4AdjointSteppingAction.hh.

86{theUserAdjointSteppingAction = anAction;}

Referenced by G4AdjointSimManager::SetAdjointSteppingAction().

◆ UserSteppingAction()

void G4AdjointSteppingAction::UserSteppingAction ( const G4Step aStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 57 of file G4AdjointSteppingAction.cc.

58{
59 //Apply first the user adjoint stepping action
60 //---------------------------
61 if (theUserAdjointSteppingAction) theUserAdjointSteppingAction->UserSteppingAction(aStep);
62
63 G4Track* aTrack =aStep->GetTrack();
64 G4double nb_nuc=1.;
65 G4ParticleDefinition* thePartDef = aTrack->GetDefinition();
66
67 if (thePartDef->GetParticleType() == "adjoint_nucleus"){
68 nb_nuc=double(thePartDef->GetBaryonNumber());
69 }
70 //Kill conditions for adjoint particles reaching the maximum energy
71 //-----------------------------------------------------------------
72 if(aTrack->GetKineticEnergy() >= ext_sourceEMax*nb_nuc){
74 did_adj_part_reach_ext_source=false;
75 return;
76 }
77
78 G4double weight_factor = aTrack->GetWeight()/prim_weight;
79
80 if ( (weight_factor>0 && weight_factor<=0) || weight_factor<= 1e-290 || weight_factor>1.e200)
81 {
82 //std::cout<<"Weight_factor problem! Value = "<<weight_factor<<std::endl;
84 did_adj_part_reach_ext_source=false;
85 return;
86 }
87
88
89 //Kill conditions for surface crossing
90 //--------------------------------------
91
92 G4String surface_name;
93 G4double cos_to_surface;
94 G4bool GoingIn;
95 G4ThreeVector crossing_pos;
96 if (theG4AdjointCrossSurfChecker->CrossingOneOfTheRegisteredSurface(aStep, surface_name, crossing_pos, cos_to_surface, GoingIn) ){
97
98 //G4cout<<"Test_step11"<<std::endl;
99 //G4cout<<surface_name<<std::endl;
100 if (surface_name == "ExternalSource") {
101 //Registering still needed
102 did_adj_part_reach_ext_source=true;
104 //now register the adjoint particles reaching the external surface
105 last_momentum =aTrack->GetMomentum();
106 last_ekin=aTrack->GetKineticEnergy();
107 last_weight = aTrack->GetWeight();
108 last_part_def = aTrack->GetDefinition();
109 last_pos = crossing_pos;
110 return;
111 }
112 else if (surface_name == "AdjointSource" && GoingIn) {
113 did_adj_part_reach_ext_source=false;
115 return;
116 }
117 }
118 //Check for reaching out of world
119 if (aStep->GetPostStepPoint()->GetStepStatus() == fWorldBoundary) {
120 did_adj_part_reach_ext_source=true;
121 last_momentum =aTrack->GetMomentum();
122 last_ekin=aTrack->GetKineticEnergy();
123 last_weight = aTrack->GetWeight();
124 last_part_def = aTrack->GetDefinition();
125 last_pos = crossing_pos;
126 return;
127
128 }
129
130}
@ fWorldBoundary
Definition: G4StepStatus.hh:52
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4bool CrossingOneOfTheRegisteredSurface(const G4Step *aStep, G4String &surface_name, G4ThreeVector &crossing_pos, G4double &cos_to_surface, G4bool &GoingIn)
const G4String & GetParticleType() const
G4StepStatus GetStepStatus() const
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4double GetWeight() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
virtual void UserSteppingAction(const G4Step *)

The documentation for this class was generated from the following files: