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

#include <G4DNAIndependentReactionTimeStepper.hh>

+ Inheritance diagram for G4DNAIndependentReactionTimeStepper:

Public Member Functions

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

Additional Inherited Members

- Static Public Member Functions inherited from G4VITTimeStepComputer
static void SetTimes (const G4double &, const G4double &)
 
- Protected Attributes inherited from G4VITTimeStepComputer
G4double fSampledMinTimeStep
 
G4TrackVectorHandle fReactants
 
const G4ITReactionTablefpReactionTable
 
- Static Protected Attributes inherited from G4VITTimeStepComputer
static G4ThreadLocal G4double fCurrentGlobalTime = -1
 
static G4ThreadLocal G4double fUserMinTimeStep = -1
 

Detailed Description

Definition at line 51 of file G4DNAIndependentReactionTimeStepper.hh.

Constructor & Destructor Documentation

◆ G4DNAIndependentReactionTimeStepper() [1/2]

G4DNAIndependentReactionTimeStepper::G4DNAIndependentReactionTimeStepper ( )

Definition at line 54 of file G4DNAIndependentReactionTimeStepper.cc.

55{
56 fReactionSet->SortByTime();
57}

Referenced by G4DNAIndependentReactionTimeStepper(), and operator=().

◆ ~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 357 of file G4DNAIndependentReactionTimeStepper.cc.

359{
360 G4double fTSTimeStep = DBL_MAX;
361 fCheckedTracks.clear();
362
363 for (auto pTrack : *fpTrackContainer->GetMainList()) {
364 if (pTrack == nullptr) {
366 msg << "No track found.";
367 G4Exception("G4DNAIndependentReactionTimeStepper::CalculateMinTimeStep",
368 "G4DNAIndependentReactionTimeStepper006", FatalErrorInArgument,
369 msg);
370 continue;
371 }
372
373 G4TrackStatus trackStatus = pTrack->GetTrackStatus();
374 if (trackStatus == fStopAndKill || trackStatus == fStopButAlive) {
375 continue;
376 }
377
378 G4double sampledMinTimeStep = CalculateStep(*pTrack, definedMinTimeStep);
379 G4TrackVectorHandle reactants = GetReactants();
380
381 if (sampledMinTimeStep < fTSTimeStep) {
382 fTSTimeStep = sampledMinTimeStep;
383 if (reactants) {
384 fReactionSet->AddReactions(fTSTimeStep, const_cast<G4Track*>(pTrack), std::move(reactants));
385
386 fSampledPositions[pTrack->GetTrackID()] = pTrack->GetPosition();
387 for (const auto& it : *fReactants) {
388 auto pTrackB = it;
389 fSampledPositions[pTrackB->GetTrackID()] = pTrackB->GetPosition();
390 }
392 }
393 }
394 else if (fTSTimeStep == sampledMinTimeStep && G4bool(reactants)) {
395 fReactionSet->AddReactions(fTSTimeStep, const_cast<G4Track*>(pTrack), std::move(reactants));
396
397 fSampledPositions[pTrack->GetTrackID()] = pTrack->GetPosition();
398 for (const auto& it : *fReactants) {
399 auto pTrackB = it;
400 fSampledPositions[pTrackB->GetTrackID()] = pTrackB->GetPosition();
401 }
403 }
404 else if (reactants) {
406 }
407 }
408 return fTSTimeStep;
409}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4shared_ptr< std::vector< G4Track * > > G4TrackVectorHandle
G4TrackStatus
@ fStopAndKill
@ fStopButAlive
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double CalculateStep(const G4Track &, const G4double &) override
G4TrackVectorHandle fReactants
G4TrackVectorHandle GetReactants()
#define DBL_MAX
Definition templates.hh:62

◆ CalculateStep()

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

Implements G4VITTimeStepComputer.

Definition at line 75 of file G4DNAIndependentReactionTimeStepper.cc.

