Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4AdjointTrackingAction Class Reference

#include <G4AdjointTrackingAction.hh>

Inheritance diagram for G4AdjointTrackingAction:
G4UserTrackingAction

Public Member Functions

void ClearEndOfAdjointTrackInfoVectors ()
 
 G4AdjointTrackingAction (G4AdjointSteppingAction *anAction)
 
G4double GetCosthAtEndOfLastAdjointTrack (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)
 
const G4StringGetFwdParticleNameAtEndOfLastAdjointTrack ()
 
G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4bool GetIsAdjointTrackingMode ()
 
G4int GetLastFwdParticleIndex (std::size_t i=0)
 
std::size_t GetNbOfAdointTracksReachingTheExternalSurface ()
 
G4ThreeVector GetPositionAtEndOfLastAdjointTrack (std::size_t i=0)
 
G4double GetWeightAtEndOfLastAdjointTrack (std::size_t i=0)
 
virtual void PostUserTrackingAction (const G4Track *)
 
virtual void PreUserTrackingAction (const G4Track *)
 
void RegisterAtEndOfAdjointTrack ()
 
void SetListOfPrimaryFwdParticles (std::vector< G4ParticleDefinition * > *aListOfParticles)
 
virtual void SetTrackingManagerPointer (G4TrackingManager *pValue)
 
void SetUserForwardTrackingAction (G4UserTrackingAction *anAction)
 
virtual ~G4AdjointTrackingAction ()
 

Protected Attributes

G4TrackingManagerfpTrackingManager = nullptr
 

Private Attributes

G4bool is_adjoint_tracking_mode = false
 
G4double last_cos_th = 0.0
 
std::vector< G4doublelast_cos_th_vec
 
G4ThreeVector last_direction
 
std::vector< G4ThreeVectorlast_direction_vec
 
G4double last_ekin = 0.0
 
G4double last_ekin_nuc = 0.0
 
std::vector< G4doublelast_ekin_nuc_vec
 
std::vector< G4doublelast_ekin_vec
 
G4int last_fwd_part_index = 0
 
std::vector< G4intlast_fwd_part_index_vec
 
G4String last_fwd_part_name
 
G4int last_fwd_part_PDGEncoding = 0
 
std::vector< G4intlast_fwd_part_PDGEncoding_vec
 
G4ThreeVector last_pos
 
std::vector< G4ThreeVectorlast_pos_vec
 
G4double last_weight = 0.0
 
std::vector< G4doublelast_weight_vec
 
std::vector< G4ParticleDefinition * > * pListOfPrimaryFwdParticles = nullptr
 
G4AdjointSteppingActiontheAdjointSteppingAction = nullptr
 
G4UserTrackingActiontheUserFwdTrackingAction = 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}
G4AdjointSteppingAction * theAdjointSteppingAction

◆ ~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();
148}
std::vector< G4ThreeVector > last_direction_vec
std::vector< G4double > last_ekin_nuc_vec
std::vector< G4double > last_cos_th_vec
std::vector< G4ThreeVector > last_pos_vec
std::vector< G4double > last_weight_vec
std::vector< G4int > last_fwd_part_index_vec
std::vector< G4int > last_fwd_part_PDGEncoding_vec
std::vector< G4double > last_ekin_vec

References last_cos_th_vec, last_direction_vec, last_ekin_nuc_vec, last_ekin_vec, last_fwd_part_index_vec, last_fwd_part_PDGEncoding_vec, last_pos_vec, and last_weight_vec.

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]; }

References last_cos_th_vec.

Referenced by G4AdjointSimManager::GetCosthAtEndOfLastAdjointTrack().

◆ GetDirectionAtEndOfLastAdjointTrack()

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

◆ GetEkinAtEndOfLastAdjointTrack()

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

Definition at line 68 of file G4AdjointTrackingAction.hh.

69 { return last_ekin_vec[i]; }

References last_ekin_vec.

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]; }

References last_ekin_nuc_vec.

