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

#include <G4AdjointTrackingAction.hh>

+ Inheritance diagram for G4AdjointTrackingAction:

Public Member Functions

 G4AdjointTrackingAction (G4AdjointSteppingAction *anAction)
 
virtual ~G4AdjointTrackingAction ()
 
virtual void PreUserTrackingAction (const G4Track *)
 
virtual void PostUserTrackingAction (const G4Track *)
 
void RegisterAtEndOfAdjointTrack ()
 
void ClearEndOfAdjointTrackInfoVectors ()
 
void SetUserForwardTrackingAction (G4UserTrackingAction *anAction)
 
G4ThreeVector GetPositionAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4ThreeVector GetDirectionAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4double GetEkinAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4double GetEkinNucAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4double GetWeightAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4double GetCosthAtEndOfLastAdjointTrack (std::size_t i=0)
 
const G4StringGetFwdParticleNameAtEndOfLastAdjointTrack ()
 
G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4bool GetIsAdjointTrackingMode ()
 
G4int GetLastFwdParticleIndex (std::size_t i=0)
 
std::size_t GetNbOfAdointTracksReachingTheExternalSurface ()
 
void SetListOfPrimaryFwdParticles (std::vector< G4ParticleDefinition * > *aListOfParticles)
 
- Public Member Functions inherited from G4UserTrackingAction
 G4UserTrackingAction ()
 
virtual ~G4UserTrackingAction ()
 
virtual void SetTrackingManagerPointer (G4TrackingManager *pValue)
 
virtual void PreUserTrackingAction (const G4Track *)
 
virtual void PostUserTrackingAction (const G4Track *)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserTrackingAction
G4TrackingManagerfpTrackingManager = nullptr
 

Detailed Description

Definition at line 49 of file G4AdjointTrackingAction.hh.

Constructor & Destructor Documentation

◆ G4AdjointTrackingAction()

G4AdjointTrackingAction::G4AdjointTrackingAction ( G4AdjointSteppingAction anAction)

Definition at line 40 of file G4AdjointTrackingAction.cc.

42 : theAdjointSteppingAction(anAction)
43{
44}

◆ ~G4AdjointTrackingAction()

G4AdjointTrackingAction::~G4AdjointTrackingAction ( )
virtual

Definition at line 48 of file G4AdjointTrackingAction.cc.

49{
50}

Member Function Documentation

◆ ClearEndOfAdjointTrackInfoVectors()

void G4AdjointTrackingAction::ClearEndOfAdjointTrackInfoVectors ( )

Definition at line 138 of file G4AdjointTrackingAction.cc.

139{
140 last_pos_vec.clear();
141 last_direction_vec.clear();
142 last_ekin_vec.clear();
143 last_ekin_nuc_vec.clear();
144 last_cos_th_vec.clear();
145 last_weight_vec.clear();
146 last_fwd_part_PDGEncoding_vec.clear();
147 last_fwd_part_index_vec.clear();
148}

Referenced by G4AdjointSimManager::ClearEndOfAdjointTrackInfoVectors().

◆ GetCosthAtEndOfLastAdjointTrack()

G4double G4AdjointTrackingAction::GetCosthAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 74 of file G4AdjointTrackingAction.hh.

75 { return last_cos_th_vec[i]; }

Referenced by G4AdjointSimManager::GetCosthAtEndOfLastAdjointTrack().

◆ GetDirectionAtEndOfLastAdjointTrack()

G4ThreeVector G4AdjointTrackingAction::GetDirectionAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 66 of file G4AdjointTrackingAction.hh.

67 { return last_direction_vec[i]; }

Referenced by G4AdjointSimManager::GetDirectionAtEndOfLastAdjointTrack().

◆ GetEkinAtEndOfLastAdjointTrack()

G4double G4AdjointTrackingAction::GetEkinAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 68 of file G4AdjointTrackingAction.hh.

69 { return last_ekin_vec[i]; }

Referenced by G4AdjointSimManager::GetEkinAtEndOfLastAdjointTrack().

◆ GetEkinNucAtEndOfLastAdjointTrack()

G4double G4AdjointTrackingAction::GetEkinNucAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 70 of file G4AdjointTrackingAction.hh.

71 { return last_ekin_nuc_vec[i]; }

Referenced by G4AdjointSimManager::GetEkinNucAtEndOfLastAdjointTrack().

◆ GetFwdParticleNameAtEndOfLastAdjointTrack()

const G4String & G4AdjointTrackingAction::GetFwdParticleNameAtEndOfLastAdjointTrack ( )
inline

Definition at line 76 of file G4AdjointTrackingAction.hh.

77 { return last_fwd_part_name; }

Referenced by G4AdjointSimManager::GetFwdParticleNameAtEndOfLastAdjointTrack().

◆ GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack()

G4int G4AdjointTrackingAction::GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 78 of file G4AdjointTrackingAction.hh.

79 { return last_fwd_part_PDGEncoding_vec[i]; }

Referenced by G4AdjointSimManager::GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack().

◆ GetIsAdjointTrackingMode()

G4bool G4AdjointTrackingAction::GetIsAdjointTrackingMode ( )
inline

Definition at line 80 of file G4AdjointTrackingAction.hh.

81 { return is_adjoint_tracking_mode; }

Referenced by G4AdjointSimManager::GetAdjointTrackingMode().

◆ GetLastFwdParticleIndex()

G4int G4AdjointTrackingAction::GetLastFwdParticleIndex ( std::size_t  i = 0)
inline

Definition at line 82 of file G4AdjointTrackingAction.hh.

83 { return last_fwd_part_index_vec[i]; }

Referenced by G4AdjointSimManager::GetFwdParticleIndexAtEndOfLastAdjointTrack().