77{
78 auto pMoleculeA = GetMolecule(trackA);
79 InitializeForNewTrack();
80 fUserMinTimeStep = userMinTimeStep;
81 fCheckedTracks.insert(trackA.GetTrackID());
82
83#ifdef G4VERBOSE
84 if (fVerbose != 0) {
85 G4cout << "________________________________________________________________"
86 "_______"
87 << G4endl;
88 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep" << G4endl;
89 G4cout << "Check done for molecule : " << pMoleculeA->GetName() << " (" << trackA.GetTrackID()
90 << ") " << G4endl;
91 }
92#endif
93
94 auto pMolConfA = pMoleculeA->GetMolecularConfiguration();
95
96 const auto pReactantList = fMolecularReactionTable->CanReactWith(pMolConfA);
97
98 if (pReactantList == nullptr) {
99 if(fVerbose > 1) {
101 msg << "G4DNAIndependentReactionTimeStepper::CalculateStep will return infinity "
102 "for the reaction because the molecule "
103 << pMoleculeA->GetName() << " does not have any reactants given in the reaction table."
104 << G4endl;
105 G4Exception("G4DNAIndependentReactionTimeStepper::CalculateStep",
106 "G4DNAIndependentReactionTimeStepper03", JustWarning, msg);
107 }
108 return DBL_MAX;
109 }
110
111 auto nbReactives = (G4int)pReactantList->size();
112
113 if (nbReactives == 0) {
114 if(fVerbose != 0){
116 msg << "G4DNAIndependentReactionTimeStepper::CalculateStep will "
117 "return infinity "
118 "for the reaction because the molecule "
119 << pMoleculeA->GetName() << " does not have any reactants given in the reaction table."
120 << "This message can also result from a wrong implementation of "
121 "the reaction table."
122 << G4endl;
123 G4Exception("G4DNAIndependentReactionTimeStepper::CalculateStep",
124 "G4DNAIndependentReactionTimeStepper04", JustWarning, msg);
125 }
126 return DBL_MAX;
127 }
128 fReactants = std::make_shared<vector<G4Track*>>();
129 fReactionModel->Initialise(pMolConfA, trackA);
130 for (G4int i = 0; i < nbReactives; ++i) {
131 auto pMoleculeB = (*pReactantList)[i];
132 G4int key = pMoleculeB->GetMoleculeID();
133
134 // fRCutOff = G4IRTUtils::GetRCutOff(1 * ps);
135 fRCutOff = G4IRTUtils::GetRCutOff();
136 //______________________________________________________________
137 // Retrieve reaction range
138 const G4double Reff = fReactionModel->GetReactionRadius(i);
139 std::vector<std::pair<G4TrackList::iterator, G4double>> resultIndices;
140 resultIndices.clear();
141 G4ChemicalMoleculeFinder::Instance()->FindNearestInRange(trackA, key, fRCutOff, resultIndices);
142
143 if (resultIndices.empty()) {
144 continue;
145 }
146 for (auto& it : resultIndices) {
147 G4Track* pTrackB = *(std::get<0>(it));
148
149 if (pTrackB == &trackA) {
150 continue;
151 }
152 if (pTrackB == nullptr) {
153 G4ExceptionDescription exceptionDescription;
154 exceptionDescription << "No trackB no valid";
156 "G4DNAIndependentReactionTimeStepper"
157 "::CalculateStep()",
158 "G4DNAIndependentReactionTimeStepper007", FatalException, exceptionDescription);
159 }
160 else {
161 if (fCheckedTracks.find(pTrackB->GetTrackID()) != fCheckedTracks.end()) {
162 continue;
163 }
164
165 Utils utils(trackA, *pTrackB);
166
167 auto pMolB = GetMolecule(pTrackB);
168 auto pMolConfB = pMolB->GetMolecularConfiguration();
169 G4double distance = (trackA.GetPosition() - pTrackB->GetPosition()).mag();
170 if (distance * distance < Reff * Reff) {
171 auto reactionData = fMolecularReactionTable->GetReactionData(pMolConfA, pMolConfB);
172 if (G4Scheduler::Instance()->GetGlobalTime() == G4Scheduler::Instance()->GetStartTime()) {
173 if (reactionData->GetProbability() > G4UniformRand()) {
174 if (!fHasAlreadyReachedNullTime) {
175 fReactants->clear();
176 fHasAlreadyReachedNullTime = true;
177 }
179 CheckAndRecordResults(utils);
180 }
181 }
182 }
183 else {
184 G4double tempMinET = GetTimeToEncounter(trackA, *pTrackB);
185 if (tempMinET < 0 || tempMinET > G4Scheduler::Instance()->GetEndTime()) {
186 continue;
187 }
188 if (tempMinET >= fSampledMinTimeStep) {
189 continue;
190 }
191 fSampledMinTimeStep = tempMinET;
192 fReactants->clear();
193 CheckAndRecordResults(utils);
194 }
195 }
196 }
197 }
198
199#ifdef G4VERBOSE
200 if (fVerbose != 0) {
201 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will finally "
202 "return :"
204
205 if (fVerbose > 1) {
206 G4cout << "Selected reactants for trackA: " << pMoleculeA->GetName() << " ("
207 << trackA.GetTrackID() << ") are: ";
208
209 vector<G4Track*>::iterator it;
210 for (it = fReactants->begin(); it != fReactants->end(); it++) {
211 G4Track* trackB = *it;
212 G4cout << GetMolecule(trackB)->GetName() << " (" << trackB->GetTrackID() << ") \t ";
213 }
214 G4cout << G4endl;
215 }
216 }
217#endif
218 return fSampledMinTimeStep;
219}
@ JustWarning
@ FatalException
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
#define G4BestUnit(a, b)
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
static G4double GetRCutOff()
Definition G4IRTUtils.cc:39
const G4String & GetName() const override
static G4OctreeFinder * Instance()
void FindNearestInRange(const G4Track &track, const int &key, G4double R, std::vector< std::pair< typename CONTAINER::iterator, G4double > > &result, G4bool isSort=false) const
static G4Scheduler * Instance()
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
static G4ThreadLocal G4double fUserMinTimeStep

