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

#include <G4ITReaction.hh>

+ Inheritance diagram for G4ITReaction:

Public Member Functions

virtual ~G4ITReaction ()
 
G4TrackGetReactant (G4Track *trackA)
 
std::pair< G4Track *, G4Track * > GetReactants () const
 
std::size_t GetHash () const
 
double GetTime ()
 
void RemoveMe ()
 
void AddIterator (G4ITReactionPerTrackPtr reactionPerTrack, G4ITReactionList::iterator it)
 
void AddIterator (G4ITReactionPerTimeIt it)
 

Static Public Member Functions

static G4ITReactionPtr New (double time, G4Track *trackA, G4Track *trackB)
 

Public Attributes

double fTime
 
std::pair< G4Track *, G4Track * > fReactants
 
G4ReactionPerTrackIt fReactionPerTrack
 
G4ITReactionPerTimeItfReactionPerTimeIt
 

Detailed Description

Definition at line 79 of file G4ITReaction.hh.

Constructor & Destructor Documentation

◆ ~G4ITReaction()

G4ITReaction::~G4ITReaction ( )
virtual

Definition at line 80 of file G4ITReaction.cc.

81{
82 //gAll->erase(this);
84}
G4ITReactionPerTimeIt * fReactionPerTimeIt

Member Function Documentation

◆ AddIterator() [1/2]

void G4ITReaction::AddIterator ( G4ITReactionPerTimeIt  it)
inline

Definition at line 107 of file G4ITReaction.hh.

108 {
110 }
std::multiset< G4ITReactionPtr, compReactionPerTime >::iterator G4ITReactionPerTimeIt
Definition: G4ITReaction.hh:77

◆ AddIterator() [2/2]

void G4ITReaction::AddIterator ( G4ITReactionPerTrackPtr  reactionPerTrack,
G4ITReactionList::iterator  it 
)
inline

Definition at line 101 of file G4ITReaction.hh.

103 {
104 fReactionPerTrack.push_back(std::make_pair(reactionPerTrack, it));
105 }
G4ReactionPerTrackIt fReactionPerTrack

◆ GetHash()

std::size_t G4ITReaction::GetHash ( ) const

Definition at line 61 of file G4ITReaction.cc.

62{
63 std::size_t hash = 0;
64 hash_combine(hash, fReactants.first->GetTrackID());
65 hash_combine(hash, fReactants.first->GetTrackID());
66 return hash;
67}
void hash_combine(std::size_t &seed, const T &v)
Definition: G4ITReaction.cc:54
std::pair< G4Track *, G4Track * > fReactants

◆ GetReactant()

G4Track * G4ITReaction::GetReactant ( G4Track trackA)
inline

Definition at line 89 of file G4ITReaction.hh.

90 {
91 if(fReactants.first != trackA) return fReactants.first;
92 return fReactants.second;
93 }

◆ GetReactants()

std::pair< G4Track *, G4Track * > G4ITReaction::GetReactants ( ) const
inline

Definition at line 95 of file G4ITReaction.hh.

95{return fReactants;}

◆ GetTime()

double G4ITReaction::GetTime ( )
inline

Definition at line 97 of file G4ITReaction.hh.

97{ return fTime; }

◆ New()

static G4ITReactionPtr G4ITReaction::New ( double  time,
G4Track trackA,
G4Track trackB 
)
inlinestatic

Definition at line 83 of file G4ITReaction.hh.

84 {
85 return G4ITReactionPtr(new G4ITReaction(time, trackA, trackB));
86 }
G4shared_ptr< G4ITReaction > G4ITReactionPtr
Definition: G4ITReaction.hh:60

Referenced by G4ITReactionSet::AddReaction().

◆ RemoveMe()

void G4ITReaction::RemoveMe ( )

Definition at line 86 of file G4ITReaction.cc.

87{
88 G4ITReactionPtr backMeUp = this->shared_from_this();
89 for(G4ReactionPerTrackIt::iterator it = fReactionPerTrack.begin() ;
90 it != fReactionPerTrack.end() ; ++it)
91 {
92 // G4cout << it->first.get() << G4endl;
93 // assert(it->first.get() != 0);
94 it->first->RemoveThisReaction(it->second);
95 }
96 fReactionPerTrack.clear();
97
99 {
101 delete fReactionPerTimeIt;
103 }
104}
static G4ITReactionSet * Instance()
G4ITReactionPerTime & GetReactionsPerTime()

Member Data Documentation

◆ fReactants

std::pair<G4Track*, G4Track*> G4ITReaction::fReactants

Definition at line 113 of file G4ITReaction.hh.

Referenced by GetHash(), GetReactant(), and GetReactants().

◆ fReactionPerTimeIt

G4ITReactionPerTimeIt* G4ITReaction::fReactionPerTimeIt

Definition at line 115 of file G4ITReaction.hh.

Referenced by AddIterator(), RemoveMe(), and ~G4ITReaction().

◆ fReactionPerTrack

G4ReactionPerTrackIt G4ITReaction::fReactionPerTrack

Definition at line 114 of file G4ITReaction.hh.

Referenced by AddIterator(), and RemoveMe().

◆ fTime

double G4ITReaction::fTime

Definition at line 112 of file G4ITReaction.hh.

Referenced by GetTime().


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