Geant4-11
Data Structures | Public Member Functions | Static Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes
G4DNAIndependentReactionTimeStepper Class Reference

#include <G4DNAIndependentReactionTimeStepper.hh>

Inheritance diagram for G4DNAIndependentReactionTimeStepper:
G4VITTimeStepComputer

Data Structures

class  Utils
 

Public Member Functions

G4double CalculateMinTimeStep (G4double, G4double) override
 
G4double CalculateStep (const G4Track &, const G4double &) override
 
std::unique_ptr< G4ITReactionChangeFindReaction (G4ITReactionSet *pReactionSet, const G4double &currentStepTime=0, const G4double &previousStepTime=0, const G4bool &reachedUserStepTimeLimit=false)
 
 G4DNAIndependentReactionTimeStepper ()
 
 G4DNAIndependentReactionTimeStepper (const G4DNAIndependentReactionTimeStepper &)=delete
 
G4TrackVectorHandle GetReactants ()
 
G4VDNAReactionModelGetReactionModel ()
 
const G4ITReactionTableGetReactionTable ()
 
G4double GetSampledMinTimeStep ()
 
virtual void Initialize ()
 
G4DNAIndependentReactionTimeStepperoperator= (const G4DNAIndependentReactionTimeStepper &)=delete
 
void Prepare () override
 
virtual void ResetReactants ()
 
void SetReactionModel (G4VDNAReactionModel *)
 
void SetReactionProcess (G4VITReactionProcess *pReactionProcess)
 
void SetReactionTable (const G4ITReactionTable *)
 
void SetReactionTypeManager (G4VReactionTypeManager *typeManager)
 
void SetVerbose (G4int)
 
 ~G4DNAIndependentReactionTimeStepper () override=default
 

Static Public Member Functions

static void SetTimes (const G4double &, const G4double &)
 

Protected Attributes

const G4ITReactionTablefpReactionTable
 
G4TrackVectorHandle fReactants
 
G4double fSampledMinTimeStep
 

Static Protected Attributes

static G4ThreadLocal G4double fCurrentGlobalTime = -1
 
static G4ThreadLocal G4double fUserMinTimeStep = -1
 

Private Member Functions

void CheckAndRecordResults (const Utils &utils)
 
ReactionType GetReactionType (const G4Track &trackA, const G4Track &trackB)
 
G4double GetTimeToEncounter (const G4Track &trackA, const G4Track &trackB)
 
void InitializeForNewTrack ()
 

Private Attributes

G4bool fHasAlreadyReachedNullTime
 
const G4DNAMolecularReactionTable *& fMolecularReactionTable
 
G4VITReactionProcessfpReactionProcess
 
G4ITTrackHolderfpTrackContainer
 
G4double fRCutOff
 
G4VDNAReactionModelfReactionModel
 
G4ITReactionSetfReactionSet
 
G4DNAReactionTypeManagerfReactionTypeManager
 
std::map< G4int, G4ThreeVectorfSampledPositions
 
G4int fVerbose
 

Detailed Description

Definition at line 51 of file G4DNAIndependentReactionTimeStepper.hh.

Constructor & Destructor Documentation

◆ G4DNAIndependentReactionTimeStepper() [1/2]

G4DNAIndependentReactionTimeStepper::G4DNAIndependentReactionTimeStepper ( )

Definition at line 57 of file G4DNAIndependentReactionTimeStepper.cc.

60 , fMolecularReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
61 , fReactionModel(nullptr)
64 , fVerbose(0)
66 , fReactionTypeManager(nullptr)
67 , fpReactionProcess(nullptr)
68{
70}
const G4DNAMolecularReactionTable *& fMolecularReactionTable
static G4double GetRCutOff()
Definition: G4IRTUtils.cc:39
static G4ITReactionSet * Instance()
static G4ITTrackHolder * Instance()
const G4ITReactionTable * fpReactionTable

References fReactionSet, and G4ITReactionSet::SortByTime().

◆ ~G4DNAIndependentReactionTimeStepper()

G4DNAIndependentReactionTimeStepper::~G4DNAIndependentReactionTimeStepper ( )
overridedefault

◆ G4DNAIndependentReactionTimeStepper() [2/2]