Referenced by G4AdjointSimManager::GetEkinNucAtEndOfLastAdjointTrack().

◆ GetFwdParticleNameAtEndOfLastAdjointTrack()

const G4String & G4AdjointTrackingAction::GetFwdParticleNameAtEndOfLastAdjointTrack ( )
inline

◆ GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack()

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

◆ GetIsAdjointTrackingMode()

G4bool G4AdjointTrackingAction::GetIsAdjointTrackingMode ( )
inline

◆ GetLastFwdParticleIndex()

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

◆ 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]; }

References last_pos_vec.

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]; }

References last_weight_vec.

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 //
81
83 {
84 if (theUserFwdTrackingAction != nullptr)
85 {
87 }
88 }
90 {
97 last_fwd_part_name.erase(0,4);
102 if (aPartDef->GetParticleType() == "adjoint_nucleus")
103 {
104 G4double nb_nuc = G4double(aPartDef->GetBaryonNumber());
105 last_ekin_nuc /= nb_nuc;
106 }
107
109 std::size_t i = 0;
110 while(i< pListOfPrimaryFwdParticles->size() && last_fwd_part_index<0)
111 {
112 if ((*pListOfPrimaryFwdParticles)[i]->GetParticleName()
114 {
116 }
117 ++i;
118 }
119
120 // Fill the vectors
121 //
122 last_pos_vec.push_back(last_pos);
124 last_ekin_vec.push_back(last_ekin);
126 last_cos_th_vec.push_back(last_cos_th);
127 last_weight_vec.push_back(last_weight);
130 }
131 else
132 {
133 }
134}
double G4double
Definition: G4Types.hh:83
double z() const
double mag() const
G4ParticleDefinition * GetLastPartDef()
G4UserTrackingAction * theUserFwdTrackingAction
std::vector< G4ParticleDefinition * > * pListOfPrimaryFwdParticles
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
virtual void PostUserTrackingAction(const G4Track *)

