Geant4 11.2.2
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
void AddReaction(G4double time, G4Track *trackA, G4Track *trackB)
G4bool CanAddThisReaction(G4Track *trackA, G4Track *trackB)
static G4ITReactionPtr New(G4double time, G4Track *trackA, G4Track *trackB)

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

◆ AddReaction() [2/2]

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

Definition at line 318 of file G4ITReaction.hh.

319 {
320 auto it = fReactionPerTrack.find(track);
321
322 G4ITReactionPerTrackPtr reactionPerTrack;
323
324 if(it == fReactionPerTrack.end())
325 {
326 reactionPerTrack = G4ITReactionPerTrack::New();
327 std::pair< G4ITReactionPerTrackMap::iterator,bool> pos =
328 fReactionPerTrack.insert(std::make_pair(track, reactionPerTrack));
329 reactionPerTrack->AddIterator(pos.first);
330 }
331 else
332 {
333 reactionPerTrack = it->second;
334 }
335
336 reactionPerTrack->AddReaction(reaction);
337 }
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
static G4ITReactionPerTrackPtr New()

◆ AddReactions()

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

Definition at line 244 of file G4ITReaction.hh.

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

Referenced by G4DNAIndependentReactionTimeStepper::CalculateMinTimeStep(), and G4DNAMoleculeEncounterStepper::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
231 reactionPerTrack = it_track->second;
232 auto list = reactionPerTrack->GetReactionList();
233 //for(auto it_list = list.begin(); it_list != list.end(); ++it_list)
234 for(const auto& it_list:list)
235 {
236 if ((*it_list).GetReactant(trackA)->GetTrackID() == trackB->GetTrackID())
237 {
238 return false;
239 }
240 }
241 return true;
242 }
G4int GetTrackID() const

Referenced by AddReaction().

◆ CleanAllReaction()

void G4ITReactionSet::CleanAllReaction ( )
inline

◆ Empty()

◆ GetReactionMap()

G4ITReactionPerTrackMap & G4ITReactionSet::GetReactionMap ( )
inline

Definition at line 276 of file G4ITReaction.hh.

277 {
278 return fReactionPerTrack;
279 }

Referenced by G4DNAMolecularReaction::FindReaction().

◆ GetReactionsPerTime()

◆ Instance()

◆ RemoveReactionPerTrack()

void G4ITReactionSet::RemoveReactionPerTrack ( G4ITReactionPerTrackPtr reactionPerTrack)
inline

Definition at line 281 of file G4ITReaction.hh.

282 {
283 for(auto & it : reactionPerTrack->GetListOfIterators())
284 {
285 fReactionPerTrack.erase(it);
286 }
287 reactionPerTrack->GetListOfIterators().clear();
288 reactionPerTrack->GetReactionList().clear();
289 }

Referenced by G4ITReactionPerTrack::RemoveThisReaction().

◆ RemoveReactionSet()

void G4ITReactionSet::RemoveReactionSet ( G4Track * track)
inline

Definition at line 253 of file G4ITReaction.hh.

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

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

◆ SelectThisReaction()

void G4ITReactionSet::SelectThisReaction ( G4ITReactionPtr reaction)
inline

Definition at line 269 of file G4ITReaction.hh.

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

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

◆ SortByTime()

void G4ITReactionSet::SortByTime ( )
inline

Member Data Documentation

◆ fpInstance

G4ThreadLocal G4ITReactionSet * G4ITReactionSet::fpInstance
staticprotected

Definition at line 342 of file G4ITReaction.hh.

Referenced by Instance().

◆ fReactionPerTime

G4ITReactionPerTime G4ITReactionSet::fReactionPerTime
protected

◆ fReactionPerTrack

◆ fSortByTime

G4bool G4ITReactionSet::fSortByTime
protected

Definition at line 341 of file G4ITReaction.hh.

Referenced by AddReaction(), and SortByTime().


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