G4DNAIndependentReactionTimeStepper::G4DNAIndependentReactionTimeStepper ( const G4DNAIndependentReactionTimeStepper )
delete

Member Function Documentation

◆ CalculateMinTimeStep()

G4double G4DNAIndependentReactionTimeStepper::CalculateMinTimeStep ( G4double  ,
G4double  definedMinTimeStep 
)
overridevirtual

Implements G4VITTimeStepComputer.

Definition at line 439 of file G4DNAIndependentReactionTimeStepper.cc.

440{
441 G4double fTSTimeStep = DBL_MAX;
442
443 for (auto pTrack : *fpTrackContainer->GetMainList())
444 {
445 if (pTrack == nullptr)
446 {
447 G4ExceptionDescription exceptionDescription;
448 exceptionDescription << "No track found.";
449 G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler006",
450 FatalErrorInArgument, exceptionDescription);
451 continue;
452 }
453
454 G4TrackStatus trackStatus = pTrack->GetTrackStatus();
455 if (trackStatus == fStopAndKill || trackStatus == fStopButAlive)
456 {
457 continue;
458 }
459
460 G4double sampledMinTimeStep = CalculateStep(*pTrack, definedMinTimeStep);
461 G4TrackVectorHandle reactants = GetReactants();
462
463 if (sampledMinTimeStep < fTSTimeStep)
464 {
465 fTSTimeStep = sampledMinTimeStep;
467 if (reactants)
468 {
469 fReactionSet->AddReactions(fTSTimeStep,
470 const_cast<G4Track*>(pTrack),
471 reactants);
473 }
474 }
475 else if (fTSTimeStep == sampledMinTimeStep && G4bool(reactants))
476 {
477 fReactionSet->AddReactions(fTSTimeStep,
478 const_cast<G4Track*>(pTrack),
479 reactants);
481 }
482 else if (reactants)
483 {
485 }
486 }
487
488 return fTSTimeStep;
489}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4shared_ptr< std::vector< G4Track * > > G4TrackVectorHandle
Definition: G4ITReaction.hh:44
G4TrackStatus
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4double CalculateStep(const G4Track &, const G4double &) override
void CleanAllReaction()
void AddReactions(G4double time, G4Track *trackA, G4TrackVectorHandle reactants)
G4TrackList * GetMainList(Key)
G4TrackVectorHandle GetReactants()
#define DBL_MAX
Definition: templates.hh:62

References G4ITReactionSet::AddReactions(), CalculateStep(), G4ITReactionSet::CleanAllReaction(), DBL_MAX, FatalErrorInArgument, fpTrackContainer, fReactionSet, fStopAndKill, fStopButAlive, G4Exception(), G4ITTrackHolder::GetMainList(), G4VITTimeStepComputer::GetReactants(), and G4VITTimeStepComputer::ResetReactants().

◆ CalculateStep()

G4double G4DNAIndependentReactionTimeStepper::CalculateStep ( const G4Track trackA,
const G4double userMinTimeStep 
)
overridevirtual

Implements G4VITTimeStepComputer.

Definition at line 96 of file G4DNAIndependentReactionTimeStepper.cc.

