Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4AdjointSteppingAction Class Reference

#include <G4AdjointSteppingAction.hh>

Inheritance diagram for G4AdjointSteppingAction:
Collaboration 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)
 
void SetUserForwardSteppingAction (G4UserSteppingAction *anAction)
 
void SetAdjointTrackingMode (G4bool aBool)
 
void ResetDidOneAdjPartReachExtSourceDuringEvent ()
 
void SetAdjointGeantinoTrackingMode (G4bool aBool)
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
virtual void SetSteppingManagerPointer (G4SteppingManager *pValue)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager
 

Detailed Description

Definition at line 69 of file G4AdjointSteppingAction.hh.

Constructor & Destructor Documentation

G4AdjointSteppingAction::G4AdjointSteppingAction ( )

Definition at line 44 of file G4AdjointSteppingAction.cc.

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

Here is the call graph for this function:

G4AdjointSteppingAction::~G4AdjointSteppingAction ( )

Definition at line 54 of file G4AdjointSteppingAction.cc.

55 {;}

Member Function Documentation

G4bool G4AdjointSteppingAction::GetDidAdjParticleReachTheExtSource ( )
inline

Definition at line 79 of file G4AdjointSteppingAction.hh.

79 {return did_adj_part_reach_ext_source;}

Here is the caller graph for this function:

G4double G4AdjointSteppingAction::GetLastEkin ( )
inline

Definition at line 82 of file G4AdjointSteppingAction.hh.

82 {return last_ekin;}

Here is the caller graph for this function:

G4ThreeVector G4AdjointSteppingAction::GetLastMomentum ( )
inline

Definition at line 80 of file G4AdjointSteppingAction.hh.

80 {return last_momentum;}

Here is the caller graph for this function:

G4ParticleDefinition* G4AdjointSteppingAction::GetLastPartDef ( )
inline

Definition at line 85 of file G4AdjointSteppingAction.hh.

85 {return last_part_def;}

Here is the caller graph for this function:

G4ThreeVector G4AdjointSteppingAction::GetLastPosition ( )
inline

Definition at line 81 of file G4AdjointSteppingAction.hh.

81 {return last_pos;}

Here is the caller graph for this function:

G4double G4AdjointSteppingAction::GetLastWeight ( )
inline

Definition at line 83 of file G4AdjointSteppingAction.hh.

83 {return last_weight;}

Here is the caller graph for this function:

void G4AdjointSteppingAction::ResetDidOneAdjPartReachExtSourceDuringEvent ( )
inline

Definition at line 91 of file G4AdjointSteppingAction.hh.

92  {did_one_adj_part_reach_ext_source_during_event =false;}

Here is the caller graph for this function:

void G4AdjointSteppingAction::SetAdjointGeantinoTrackingMode ( G4bool  aBool)
inline

Definition at line 93 of file G4AdjointSteppingAction.hh.

93 {is_adjoint_geantino_tracking_mode =aBool;}
void G4AdjointSteppingAction::SetAdjointTrackingMode ( G4bool  aBool)
inline

Definition at line 90 of file G4AdjointSteppingAction.hh.

90 {is_adjoint_tracking_mode =aBool;}

Here is the caller graph for this function:

void G4AdjointSteppingAction::SetExtSourceEMax ( G4double  Emax)
inline

Definition at line 77 of file G4AdjointSteppingAction.hh.

77 {ext_sourceEMax=Emax;}
static const G4double Emax

Here is the caller graph for this function:

void G4AdjointSteppingAction::SetPrimWeight ( G4double  weight)
inline

Definition at line 84 of file G4AdjointSteppingAction.hh.

84 {prim_weight=weight;}

Here is the caller graph for this function:

void G4AdjointSteppingAction::SetStartEvent ( G4bool  aBool)
inline

Definition at line 78 of file G4AdjointSteppingAction.hh.

78 {start_event =aBool;}
void G4AdjointSteppingAction::SetUserAdjointSteppingAction ( G4UserSteppingAction anAction)
inline

Definition at line 86 of file G4AdjointSteppingAction.hh.

86  {
87  theUserAdjointSteppingAction = anAction;}

Here is the caller graph for this function:

void G4AdjointSteppingAction::SetUserForwardSteppingAction ( G4UserSteppingAction anAction)
inline

