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

#include <G4ITReaction.hh>

Public Member Functions

virtual ~G4ITReactionSet ()
 
void AddReaction (double time, G4Track *trackA, G4Track *trackB)
 
void AddReactions (double time, G4Track *trackA, G4TrackVectorHandle reactants)
 
void RemoveReactionSet (G4Track *track)
 
void SelectThisReaction (G4ITReactionPtr reaction)
 
G4ITReactionPerTrackMapGetReactionMap ()
 
void RemoveReactionPerTrack (G4ITReactionPerTrackPtr reactionPerTrack)
 
void CleanAllReaction ()
 
bool Empty ()
 
G4ITReactionPerTimeGetReactionsPerTime ()
 
void SortByTime ()
 

Static Public Member Functions

static G4ITReactionSetInstance ()
 

Protected Member Functions

void AddReaction (G4Track *track, G4ITReactionPtr reaction)
 

Protected Attributes

G4ITReactionPerTrackMap fReactionPerTrack
 
G4ITReactionPerTime fReactionPerTime
 
bool fSortByTime
 

Static Protected Attributes

static G4ThreadLocal G4ITReactionSetfpInstance
 

Detailed Description

Definition at line 178 of file G4ITReaction.hh.

Constructor & Destructor Documentation

◆ ~G4ITReactionSet()

virtual G4ITReactionSet::~G4ITReactionSet ( )
inlinevirtual

Definition at line 186 of file G4ITReaction.hh.

187 {
188 fReactionPerTrack.clear();
189 fReactionPerTime.clear();
190 }
G4ITReactionPerTime fReactionPerTime
G4ITReactionPerTrackMap fReactionPerTrack

Member Function Documentation

◆ AddReaction() [1/2]

void G4ITReactionSet::AddReaction ( double  time,
G4Track trackA,
G4Track trackB 
)
inline

Definition at line 201 of file G4ITReaction.hh.

202 {
203 G4ITReactionPtr reaction(G4ITReaction::New(time, trackA, trackB));
204 AddReaction(trackA, reaction);
205 AddReaction(trackB, reaction);
206
207 if(fSortByTime)
208 {
209 G4ITReactionPerTime::iterator it = fReactionPerTime.insert(reaction);
210 reaction->AddIterator(it);
211 }
212 }
G4shared_ptr< G4ITReaction > G4ITReactionPtr
Definition: G4ITReaction.hh:60
void AddReaction(double time, G4Track *trackA, G4Track *trackB)
static G4ITReactionPtr New(double time, G4Track *trackA, G4Track *trackB)
Definition: G4ITReaction.hh:83

Referenced by AddReaction(), AddReactions(), and G4DNAIRT::Sampling().

◆ AddReaction() [2/2]

void G4ITReactionSet::AddReaction ( G4Track track,
G4ITReactionPtr  reaction 
)
inlineprotected

Definition at line 291 of file G4ITReaction.hh.

292 {
293 G4ITReactionPerTrackMap::iterator it = fReactionPerTrack.find(track);
294
295 G4ITReactionPerTrackPtr reactionPerTrack;
296
297 if(it == fReactionPerTrack.end())
298 {
299 reactionPerTrack = G4ITReactionPerTrack::New();
300 std::pair< G4ITReactionPerTrackMap::iterator,bool> pos =
301 fReactionPerTrack.insert(std::make_pair(track, reactionPerTrack));
302 reactionPerTrack->AddIterator(pos.first);
303 }
304 else
305 {
306 reactionPerTrack = it->second;
307 }
308
309 reactionPerTrack->AddReaction(reaction);
310 }
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
Definition: G4ITReaction.hh:61
static G4ITReactionPerTrackPtr New()

◆ AddReactions()

void G4ITReactionSet::AddReactions ( double  time,
G4Track trackA,
G4TrackVectorHandle  reactants 
)
inline

Definition at line 214 of file G4ITReaction.hh.

215 {
216 std::vector<G4Track*>::iterator it = reactants->begin();
217 for(;it != reactants->end() ; ++it)
218 {
219 AddReaction(time, trackA, *it);
220 }
221 }

Referenced by G4DNAMoleculeEncounterStepper::CalculateMinTimeStep().

◆ CleanAllReaction()