98{
99 auto pMoleculeA = GetMolecule(trackA);
101 fUserMinTimeStep = userMinTimeStep;
102
103#ifdef G4VERBOSE
104 if (fVerbose)
105 {
106 G4cout
107 << "_______________________________________________________________________"
108 << G4endl;
109 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep" << G4endl;
110 G4cout << "Check done for molecule : " << pMoleculeA->GetName()
111 << " (" << trackA.GetTrackID() << ") "
112 << G4endl;
113 }
114#endif
115
116 auto pMolConfA = pMoleculeA->GetMolecularConfiguration();
117
118 const auto pReactantList = fMolecularReactionTable->CanReactWith(pMolConfA);
119
120 if (!pReactantList)
121 {
122#ifdef G4VERBOSE
123 if (fVerbose > 1)
124 {
125 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
126 G4cout << "!!! WARNING" << G4endl;
127 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will return infinity "
128 "for the reaction because the molecule "
129 << pMoleculeA->GetName()
130 << " does not have any reactants given in the reaction table."
131 << G4endl;
132 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
133 }
134#endif
135 return DBL_MAX;
136 }
137
138 G4int nbReactives = pReactantList->size();
139
140 if (nbReactives == 0)
141 {
142#ifdef G4VERBOSE
143 // DEBUG
144 if (fVerbose)
145 {
146 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
147 G4cout << "!!! WARNING" << G4endl;
148 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will return infinity "
149 "for the reaction because the molecule "
150 << pMoleculeA->GetName()
151 << " does not have any reactants given in the reaction table."
152 << "This message can also result from a wrong implementation of the reaction table."
153 << G4endl;
154 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
155 }
156#endif
157 return DBL_MAX;
158 }
159 fReactants.reset(new vector<G4Track*>());
160 fReactionModel->Initialise(pMolConfA, trackA);
161 for (G4int i = 0; i < nbReactives; i++)
162 {
163 auto pMoleculeB = (*pReactantList)[i];
164 G4int key = pMoleculeB->GetMoleculeID();
165
166 //fRCutOff = G4IRTUtils::GetRCutOff(1 * ps);
168 //______________________________________________________________
169 // Retrieve reaction range
171 std::vector<std::pair<G4TrackList::iterator,G4double>> resultIndices;
172 resultIndices.clear();
174 FindNearestInRange(trackA,
175 key,
176 fRCutOff,
177 resultIndices);
178
179 if(resultIndices.empty())
180 {
181 continue;
182 }
183 for(auto& it : resultIndices)
184 {
185 G4Track* pTrackB = *(std::get<0>(it));
186
187 if(pTrackB == &trackA)
188 {
189 continue;
190 }
191 if(pTrackB == nullptr)
192 {
193 G4ExceptionDescription exceptionDescription;
194 exceptionDescription << "No trackB no valid";
195 G4Exception("G4DNAIndependentReactionTimeModel"
196 "::BuildReactionMap()", "NO_TRACK02",
197 FatalException, exceptionDescription);
198 }
199 Utils utils(trackA, *pTrackB);
200
201 auto pMolB = GetMolecule(pTrackB);
202 auto pMolConfB = pMolB->GetMolecularConfiguration();
203 G4double distance = (trackA.GetPosition() - pTrackB->GetPosition()).mag();
204 if(distance * distance < Reff * Reff)
205 {
206 auto processTable = *(fReactionTypeManager->GetReactionTypeTable());
207 auto typeOfReaction = (G4int)GetReactionType(trackA, *pTrackB);
208 if(processTable[typeOfReaction]->
209 GeminateRecombinationProbability(pMolConfA, pMolConfB))
210 {
212 {
213 fReactants->clear();
215 }
218 }
219 }
220 else
221 {
222 G4double tempMinET = GetTimeToEncounter(trackA, *pTrackB);
223 if(tempMinET < 0 ||
224 tempMinET > G4Scheduler::Instance()->GetEndTime())
225 {
226 continue;
227 }
228 if (tempMinET >= fSampledMinTimeStep)
229 {
230 continue;
231 }
232 fSampledMinTimeStep = tempMinET;
233 fReactants->clear();
235 }
236 }
237 }
238
239#ifdef G4VERBOSE
240 if (fVerbose)
241 {
242 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will finally return :"
244
245 if (fVerbose > 1)
246 {
247 G4cout << "Selected reactants for trackA: " << pMoleculeA->GetName()
248 << " (" << trackA.GetTrackID() << ") are: ";
249
250 vector<G4Track*>::iterator it;
251 for (it = fReactants->begin(); it != fReactants->end(); it++)
252 {
253 G4Track* trackB = *it;
254 G4cout << GetMolecule(trackB)->GetName() << " ("
255 << trackB->GetTrackID() << ") \t ";
256 }
257 G4cout << G4endl;
258 }
259 }
260#endif
261
262 for (const auto& it : *fReactants)
263 {
264 auto pTrackB = it;
265 fSampledPositions[pTrackB->GetTrackID()] = pTrackB->GetPosition();
266 }
267 return fSampledMinTimeStep;
268}
@ FatalException
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
#define G4BestUnit(a, b)
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetTimeToEncounter(const G4Track &trackA, const G4Track &trackB)
ReactionType GetReactionType(const G4Track &trackA, const G4Track &trackB)
const ReactantList * CanReactWith(Reactant *) const
const ReactionTypeTable * GetReactionTypeTable() const override
const G4String & GetName() const
Definition: G4Molecule.cc:338
static G4OctreeFinder * Instance()
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
virtual void Initialise(const G4MolecularConfiguration *, const G4Track &)
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0
G4TrackVectorHandle fReactants
static G4ThreadLocal G4double fUserMinTimeStep