◆ GetNbOfAdointTracksReachingTheExternalSurface()

std::size_t G4AdjointTrackingAction::GetNbOfAdointTracksReachingTheExternalSurface ( )
inline

◆ GetPositionAtEndOfLastAdjointTrack()

G4ThreeVector G4AdjointTrackingAction::GetPositionAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 64 of file G4AdjointTrackingAction.hh.

65 { return last_pos_vec[i]; }

Referenced by G4AdjointSimManager::GetPositionAtEndOfLastAdjointTrack().

◆ GetWeightAtEndOfLastAdjointTrack()

G4double G4AdjointTrackingAction::GetWeightAtEndOfLastAdjointTrack ( std::size_t  i = 0)
inline

Definition at line 72 of file G4AdjointTrackingAction.hh.

73 { return last_weight_vec[i]; }

Referenced by G4AdjointSimManager::GetWeightAtEndOfLastAdjointTrack().

◆ PostUserTrackingAction()

void G4AdjointTrackingAction::PostUserTrackingAction ( const G4Track aTrack)
virtual

Reimplemented from G4UserTrackingAction.

Definition at line 75 of file G4AdjointTrackingAction.cc.

76{
77 // important to have it here !
78 //
79 last_weight = theAdjointSteppingAction->GetLastWeight();
80 last_ekin = theAdjointSteppingAction->GetLastEkin();
81
82 if(!is_adjoint_tracking_mode)
83 {
84 if (theUserFwdTrackingAction != nullptr)
85 {
86 theUserFwdTrackingAction->PostUserTrackingAction(aTrack);
87 }
88 }
89 else if (theAdjointSteppingAction->GetDidAdjParticleReachTheExtSource())
90 {
91 last_pos = theAdjointSteppingAction->GetLastPosition();
92 last_direction = theAdjointSteppingAction->GetLastMomentum();
93 last_direction /=last_direction.mag();
94 last_cos_th = last_direction.z();
95 G4ParticleDefinition* aPartDef= theAdjointSteppingAction->GetLastPartDef();
96 last_fwd_part_name= aPartDef->GetParticleName();
97 last_fwd_part_name.erase(0,4);
98 last_fwd_part_PDGEncoding=G4ParticleTable::GetParticleTable()
99 ->FindParticle(last_fwd_part_name)->GetPDGEncoding();
100 last_ekin = theAdjointSteppingAction->GetLastEkin();
101 last_ekin_nuc = last_ekin;
102 if (aPartDef->GetParticleType() == "adjoint_nucleus")
103 {
104 G4double nb_nuc = G4double(aPartDef->GetBaryonNumber());
105 last_ekin_nuc /= nb_nuc;
106 }
107
108 last_fwd_part_index = -1;
109 std::size_t i = 0;
110 while(i< pListOfPrimaryFwdParticles->size() && last_fwd_part_index<0)
111 {
112 if ((*pListOfPrimaryFwdParticles)[i]->GetParticleName()
113 == last_fwd_part_name)
114 {
115 last_fwd_part_index=(G4int)i;
116 }
117 ++i;
118 }
119
120 // Fill the vectors
121 //
122 last_pos_vec.push_back(last_pos);
123 last_direction_vec.push_back(last_direction);
124 last_ekin_vec.push_back(last_ekin);
125 last_ekin_nuc_vec.push_back(last_ekin_nuc);
126 last_cos_th_vec.push_back(last_cos_th);
127 last_weight_vec.push_back(last_weight);
128 last_fwd_part_PDGEncoding_vec.push_back(last_fwd_part_PDGEncoding);
129 last_fwd_part_index_vec.push_back(last_fwd_part_index);
130 }
131 else
132 {
133 }
134}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double z() const
double mag() const
G4ParticleDefinition * GetLastPartDef()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
virtual void PostUserTrackingAction(const G4Track *)

◆ PreUserTrackingAction()

void G4AdjointTrackingAction::PreUserTrackingAction ( const G4Track aTrack)
virtual

Reimplemented from G4UserTrackingAction.

Definition at line 54 of file G4AdjointTrackingAction.cc.

55{
56 G4String partType = aTrack->GetParticleDefinition()->GetParticleType();
57 if (G4StrUtil::contains(partType, "adjoint"))
58 {
59 is_adjoint_tracking_mode = true;
60 theAdjointSteppingAction->SetPrimWeight(aTrack->GetWeight());
61 }
62 else
63 {
64 is_adjoint_tracking_mode = false;
65 if (theUserFwdTrackingAction)
66 {
67 theUserFwdTrackingAction->PreUserTrackingAction(aTrack);
68 }
69 }
70 theAdjointSteppingAction->SetAdjointTrackingMode(is_adjoint_tracking_mode);
71}
void SetAdjointTrackingMode(G4bool aBool)
void SetPrimWeight(G4double weight)
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetWeight() const
virtual void PreUserTrackingAction(const G4Track *)
G4bool contains(const G4String &str, std::string_view ss)
Check if a string contains a given substring.

◆ RegisterAtEndOfAdjointTrack()

void G4AdjointTrackingAction::RegisterAtEndOfAdjointTrack ( )

◆ SetListOfPrimaryFwdParticles()

void G4AdjointTrackingAction::SetListOfPrimaryFwdParticles ( std::vector< G4ParticleDefinition * > *  aListOfParticles)
inline

Definition at line 86 of file G4AdjointTrackingAction.hh.

87 { pListOfPrimaryFwdParticles = aListOfParticles; }

◆ SetUserForwardTrackingAction()

void G4AdjointTrackingAction::SetUserForwardTrackingAction ( G4UserTrackingAction anAction)
inline

Definition at line 62 of file G4AdjointTrackingAction.hh.

63 { theUserFwdTrackingAction = anAction; }

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