Referenced by CalculateMinTimeStep().

◆ FindReaction()

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

Definition at line 266 of file G4DNAIndependentReactionTimeStepper.cc.

269{
270 if (pReactionSet == nullptr) {
271 return nullptr;
272 }
273
274 G4ITReactionPerTime& reactionPerTime = pReactionSet->GetReactionsPerTime();
275 if (reactionPerTime.empty()) {
276 return nullptr;
277 }
278 for (auto reaction_i = reactionPerTime.begin(); reaction_i != reactionPerTime.end();
279 reaction_i = reactionPerTime.begin())
280 {
281 if ((*reaction_i)->GetTime() > currentStepTime) {
282 fReactionSet->CleanAllReaction();
283 return nullptr;
284 }
285
286 G4Track* pTrackA = (*reaction_i)->GetReactants().first;
287 if (pTrackA->GetTrackStatus() == fStopAndKill) {
288 continue;
289 }
290 G4Track* pTrackB = (*reaction_i)->GetReactant(pTrackA);
291 if (pTrackB->GetTrackStatus() == fStopAndKill) {
292 continue;
293 }
294
295 if (pTrackB == pTrackA) {
297 msg << "The IT reaction process sent back a reaction "
298 "between trackA and trackB. ";
299 msg << "The problem is trackA == trackB";
300 G4Exception("G4DNAIndependentReactionTimeStepper::FindReaction",
301 "G4DNAIndependentReactionTimeStepper02", FatalErrorInArgument,
302 msg);
303 }
304 pReactionSet->SelectThisReaction(*reaction_i);
305 if (fpReactionProcess != nullptr
306 && fpReactionProcess->TestReactibility(*pTrackA, *pTrackB, currentStepTime, false))
307 {
308 if ((fSampledPositions.find(pTrackA->GetTrackID()) == fSampledPositions.end()
309 && (fSampledPositions.find(pTrackB->GetTrackID()) == fSampledPositions.end())))
310 {
312 msg << "The positions of trackA and trackB have no counted ";
313 G4Exception("G4DNAIndependentReactionTimeStepper::FindReaction",
314 "G4DNAIndependentReactionTimeStepper0001", FatalErrorInArgument,
315 msg);
316 }
317
318 pTrackA->SetPosition(fSampledPositions[pTrackA->GetTrackID()]);
319 pTrackB->SetPosition(fSampledPositions[pTrackB->GetTrackID()]);
320 auto pReactionChange = fpReactionProcess->MakeReaction(*pTrackA, *pTrackB);
321 if (pReactionChange == nullptr) {
322 return nullptr;
323 }
324 return pReactionChange;
325 }
326 }
327 return nullptr;
328}
std::multiset< G4ITReactionPtr, compReactionPerTime > G4ITReactionPerTime
void SelectThisReaction(G4ITReactionPtr reaction)
G4ITReactionPerTime & GetReactionsPerTime()
G4TrackStatus GetTrackStatus() const
void SetPosition(const G4ThreeVector &aValue)

◆ GetReactionModel()

G4VDNAReactionModel * G4DNAIndependentReactionTimeStepper::GetReactionModel ( )

Definition at line 335 of file G4DNAIndependentReactionTimeStepper.cc.

336{
337 return fReactionModel;
338}

◆ operator=()

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

◆ Prepare()

void G4DNAIndependentReactionTimeStepper::Prepare ( )
overridevirtual

Reimplemented from G4VITTimeStepComputer.

Definition at line 59 of file G4DNAIndependentReactionTimeStepper.cc.

60{
62 fSampledPositions.clear();
64}
#define BuildChemicalMoleculeFinder()

◆ SetReactionModel()

void G4DNAIndependentReactionTimeStepper::SetReactionModel ( G4VDNAReactionModel * pReactionModel)

Definition at line 330 of file G4DNAIndependentReactionTimeStepper.cc.

331{
332 fReactionModel = pReactionModel;
333}

◆ SetReactionProcess()

void G4DNAIndependentReactionTimeStepper::SetReactionProcess ( G4VITReactionProcess * pReactionProcess)

Definition at line 353 of file G4DNAIndependentReactionTimeStepper.cc.

354{
355 fpReactionProcess = pReactionProcess;
356}

◆ SetVerbose()

void G4DNAIndependentReactionTimeStepper::SetVerbose ( G4int flag)

Definition at line 340 of file G4DNAIndependentReactionTimeStepper.cc.

341{
342 fVerbose = flag;
343}

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