References G4DNAMolecularReactionTable::CanReactWith(), CheckAndRecordResults(), DBL_MAX, FatalException, fHasAlreadyReachedNullTime, fMolecularReactionTable, fRCutOff, G4VITTimeStepComputer::fReactants, fReactionModel, fReactionTypeManager, G4VITTimeStepComputer::fSampledMinTimeStep, fSampledPositions, G4VITTimeStepComputer::fUserMinTimeStep, fVerbose, G4BestUnit, G4cout, G4endl, G4Exception(), GetMolecule(), G4Molecule::GetName(), G4Track::GetPosition(), G4IRTUtils::GetRCutOff(), G4VDNAReactionModel::GetReactionRadius(), GetReactionType(), G4DNAReactionTypeManager::GetReactionTypeTable(), GetTimeToEncounter(), G4Track::GetTrackID(), G4VDNAReactionModel::Initialise(), InitializeForNewTrack(), G4OctreeFinder< T, CONTAINER >::Instance(), and G4Scheduler::Instance().

Referenced by CalculateMinTimeStep().

◆ CheckAndRecordResults()

void G4DNAIndependentReactionTimeStepper::CheckAndRecordResults ( const Utils utils)
private

Definition at line 270 of file G4DNAIndependentReactionTimeStepper.cc.

271{
272 if (utils.fTrackB.GetTrackStatus() != fAlive)
273 {
274 return;
275 }
276
277 if (&utils.fTrackB == &utils.fTrackA)
278 {
279 G4ExceptionDescription exceptionDescription;
280 exceptionDescription<< "A track is reacting with itself"
281 " (which is impossible) ie fpTrackA == trackB"<< G4endl;
282 exceptionDescription << "Molecule A is of type : "
283 << utils.fpMoleculeA->GetName() << " with trackID : "
284 << utils.fTrackA.GetTrackID()<<" and B : "
285 << utils.fpMoleculeB->GetName() << " with trackID : "
286 << utils.fTrackB.GetTrackID() << G4endl;
287 G4Exception("G4DNAIndependentReactionTimeStepper::RetrieveResults",
288 "G4DNAIndependentReactionTimeStepper003", FatalErrorInArgument,
289 exceptionDescription);
290 }
291
292 if (fabs(utils.fTrackB.GetGlobalTime() - utils.fTrackA.GetGlobalTime())
293 > utils.fTrackA.GetGlobalTime() * (1. - 1. / 100))
294 {
295 // DEBUG
296 G4ExceptionDescription exceptionDescription;
297 exceptionDescription
298 << "The interacting tracks are not synchronized in time" << G4endl;
299 exceptionDescription
300 << "trackB->GetGlobalTime() != fpTrackA.GetGlobalTime()" << G4endl;
301
302 exceptionDescription << "fpTrackA : trackID : " << utils.fTrackA.GetTrackID()
303 << "\t Name :" << utils.fpMoleculeA->GetName()
304 << "\t fpTrackA->GetGlobalTime() = "
305 << G4BestUnit(utils.fTrackA.GetGlobalTime(), "Time") << G4endl;
306
307 exceptionDescription << "trackB : trackID : " << utils.fTrackB.GetTrackID()
308 << "\t Name :" << utils.fpMoleculeB->GetName()
309 << "\t trackB->GetGlobalTime() = "
310 << G4BestUnit(utils.fTrackB.GetGlobalTime(), "Time") << G4endl;
311
312 G4Exception("G4DNAIndependentReactionTimeStepper::RetrieveResults",
313 "G4DNAIndependentReactionTimeStepper004", FatalErrorInArgument,
314 exceptionDescription);
315 }
316 fReactants->push_back(const_cast<G4Track*>(&utils.fTrackB));
317}
@ fAlive

