Geant4 11.2.2
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)
 
 ~G4AdjointTrackingAction () override=default
 
void PreUserTrackingAction (const G4Track *) override
 
void PostUserTrackingAction (const G4Track *) override
 
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 ()=default
 
virtual void SetTrackingManagerPointer (G4TrackingManager *pValue)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserTrackingAction
G4TrackingManagerfpTrackingManager = nullptr
 

Detailed Description

Definition at line 50 of file G4AdjointTrackingAction.hh.

Constructor & Destructor Documentation

◆ G4AdjointTrackingAction()

G4AdjointTrackingAction::G4AdjointTrackingAction ( G4AdjointSteppingAction * anAction)

Definition at line 41 of file G4AdjointTrackingAction.cc.

42 : theAdjointSteppingAction(anAction)
43{}

◆ ~G4AdjointTrackingAction()

G4AdjointTrackingAction::~G4AdjointTrackingAction ( )
overridedefault

Member Function Documentation

◆ ClearEndOfAdjointTrackInfoVectors()

void G4AdjointTrackingAction::ClearEndOfAdjointTrackInfoVectors ( )

Definition at line 120 of file G4AdjointTrackingAction.cc.

121{
122 last_pos_vec.clear();
123 last_direction_vec.clear();
124 last_ekin_vec.clear();
125 last_ekin_nuc_vec.clear();
126 last_cos_th_vec.clear();
127 last_weight_vec.clear();
128 last_fwd_part_PDGEncoding_vec.clear();
129 last_fwd_part_index_vec.clear();
130}

Referenced by G4AdjointSimManager::ClearEndOfAdjointTrackInfoVectors().

◆ GetCosthAtEndOfLastAdjointTrack()

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

Definition at line 79 of file G4AdjointTrackingAction.hh.

79{ return last_cos_th_vec[i]; }

Referenced by G4AdjointSimManager::GetCosthAtEndOfLastAdjointTrack().

◆ GetDirectionAtEndOfLastAdjointTrack()

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

Definition at line 69 of file G4AdjointTrackingAction.hh.

70 {
71 return last_direction_vec[i];
72 }

Referenced by G4AdjointSimManager::GetDirectionAtEndOfLastAdjointTrack().

◆ GetEkinAtEndOfLastAdjointTrack()

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

Definition at line 73 of file G4AdjointTrackingAction.hh.

73{ return last_ekin_vec[i]; }

Referenced by G4AdjointSimManager::GetEkinAtEndOfLastAdjointTrack().

◆ GetEkinNucAtEndOfLastAdjointTrack()

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

Definition at line 74 of file G4AdjointTrackingAction.hh.

75 {
76 return last_ekin_nuc_vec[i];
77 }

Referenced by G4AdjointSimManager::GetEkinNucAtEndOfLastAdjointTrack().

◆ GetFwdParticleNameAtEndOfLastAdjointTrack()

const G4String & G4AdjointTrackingAction::GetFwdParticleNameAtEndOfLastAdjointTrack ( )
inline

Definition at line 80 of file G4AdjointTrackingAction.hh.

80{ return last_fwd_part_name; }

Referenced by G4AdjointSimManager::GetFwdParticleNameAtEndOfLastAdjointTrack().

◆ GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack()

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

Definition at line 81 of file G4AdjointTrackingAction.hh.

82 {
83 return last_fwd_part_PDGEncoding_vec[i];
84 }

Referenced by G4AdjointSimManager::GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack().

◆ GetIsAdjointTrackingMode()

G4bool G4AdjointTrackingAction::GetIsAdjointTrackingMode ( )
inline

Definition at line 85 of file G4AdjointTrackingAction.hh.

85{ return is_adjoint_tracking_mode; }

Referenced by G4AdjointSimManager::GetAdjointTrackingMode().

◆ GetLastFwdParticleIndex()

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

Definition at line 86 of file G4AdjointTrackingAction.hh.