References G4ParticleTable::FindParticle(), G4ParticleDefinition::GetBaryonNumber(), G4AdjointSteppingAction::GetDidAdjParticleReachTheExtSource(), G4AdjointSteppingAction::GetLastEkin(), G4AdjointSteppingAction::GetLastMomentum(), G4AdjointSteppingAction::GetLastPartDef(), G4AdjointSteppingAction::GetLastPosition(), G4AdjointSteppingAction::GetLastWeight(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGEncoding(), is_adjoint_tracking_mode, last_cos_th, last_cos_th_vec, last_direction, last_direction_vec, last_ekin, last_ekin_nuc, last_ekin_nuc_vec, last_ekin_vec, last_fwd_part_index, last_fwd_part_index_vec, last_fwd_part_name, last_fwd_part_PDGEncoding, last_fwd_part_PDGEncoding_vec, last_pos, last_pos_vec, last_weight, last_weight_vec, CLHEP::Hep3Vector::mag(), pListOfPrimaryFwdParticles, G4UserTrackingAction::PostUserTrackingAction(), theAdjointSteppingAction, theUserFwdTrackingAction, and CLHEP::Hep3Vector::z().

◆ 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 {
61 }
62 else
63 {
66 {
68 }
69 }
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.

References G4StrUtil::contains(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleType(), G4Track::GetWeight(), is_adjoint_tracking_mode, G4UserTrackingAction::PreUserTrackingAction(), G4AdjointSteppingAction::SetAdjointTrackingMode(), G4AdjointSteppingAction::SetPrimWeight(), theAdjointSteppingAction, and theUserFwdTrackingAction.

◆ RegisterAtEndOfAdjointTrack()

void G4AdjointTrackingAction::RegisterAtEndOfAdjointTrack ( )

◆ SetListOfPrimaryFwdParticles()

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

Definition at line 86 of file G4AdjointTrackingAction.hh.

87 { pListOfPrimaryFwdParticles = aListOfParticles; }

References pListOfPrimaryFwdParticles.

Referenced by G4AdjointSimManager::G4AdjointSimManager().

◆ SetTrackingManagerPointer()

void G4UserTrackingAction::SetTrackingManagerPointer ( G4TrackingManager pValue)
virtualinherited

Reimplemented in G4MultiTrackingAction.

Definition at line 63 of file G4UserTrackingAction.cc.

65{
66 fpTrackingManager = pValue;
67}
G4TrackingManager * fpTrackingManager

References G4UserTrackingAction::fpTrackingManager.

Referenced by G4TrackingManager::SetUserAction().

◆ SetUserForwardTrackingAction()

void G4AdjointTrackingAction::SetUserForwardTrackingAction ( G4UserTrackingAction anAction)
inline

Definition at line 62 of file G4AdjointTrackingAction.hh.

63 { theUserFwdTrackingAction = anAction; }

References theUserFwdTrackingAction.

Referenced by G4AdjointSimManager::SetAdjointActions().

Field Documentation

◆ fpTrackingManager

G4TrackingManager* G4UserTrackingAction::fpTrackingManager = nullptr
protectedinherited

◆ is_adjoint_tracking_mode

G4bool G4AdjointTrackingAction::is_adjoint_tracking_mode = false
private

◆ last_cos_th

G4double G4AdjointTrackingAction::last_cos_th = 0.0
private

Definition at line 101 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_cos_th_vec

std::vector<G4double> G4AdjointTrackingAction::last_cos_th_vec
private

◆ last_direction

G4ThreeVector G4AdjointTrackingAction::last_direction
private

Definition at line 98 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_direction_vec

std::vector<G4ThreeVector> G4AdjointTrackingAction::last_direction_vec
private

◆ last_ekin

G4double G4AdjointTrackingAction::last_ekin = 0.0
private

Definition at line 99 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_ekin_nuc

G4double G4AdjointTrackingAction::last_ekin_nuc = 0.0
private

Definition at line 99 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_ekin_nuc_vec

std::vector<G4double> G4AdjointTrackingAction::last_ekin_nuc_vec
private

◆ last_ekin_vec

std::vector<G4double> G4AdjointTrackingAction::last_ekin_vec
private

◆ last_fwd_part_index

G4int G4AdjointTrackingAction::last_fwd_part_index = 0
private

Definition at line 105 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_fwd_part_index_vec

std::vector<G4int> G4AdjointTrackingAction::last_fwd_part_index_vec
private

◆ last_fwd_part_name

G4String G4AdjointTrackingAction::last_fwd_part_name
private

◆ last_fwd_part_PDGEncoding

G4int G4AdjointTrackingAction::last_fwd_part_PDGEncoding = 0
private

Definition at line 103 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_fwd_part_PDGEncoding_vec

std::vector<G4int> G4AdjointTrackingAction::last_fwd_part_PDGEncoding_vec
private

◆ last_pos

G4ThreeVector G4AdjointTrackingAction::last_pos
private

Definition at line 97 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_pos_vec

std::vector<G4ThreeVector> G4AdjointTrackingAction::last_pos_vec
private

◆ last_weight

G4double G4AdjointTrackingAction::last_weight = 0.0
private

Definition at line 104 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction().

◆ last_weight_vec

std::vector<G4double> G4AdjointTrackingAction::last_weight_vec
private

◆ pListOfPrimaryFwdParticles

std::vector<G4ParticleDefinition*>* G4AdjointTrackingAction::pListOfPrimaryFwdParticles = nullptr
private

◆ theAdjointSteppingAction

G4AdjointSteppingAction* G4AdjointTrackingAction::theAdjointSteppingAction = nullptr
private

Definition at line 91 of file G4AdjointTrackingAction.hh.

Referenced by PostUserTrackingAction(), and PreUserTrackingAction().

◆ theUserFwdTrackingAction

G4UserTrackingAction* G4AdjointTrackingAction::theUserFwdTrackingAction = nullptr
private

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