References fAlive, FatalErrorInArgument, G4DNAIndependentReactionTimeStepper::Utils::fpMoleculeA, G4DNAIndependentReactionTimeStepper::Utils::fpMoleculeB, G4VITTimeStepComputer::fReactants, G4DNAIndependentReactionTimeStepper::Utils::fTrackA, G4DNAIndependentReactionTimeStepper::Utils::fTrackB, G4BestUnit, G4endl, G4Exception(), G4Track::GetGlobalTime(), G4Molecule::GetName(), G4Track::GetTrackID(), and G4Track::GetTrackStatus().

Referenced by CalculateStep().

◆ FindReaction()

std::unique_ptr< G4ITReactionChange > G4DNAIndependentReactionTimeStepper::FindReaction ( G4ITReactionSet pReactionSet,
const G4double currentStepTime = 0,
const G4double previousStepTime = 0,
const G4bool reachedUserStepTimeLimit = false 
)

Definition at line 320 of file G4DNAIndependentReactionTimeStepper.cc.

324{
325 if (pReactionSet == nullptr)
326 {
327 return nullptr;
328 }
329
330 G4ITReactionPerTime& reactionPerTime = pReactionSet->GetReactionsPerTime();
331 if(reactionPerTime.empty())
332 {
333 return nullptr;
334 }
335
336 for (auto reaction_i = reactionPerTime.begin();
337 reaction_i != reactionPerTime.end();
338 reaction_i = reactionPerTime.begin())
339 {
340 G4Track* pTrackA = (*reaction_i)->GetReactants().first;
341 if (pTrackA->GetTrackStatus() == fStopAndKill)
342 {
343 continue;
344 }
345 G4Track* pTrackB = (*reaction_i)->GetReactant(pTrackA);
346 if (pTrackB->GetTrackStatus() == fStopAndKill)
347 {
348 continue;
349 }
350
351 if (pTrackB == pTrackA)
352 {
353 G4ExceptionDescription exceptionDescription;
354 exceptionDescription
355 << "The IT reaction process sent back a reaction between trackA and trackB. ";
356 exceptionDescription << "The problem is trackA == trackB";
357 G4Exception("G4ITModelProcessor::FindReaction",
358 "ITModelProcessor005",
360 exceptionDescription);
361 }
362 pReactionSet->SelectThisReaction(*reaction_i);
363 if(fpReactionProcess != nullptr && fpReactionProcess->TestReactibility(*pTrackA,
364 *pTrackB,
365 currentStepTime,
366 false))
367 {
368 pTrackA->SetPosition(fSampledPositions[pTrackA->GetTrackID()]);
369 pTrackB->SetPosition(fSampledPositions[pTrackB->GetTrackID()]);
370 auto pReactionChange = fpReactionProcess->MakeReaction(*pTrackA, *pTrackB);
371 if (pReactionChange == nullptr)
372 {
373 return nullptr;
374 }
375 return pReactionChange;
376 }
377 }
378 return nullptr;
379}
std::multiset< G4ITReactionPtr, compReactionPerTime > G4ITReactionPerTime
Definition: G4ITReaction.hh:77
void SelectThisReaction(G4ITReactionPtr reaction)
G4ITReactionPerTime & GetReactionsPerTime()
G4TrackStatus GetTrackStatus() const
void SetPosition(const G4ThreeVector &aValue)
virtual std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &)=0
virtual G4bool TestReactibility(const G4Track &, const G4Track &, double, bool)=0

References FatalErrorInArgument, fpReactionProcess, fSampledPositions, fStopAndKill, G4Exception(), G4ITReactionSet::GetReactionsPerTime(), G4Track::GetTrackID(), G4Track::GetTrackStatus(), G4VITReactionProcess::MakeReaction(), G4ITReactionSet::SelectThisReaction(), G4Track::SetPosition(), and G4VITReactionProcess::TestReactibility().

◆ GetReactants()

G4TrackVectorHandle G4VITTimeStepComputer::GetReactants ( )
inlineinherited

◆ GetReactionModel()

G4VDNAReactionModel * G4DNAIndependentReactionTimeStepper::GetReactionModel ( )

Definition at line 386 of file G4DNAIndependentReactionTimeStepper.cc.

387{
388 return fReactionModel;
389}

References fReactionModel.

◆ GetReactionTable()

const G4ITReactionTable * G4VITTimeStepComputer::GetReactionTable ( )
inlineinherited

Definition at line 123 of file G4VITTimeStepComputer.hh.

124{
125 return fpReactionTable ;
126}

