Geant4 11.1.1
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 (G4double time, G4Track *trackA, G4Track *trackB)
 
G4bool CanAddThisReaction (G4Track *trackA, G4Track *trackB)
 
void AddReactions (G4double time, G4Track *trackA, G4TrackVectorHandle reactants)
 
void RemoveReactionSet (G4Track *track)
 
void SelectThisReaction (G4ITReactionPtr reaction)
 
G4ITReactionPerTrackMapGetReactionMap ()
 
void RemoveReactionPerTrack (G4ITReactionPerTrackPtr reactionPerTrack)
 
void CleanAllReaction ()
 
G4bool 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
 
G4bool fSortByTime
 

Static Protected Attributes

static G4ThreadLocal G4ITReactionSetfpInstance
 

Detailed Description

Definition at line 182 of file G4ITReaction.hh.

Constructor & Destructor Documentation

◆ ~G4ITReactionSet()

virtual G4ITReactionSet::~G4ITReactionSet ( )
inlinevirtual

Definition at line 190 of file G4ITReaction.hh.

191 {
192 fReactionPerTrack.clear();
193 fReactionPerTime.clear();
194 }
G4ITReactionPerTime fReactionPerTime
G4ITReactionPerTrackMap fReactionPerTrack

Member Function Documentation

◆ AddReaction() [1/2]

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

Definition at line 205 of file G4ITReaction.hh.

206 {
207 if(CanAddThisReaction(trackA, trackB))
208 {
209 G4ITReactionPtr reaction(G4ITReaction::New(time, trackA, trackB));
210 AddReaction(trackA, reaction);
211 AddReaction(trackB, reaction);
212
213 if(fSortByTime)
214 {
215 auto it = fReactionPerTime.insert(reaction);
216 reaction->AddIterator(it);
217 }
218 }
219 }
G4shared_ptr< G4ITReaction > G4ITReactionPtr
Definition: G4ITReaction.hh:61
void AddReaction(G4double time, G4Track *trackA, G4Track *trackB)
G4bool CanAddThisReaction(G4Track *trackA, G4Track *trackB)
static G4ITReactionPtr New(G4double time, G4Track *trackA, G4Track *trackB)
Definition: G4ITReaction.hh:84

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

◆ AddReaction() [2/2]

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

Definition at line 323 of file G4ITReaction.hh.

324 {
325 auto it = fReactionPerTrack.find(track);
326
327 G4ITReactionPerTrackPtr reactionPerTrack;
328
329 if(it == fReactionPerTrack.end())
330 {
331 reactionPerTrack = G4ITReactionPerTrack::New();
332 std::pair< G4ITReactionPerTrackMap::iterator,bool> pos =
333 fReactionPerTrack.insert(std::make_pair(track, reactionPerTrack));
334 reactionPerTrack->AddIterator(pos.first);
335 }
336 else
337 {
338 reactionPerTrack = it->second;
339 }
340
341 reactionPerTrack->AddReaction(reaction);
342 }
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
Definition: G4ITReaction.hh:62
static G4ITReactionPerTrackPtr New()

◆ AddReactions()

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

Definition at line 246 of file G4ITReaction.hh.

247 {
248 auto it = reactants->begin();
249 for(;it != reactants->end() ; ++it)
250 {
251 AddReaction(time, trackA, *it);
252 }
253 }

Referenced by G4DNAMoleculeEncounterStepper::CalculateMinTimeStep(), and G4DNAIndependentReactionTimeStepper::CalculateMinTimeStep().

◆ CanAddThisReaction()

G4bool G4ITReactionSet::CanAddThisReaction ( G4Track trackA,
G4Track trackB 
)
inline

Definition at line 222 of file G4ITReaction.hh.

223 {
224 auto it_track = fReactionPerTrack.find(trackA);
225 G4ITReactionPerTrackPtr reactionPerTrack;
226 if(it_track == fReactionPerTrack.end())
227 {
228 return true;
229 }
230 else
231 {
232 reactionPerTrack = it_track->second;
233 auto list = reactionPerTrack->GetReactionList();
234 //for(auto it_list = list.begin(); it_list != list.end(); ++it_list)
235 for(const auto& it_list:list)
236 {
237 if ((*it_list).GetReactant(trackA)->GetTrackID() == trackB->GetTrackID())
238 {
239 return false;
240 }
241 }
242 return true;
243 }
244 }
G4int GetTrackID() const

Referenced by AddReaction().

◆ CleanAllReaction()

void G4ITReactionSet::CleanAllReaction ( )
inline

◆ Empty()

◆ GetReactionMap()

G4ITReactionPerTrackMap & G4ITReactionSet::GetReactionMap ( )
inline

Definition at line 278 of file G4ITReaction.hh.

279 {
280 return fReactionPerTrack;
281 }

Referenced by G4DNAMolecularReaction::FindReaction().

◆ GetReactionsPerTime()

◆ Instance()

◆ RemoveReactionPerTrack()

void G4ITReactionSet::RemoveReactionPerTrack ( G4ITReactionPerTrackPtr  reactionPerTrack)
inline

Definition at line 283 of file G4ITReaction.hh.

284 {
285 for(auto it =
286 reactionPerTrack->GetListOfIterators().begin() ;
287 it != reactionPerTrack->GetListOfIterators().end() ;
288 ++it)
289 {
290 fReactionPerTrack.erase(*it);
291 }
292 reactionPerTrack->GetListOfIterators().clear();
293 reactionPerTrack->GetReactionList().clear();
294 }

Referenced by G4ITReactionPerTrack::RemoveThisReaction().

◆ RemoveReactionSet()

void G4ITReactionSet::RemoveReactionSet ( G4Track track)
inline

Definition at line 255 of file G4ITReaction.hh.

256 {
257 auto it = fReactionPerTrack.find(track);
258 if(it != fReactionPerTrack.end())
259 {
260 G4ITReactionPerTrackPtr backItUp = it->second->shared_from_this();
261 backItUp->RemoveMe();
262 //fReactionPerTrack.erase(it); // not needed : once empty ==> auto-erase
263 it = fReactionPerTrack.find(track);
264 if(it != fReactionPerTrack.end())
265 {
266 fReactionPerTrack.erase(it);
267 }
268 }
269 }

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

◆ SelectThisReaction()

void G4ITReactionSet::SelectThisReaction ( G4ITReactionPtr  reaction)
inline

Definition at line 271 of file G4ITReaction.hh.

272 {
273 reaction->RemoveMe();
274 RemoveReactionSet(reaction->GetReactants().first);
275 RemoveReactionSet(reaction->GetReactants().second);
276 }
void RemoveReactionSet(G4Track *track)

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

◆ SortByTime()

void G4ITReactionSet::SortByTime ( )
inline

Member Data Documentation

◆ fpInstance

G4ThreadLocal G4ITReactionSet * G4ITReactionSet::fpInstance
staticprotected

Definition at line 347 of file G4ITReaction.hh.

Referenced by Instance().

◆ fReactionPerTime

G4ITReactionPerTime G4ITReactionSet::fReactionPerTime
protected

◆ fReactionPerTrack

◆ fSortByTime

G4bool G4ITReactionSet::fSortByTime
protected

Definition at line 346 of file G4ITReaction.hh.

Referenced by AddReaction(), and SortByTime().


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