void G4ITReactionSet::CleanAllReaction ( )
inline

Definition at line 264 of file G4ITReaction.hh.

265 {
266 for(G4ITReactionPerTrackMap::iterator it = fReactionPerTrack.begin();
267 it != fReactionPerTrack.end() ;
268 it = fReactionPerTrack.begin())
269 {
270 it->second->RemoveMe();
271 }
272 fReactionPerTrack.clear();
273 fReactionPerTime.clear();
274 }

Referenced by G4DNAMoleculeEncounterStepper::CalculateMinTimeStep(), G4DNAMolecularReaction::FindReaction(), G4DNAIRT::Initialize(), G4Scheduler::Reset(), and G4Scheduler::Stepping().

◆ Empty()

bool G4ITReactionSet::Empty ( )
inline

◆ GetReactionMap()

G4ITReactionPerTrackMap & G4ITReactionSet::GetReactionMap ( )
inline

Definition at line 246 of file G4ITReaction.hh.

247 {
248 return fReactionPerTrack;
249 }

Referenced by G4DNAMolecularReaction::FindReaction().

◆ GetReactionsPerTime()

G4ITReactionPerTime & G4ITReactionSet::GetReactionsPerTime ( )
inline

◆ Instance()

◆ RemoveReactionPerTrack()

void G4ITReactionSet::RemoveReactionPerTrack ( G4ITReactionPerTrackPtr  reactionPerTrack)
inline

Definition at line 251 of file G4ITReaction.hh.

252 {
253 for(std::list<G4ITReactionPerTrackMap::iterator>::iterator it =
254 reactionPerTrack->GetListOfIterators().begin() ;
255 it != reactionPerTrack->GetListOfIterators().end() ;
256 ++it)
257 {
258 fReactionPerTrack.erase(*it);
259 }
260 reactionPerTrack->GetListOfIterators().clear();
261 reactionPerTrack->GetReactionList().clear();
262 }

Referenced by G4ITReactionPerTrack::RemoveThisReaction().

◆ RemoveReactionSet()

void G4ITReactionSet::RemoveReactionSet ( G4Track track)
inline

Definition at line 223 of file G4ITReaction.hh.

224 {
225 G4ITReactionPerTrackMap::iterator it = fReactionPerTrack.find(track);
226 if(it != fReactionPerTrack.end())
227 {
228 G4ITReactionPerTrackPtr backItUp = it->second->shared_from_this();
229 backItUp->RemoveMe();
230 //fReactionPerTrack.erase(it); // not needed : once empty ==> auto-erase
231 it = fReactionPerTrack.find(track);
232 if(it != fReactionPerTrack.end())
233 {
234 fReactionPerTrack.erase(it);
235 }
236 }
237 }

Referenced by G4ITStepProcessor::ExtractDoItData(), and SelectThisReaction().

◆ SelectThisReaction()

void G4ITReactionSet::SelectThisReaction ( G4ITReactionPtr  reaction)
inline

Definition at line 239 of file G4ITReaction.hh.

240 {
241 reaction->RemoveMe();
242 RemoveReactionSet(reaction->GetReactants().first);
243 RemoveReactionSet(reaction->GetReactants().second);
244 }
void RemoveReactionSet(G4Track *track)

Referenced by G4DNAIRT::FindReaction(), and G4DNAMolecularReaction::FindReaction().

◆ SortByTime()

void G4ITReactionSet::SortByTime ( )
inline

Definition at line 286 of file G4ITReaction.hh.

286 {
287 fSortByTime = true;
288 }

Referenced by G4DNAIRT::Initialize().

Member Data Documentation

◆ fpInstance

G4ThreadLocal G4ITReactionSet * G4ITReactionSet::fpInstance
staticprotected

Definition at line 315 of file G4ITReaction.hh.

Referenced by Instance().

◆ fReactionPerTime

G4ITReactionPerTime G4ITReactionSet::fReactionPerTime
protected

◆ fReactionPerTrack

G4ITReactionPerTrackMap G4ITReactionSet::fReactionPerTrack
protected

◆ fSortByTime

bool G4ITReactionSet::fSortByTime
protected

Definition at line 314 of file G4ITReaction.hh.

Referenced by AddReaction(), and SortByTime().


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