86{ 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 65 of file G4AdjointTrackingAction.hh.

66 {
67 return last_pos_vec[i];
68 }

Referenced by G4AdjointSimManager::GetPositionAtEndOfLastAdjointTrack().

◆ GetWeightAtEndOfLastAdjointTrack()

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

Definition at line 78 of file G4AdjointTrackingAction.hh.

78{ return last_weight_vec[i]; }

Referenced by G4AdjointSimManager::GetWeightAtEndOfLastAdjointTrack().

◆ PostUserTrackingAction()

void G4AdjointTrackingAction::PostUserTrackingAction ( const G4Track * aTrack)
overridevirtual

Reimplemented from G4UserTrackingAction.

Definition at line 65 of file G4AdjointTrackingAction.cc.

66{
67 // important to have it here !
68 //
69 last_weight = theAdjointSteppingAction->GetLastWeight();
70 last_ekin = theAdjointSteppingAction->GetLastEkin();
71
72 if (! is_adjoint_tracking_mode) {
73 if (theUserFwdTrackingAction != nullptr) {
74 theUserFwdTrackingAction->PostUserTrackingAction(aTrack);
75 }
76 }
77 else if (theAdjointSteppingAction->GetDidAdjParticleReachTheExtSource()) {
78 last_pos = theAdjointSteppingAction->GetLastPosition();
79 last_direction = theAdjointSteppingAction->GetLastMomentum();
80 last_direction /= last_direction.mag();
81 last_cos_th = last_direction.z();
82 G4ParticleDefinition* aPartDef = theAdjointSteppingAction->GetLastPartDef();
83 last_fwd_part_name = aPartDef->GetParticleName();
84 last_fwd_part_name.erase(0, 4);
85 last_fwd_part_PDGEncoding =
87 last_ekin = theAdjointSteppingAction->GetLastEkin();
88 last_ekin_nuc = last_ekin;
89 if (aPartDef->GetParticleType() == "adjoint_nucleus") {
90 auto nb_nuc = G4double(aPartDef->GetBaryonNumber());
91 last_ekin_nuc /= nb_nuc;
92 }
93
94 last_fwd_part_index = -1;
95 std::size_t i = 0;
96 while (i < pListOfPrimaryFwdParticles->size() && last_fwd_part_index < 0) {
97 if ((*pListOfPrimaryFwdParticles)[i]->GetParticleName() == last_fwd_part_name) {
98 last_fwd_part_index = (G4int)i;
99 }
100 ++i;
101 }
102
103 // Fill the vectors
104 //
105 last_pos_vec.push_back(last_pos);
106 last_direction_vec.push_back(last_direction);
107 last_ekin_vec.push_back(last_ekin);
108 last_ekin_nuc_vec.push_back(last_ekin_nuc);
109 last_cos_th_vec.push_back(last_cos_th);
110 last_weight_vec.push_back(last_weight);
111 last_fwd_part_PDGEncoding_vec.push_back(last_fwd_part_PDGEncoding);
112 last_fwd_part_index_vec.push_back(last_fwd_part_index);
113 }
114 else {
115 }
116}
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)
overridevirtual

Reimplemented from G4UserTrackingAction.

Definition at line 47 of file G4AdjointTrackingAction.cc.

48{
49 G4String partType = aTrack->GetParticleDefinition()->GetParticleType();
50 if (G4StrUtil::contains(partType, "adjoint")) {
51 is_adjoint_tracking_mode = true;
52 theAdjointSteppingAction->SetPrimWeight(aTrack->GetWeight());
53 }
54 else {
55 is_adjoint_tracking_mode = false;
56 if (theUserFwdTrackingAction != nullptr) {
57 theUserFwdTrackingAction->PreUserTrackingAction(aTrack);
58 }
59 }
60 theAdjointSteppingAction->SetAdjointTrackingMode(is_adjoint_tracking_mode);
61}
void SetAdjointTrackingMode(G4bool aBool)
void SetPrimWeight(G4double weight)
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetWeight() const
virtual void PreUserTrackingAction(const G4Track *)

◆ RegisterAtEndOfAdjointTrack()

void G4AdjointTrackingAction::RegisterAtEndOfAdjointTrack ( )

◆ SetListOfPrimaryFwdParticles()

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

Definition at line 88 of file G4AdjointTrackingAction.hh.

89 {
90 pListOfPrimaryFwdParticles = aListOfParticles;
91 }

◆ SetUserForwardTrackingAction()

void G4AdjointTrackingAction::SetUserForwardTrackingAction ( G4UserTrackingAction * anAction)
inline

Definition at line 61 of file G4AdjointTrackingAction.hh.

62 {
63 theUserFwdTrackingAction = anAction;
64 }

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