Definition at line 88 of file G4AdjointSteppingAction.hh.

88  {
89  theUserFwdSteppingAction = anAction;}
void G4AdjointSteppingAction::UserSteppingAction ( const G4Step aStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 59 of file G4AdjointSteppingAction.cc.

60 {
61  G4Track* aTrack =aStep->GetTrack();
62  //G4cout<<"Ste weight "<<<<std:.endl;
63  //forward tracking mode
64  if(!is_adjoint_tracking_mode){
65  if (!did_one_adj_part_reach_ext_source_during_event) {
67  return;
68  }
69  if(theUserFwdSteppingAction)
70  theUserFwdSteppingAction->UserSteppingAction(aStep);
71  return;
72  }
73  //Apply first the user adjoint stepping action
74  //---------------------------
75  did_adj_part_reach_ext_source=false;
76  if (theUserAdjointSteppingAction) theUserAdjointSteppingAction->UserSteppingAction(aStep);
77 
78 
79  G4double nb_nuc=1.;
80  G4ParticleDefinition* thePartDef = aTrack->GetDefinition();
81 
82  if (thePartDef->GetParticleType() == "adjoint_nucleus"){
83  nb_nuc=double(thePartDef->GetBaryonNumber());
84  }
85  //Kill conditions for adjoint particles reaching the maximum energy
86  //-----------------------------------------------------------------
87  if(aTrack->GetKineticEnergy() >= ext_sourceEMax*nb_nuc){
89  did_adj_part_reach_ext_source=false;
90  return;
91  }
92  //G4cout<<"Weight adjoint "<<aTrack->GetWeight()<<std::endl;
93  /*G4double weight_factor = aTrack->GetWeight()/prim_weight;
94  if ( (weight_factor>0 && weight_factor<=0) || weight_factor<= 1e-290 || weight_factor>1.e200)
95  {
96  //std::cout<<"Weight_factor problem! Value = "<<weight_factor<<std::endl;
97  aTrack->SetTrackStatus(fStopAndKill);
98  did_adj_part_reach_ext_source=false;
99  //G4cout<<thePartDef->GetParticleName()<<std::endl;
100  return;
101  }
102  */
103 
104 
105  //Kill conditions for surface crossing
106  //--------------------------------------
107 
108  G4String surface_name;
109  G4double cos_to_surface;
110  G4bool GoingIn;
111  G4ThreeVector crossing_pos;
112  if (theG4AdjointCrossSurfChecker->CrossingOneOfTheRegisteredSurface(aStep, surface_name, crossing_pos, cos_to_surface, GoingIn) ){
113 
114  //G4cout<<"Test_step11"<<std::endl;
115  if (surface_name == "ExternalSource") {
116  //Registering still needed
117  did_adj_part_reach_ext_source=true;
118  did_one_adj_part_reach_ext_source_during_event=true;
119  aTrack->SetTrackStatus(fStopAndKill);
120  //now register the adjoint particles reaching the external surface
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  else if (surface_name == "AdjointSource" && GoingIn) {
129  did_adj_part_reach_ext_source=false;
130  aTrack->SetTrackStatus(fStopAndKill);
131  return;
132  }
133  }
134  //Check for reaching out of world
135  //G4cout<<aStep->GetPostStepPoint()->GetStepStatus()<<std::endl;
136  if (aStep->GetPostStepPoint()->GetStepStatus() == fWorldBoundary) {
137  did_adj_part_reach_ext_source=true;
138  did_one_adj_part_reach_ext_source_during_event=true;
139  last_momentum =aTrack->GetMomentum();
140  last_ekin=aTrack->GetKineticEnergy();
141  last_weight = aTrack->GetWeight();
142  last_part_def = aTrack->GetDefinition();
143  last_pos = crossing_pos;
144  return;
145  }
146 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
G4StepStatus GetStepStatus() const
virtual void UserSteppingAction(const G4Step *)
G4double GetKineticEnergy() const
bool G4bool
Definition: G4Types.hh:79
const G4String & GetParticleType() const
G4ThreeVector GetMomentum() const
G4bool CrossingOneOfTheRegisteredSurface(const G4Step *aStep, G4String &surface_name, G4ThreeVector &crossing_pos, G4double &cos_to_surface, G4bool &GoingIn)
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const

Here is the call graph for this function:


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