References G4VITTimeStepComputer::fpReactionTable.

◆ GetReactionType()

ReactionType G4DNAIndependentReactionTimeStepper::GetReactionType ( const G4Track trackA,
const G4Track trackB 
)
private

Definition at line 396 of file G4DNAIndependentReactionTimeStepper.cc.

398{
399 auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
400 auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
401 auto pData = fMolecularReactionTable->GetReactionData(pMoleculeA,pMoleculeB);
402 G4int reactionID = pData->GetReactionID();
403 return fReactionTypeManager->GetReactionTypeByID(reactionID);
404}
Data * GetReactionData(Reactant *, Reactant *) const
ReactionType GetReactionTypeByID(ReactionID iD)
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532

References fMolecularReactionTable, fReactionTypeManager, G4Molecule::GetMolecularConfiguration(), GetMolecule(), G4DNAMolecularReactionTable::GetReactionData(), G4DNAMolecularReactionData::GetReactionID(), and G4DNAReactionTypeManager::GetReactionTypeByID().

Referenced by CalculateStep(), and GetTimeToEncounter().

◆ GetSampledMinTimeStep()

G4double G4VITTimeStepComputer::GetSampledMinTimeStep ( )
inlineinherited

Definition at line 134 of file G4VITTimeStepComputer.hh.

135{
136 return fSampledMinTimeStep ;
137}

References G4VITTimeStepComputer::fSampledMinTimeStep.

◆ GetTimeToEncounter()

G4double G4DNAIndependentReactionTimeStepper::GetTimeToEncounter ( const G4Track trackA,
const G4Track trackB 
)
private

Definition at line 406 of file G4DNAIndependentReactionTimeStepper.cc.

408{
409 if(fReactionTypeManager == nullptr)
410 {
411 G4ExceptionDescription exceptionDescription;
412 exceptionDescription << "fpProManager is not "
413 "initialized ";
414 G4Exception("G4DNAIndependentReactionTimeModel::"
415 "GetIndependentReactionTime()",
416 "G4DNAIndependentReactionTimeModel002",
417 FatalErrorInArgument,exceptionDescription);
418 }
419 auto processTable = *(fReactionTypeManager->GetReactionTypeTable());
420 ReactionType reactionType = GetReactionType(trackA,trackB);
421
422#ifdef DEBUG
423 G4cout<<"A: "<<GetMolecule(trackA)->GetName()<<"("<<trackA.GetTrackID()<<")"<<" + B : "
424 <<GetMolecule(trackB)->GetName()<<"("<<trackB.GetTrackID()<<")"<<G4endl;
425#endif
426
427 return processTable[(G4int)reactionType]->GetTimeToEncounter(trackA,trackB);
428}
ReactionType

References FatalErrorInArgument, fReactionTypeManager, G4cout, G4endl, G4Exception(), GetMolecule(), G4Molecule::GetName(), GetReactionType(), G4DNAReactionTypeManager::GetReactionTypeTable(), GetTimeToEncounter(), and G4Track::GetTrackID().

Referenced by CalculateStep(), and GetTimeToEncounter().

◆ Initialize()

virtual void G4VITTimeStepComputer::Initialize ( )
inlinevirtualinherited

This macro defined in AddClone_def

Definition at line 83 of file G4VITTimeStepComputer.hh.

83{;}

◆ InitializeForNewTrack()

void G4DNAIndependentReactionTimeStepper::InitializeForNewTrack ( )
private

◆ operator=()

G4DNAIndependentReactionTimeStepper & G4DNAIndependentReactionTimeStepper::operator= ( const G4DNAIndependentReactionTimeStepper )
delete

◆ Prepare()

void G4DNAIndependentReactionTimeStepper::Prepare ( )
overridevirtual

◆ ResetReactants()

virtual void G4VITTimeStepComputer::ResetReactants ( )
inlinevirtualinherited

◆ SetReactionModel()

void G4DNAIndependentReactionTimeStepper::SetReactionModel ( G4VDNAReactionModel pReactionModel)

Definition at line 381 of file G4DNAIndependentReactionTimeStepper.cc.

382{
383 fReactionModel = pReactionModel;
384}

References fReactionModel.

◆ SetReactionProcess()

void G4DNAIndependentReactionTimeStepper::SetReactionProcess ( G4VITReactionProcess pReactionProcess)

Definition at line 435 of file G4DNAIndependentReactionTimeStepper.cc.

436{
437 fpReactionProcess = pReactionProcess;
438}

References fpReactionProcess.

◆ SetReactionTable()

void G4VITTimeStepComputer::SetReactionTable ( const G4ITReactionTable table)
inlineinherited

Definition at line 118 of file G4VITTimeStepComputer.hh.

119{
120 fpReactionTable = table;
121}

References G4VITTimeStepComputer::fpReactionTable.

◆ SetReactionTypeManager()

void G4DNAIndependentReactionTimeStepper::SetReactionTypeManager ( G4VReactionTypeManager typeManager)

◆ SetTimes()

void G4VITTimeStepComputer::SetTimes ( const G4double currentGlobalTime,
const G4double userMinStepTime 
)
staticinherited

Definition at line 68 of file G4VITTimeStepComputer.cc.

70{
71 fCurrentGlobalTime = currentGlobalTime ;
72 fUserMinTimeStep = userMinStepTime ;
73}
static G4ThreadLocal G4double fCurrentGlobalTime

References G4VITTimeStepComputer::fCurrentGlobalTime, and G4VITTimeStepComputer::fUserMinTimeStep.

Referenced by G4ITModelProcessor::InitializeStepper().

◆ SetVerbose()

void G4DNAIndependentReactionTimeStepper::SetVerbose ( G4int  flag)

Definition at line 391 of file G4DNAIndependentReactionTimeStepper.cc.

392{
393 fVerbose = flag;
394}

References fVerbose.

Field Documentation

◆ fCurrentGlobalTime

G4ThreadLocal G4double G4VITTimeStepComputer::fCurrentGlobalTime = -1
staticprotectedinherited

Definition at line 106 of file G4VITTimeStepComputer.hh.

Referenced by G4VITTimeStepComputer::SetTimes().

◆ fHasAlreadyReachedNullTime

G4bool G4DNAIndependentReactionTimeStepper::fHasAlreadyReachedNullTime
private

Definition at line 87 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateStep(), and InitializeForNewTrack().

◆ fMolecularReactionTable

const G4DNAMolecularReactionTable*& G4DNAIndependentReactionTimeStepper::fMolecularReactionTable
private

Definition at line 88 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateStep(), and GetReactionType().

◆ fpReactionProcess

G4VITReactionProcess* G4DNAIndependentReactionTimeStepper::fpReactionProcess
private

Definition at line 96 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by FindReaction(), and SetReactionProcess().

◆ fpReactionTable

const G4ITReactionTable* G4VITTimeStepComputer::fpReactionTable
protectedinherited

◆ fpTrackContainer

G4ITTrackHolder* G4DNAIndependentReactionTimeStepper::fpTrackContainer
private

Definition at line 90 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateMinTimeStep().

◆ fRCutOff

G4double G4DNAIndependentReactionTimeStepper::fRCutOff
private

Definition at line 93 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateStep().

◆ fReactants

G4TrackVectorHandle G4VITTimeStepComputer::fReactants
protectedinherited

◆ fReactionModel

G4VDNAReactionModel* G4DNAIndependentReactionTimeStepper::fReactionModel
private

◆ fReactionSet

G4ITReactionSet* G4DNAIndependentReactionTimeStepper::fReactionSet
private

◆ fReactionTypeManager

G4DNAReactionTypeManager* G4DNAIndependentReactionTimeStepper::fReactionTypeManager
private

◆ fSampledMinTimeStep

G4double G4VITTimeStepComputer::fSampledMinTimeStep
protectedinherited

◆ fSampledPositions

std::map<G4int,G4ThreeVector> G4DNAIndependentReactionTimeStepper::fSampledPositions
private

Definition at line 97 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateStep(), FindReaction(), and Prepare().

◆ fUserMinTimeStep

G4ThreadLocal G4double G4VITTimeStepComputer::fUserMinTimeStep = -1
staticprotectedinherited

◆ fVerbose

G4int G4DNAIndependentReactionTimeStepper::fVerbose
private

Definition at line 92 of file G4DNAIndependentReactionTimeStepper.hh.

Referenced by CalculateStep(), and SetVerbose().


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