Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAIndependentReactionTimeStepper.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25// 20/2/2019
26// Author: HoangTRAN
27
28#include <memory>
32#include "G4UnitsTable.hh"
33#include "G4Molecule.hh"
36#include "G4DNAMakeReaction.hh"
37#include "G4ITReactionChange.hh"
38#include "G4Scheduler.hh"
39#include "G4IRTUtils.hh"
40#include "Randomize.hh"
42using namespace std;
43using namespace CLHEP;
44
45G4DNAIndependentReactionTimeStepper::Utils::Utils(const G4Track& trackA,
46 const G4Track& trackB)
47 : fTrackA(trackA)
48 , fTrackB(trackB)
49{
50 fpMoleculeA = GetMolecule(trackA);
51 fpMoleculeB = GetMolecule(trackA);
52}
53
56{
57 fReactionSet->SortByTime();
58}
59
61{
63 fSampledPositions.clear();
65}
66
67void G4DNAIndependentReactionTimeStepper::InitializeForNewTrack()
68{
69 if(fReactants != nullptr)
70 {
71 fReactants.reset();
72 }
74 fHasAlreadyReachedNullTime = false;
75}
76
78 const G4Track& trackA, const G4double& userMinTimeStep)
79{
80 auto pMoleculeA = GetMolecule(trackA);
81 InitializeForNewTrack();
82 fUserMinTimeStep = userMinTimeStep;
83 fCheckedTracks.insert(trackA.GetTrackID());
84
85#ifdef G4VERBOSE
86 if(fVerbose != 0)
87 {
88 G4cout << "________________________________________________________________"
89 "_______"
90 << G4endl;
91 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep" << G4endl;
92 G4cout << "Check done for molecule : " << pMoleculeA->GetName() << " ("
93 << trackA.GetTrackID() << ") " << G4endl;
94 }
95#endif
96
97 auto pMolConfA = pMoleculeA->GetMolecularConfiguration();
98
99 const auto pReactantList = fMolecularReactionTable->CanReactWith(pMolConfA);
100
101 if(pReactantList == nullptr)
102 {
103#ifdef G4VERBOSE
104 if(fVerbose > 1)
105 {
106 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
107 G4cout << "!!! WARNING" << G4endl;
108 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will "
109 "return infinity "
110 "for the reaction because the molecule "
111 << pMoleculeA->GetName()
112 << " does not have any reactants given in the reaction table."
113 << G4endl;
114 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
115 }
116#endif
117 return DBL_MAX;
118 }
119
120 G4int nbReactives = (G4int)pReactantList->size();
121
122 if(nbReactives == 0)
123 {
124#ifdef G4VERBOSE
125 // DEBUG
126 if(fVerbose != 0)
127 {
128 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
129 G4cout << "!!! WARNING" << G4endl;
130 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will "
131 "return infinity "
132 "for the reaction because the molecule "
133 << pMoleculeA->GetName()
134 << " does not have any reactants given in the reaction table."
135 << "This message can also result from a wrong implementation of "
136 "the reaction table."
137 << G4endl;
138 G4cout << "!!!!!!!!!!!!!!!!!!!!" << G4endl;
139 }
140#endif
141 return DBL_MAX;
142 }
143 fReactants = std::make_shared<vector<G4Track*>>();
144 fReactionModel->Initialise(pMolConfA, trackA);
145 for(G4int i = 0; i < nbReactives; ++i)
146 {
147 auto pMoleculeB = (*pReactantList)[i];
148 G4int key = pMoleculeB->GetMoleculeID();
149
150 // fRCutOff = G4IRTUtils::GetRCutOff(1 * ps);
151 fRCutOff = G4IRTUtils::GetRCutOff();
152 //______________________________________________________________
153 // Retrieve reaction range
154 const G4double Reff = fReactionModel->GetReactionRadius(i);
155 std::vector<std::pair<G4TrackList::iterator, G4double>> resultIndices;
156 resultIndices.clear();
158 trackA, key, fRCutOff, resultIndices);
159
160 if(resultIndices.empty())
161 {
162 continue;
163 }
164 for(auto& it : resultIndices)
165 {
166 G4Track* pTrackB = *(std::get<0>(it));
167
168 if(pTrackB == &trackA)
169 {
170 continue;
171 }
172 if(pTrackB == nullptr)
173 {
174 G4ExceptionDescription exceptionDescription;
175 exceptionDescription << "No trackB no valid";
176 G4Exception("G4DNAIndependentReactionTimeStepper"
177 "::CalculateStep()",
178 "G4DNAIndependentReactionTimeStepper007", FatalException,
179 exceptionDescription);
180 }else
181 {
182 if(fCheckedTracks.find(pTrackB->GetTrackID()) != fCheckedTracks.end())
183 {
184 continue;
185 }
186
187 Utils utils(trackA, *pTrackB);
188
189 auto pMolB = GetMolecule(pTrackB);
190 auto pMolConfB = pMolB->GetMolecularConfiguration();
191 G4double distance = (trackA.GetPosition() - pTrackB->GetPosition()).mag();
192 if(distance * distance < Reff * Reff)
193 {
194 auto reactionData =
195 fMolecularReactionTable->GetReactionData(pMolConfA, pMolConfB);
197 {
198 if(reactionData->GetProbability() > G4UniformRand())
199 {
200 if(!fHasAlreadyReachedNullTime)
201 {
202 fReactants->clear();
203 fHasAlreadyReachedNullTime = true;
204 }
206 CheckAndRecordResults(utils);
207 }
208 }
209 }
210 else
211 {
212 G4double tempMinET = GetTimeToEncounter(trackA, *pTrackB);
213 if(tempMinET < 0 || tempMinET > G4Scheduler::Instance()->GetEndTime())
214 {
215 continue;
216 }
217 if(tempMinET >= fSampledMinTimeStep)
218 {
219 continue;
220 }
221 fSampledMinTimeStep = tempMinET;
222 fReactants->clear();
223 CheckAndRecordResults(utils);
224 }
225 }
226 }
227 }
228
229#ifdef G4VERBOSE
230 if(fVerbose != 0)
231 {
232 G4cout << "G4DNAIndependentReactionTimeStepper::CalculateStep will finally "
233 "return :"
235
236 if(fVerbose > 1)
237 {
238 G4cout << "Selected reactants for trackA: " << pMoleculeA->GetName()
239 << " (" << trackA.GetTrackID() << ") are: ";
240
241 vector<G4Track*>::iterator it;
242 for(it = fReactants->begin(); it != fReactants->end(); it++)
243 {
244 G4Track* trackB = *it;
245 G4cout << GetMolecule(trackB)->GetName() << " (" << trackB->GetTrackID()
246 << ") \t ";
247 }
248 G4cout << G4endl;
249 }
250 }
251#endif
252 return fSampledMinTimeStep;
253}
254
255void G4DNAIndependentReactionTimeStepper::CheckAndRecordResults(
256 const Utils& utils)
257{
258 if(utils.fTrackB.GetTrackStatus() != fAlive)
259 {
260 return;
261 }
262
263 if(&utils.fTrackB == &utils.fTrackA)
264 {
265 G4ExceptionDescription exceptionDescription;
266 exceptionDescription << "A track is reacting with itself"
267 " (which is impossible) ie fpTrackA == trackB"
268 << G4endl;
269 exceptionDescription << "Molecule A is of type : "
270 << utils.fpMoleculeA->GetName()
271 << " with trackID : " << utils.fTrackA.GetTrackID()
272 << " and B : " << utils.fpMoleculeB->GetName()
273 << " with trackID : " << utils.fTrackB.GetTrackID()
274 << G4endl;
275 G4Exception("G4DNAIndependentReactionTimeStepper::RetrieveResults",
276 "G4DNAIndependentReactionTimeStepper003", FatalErrorInArgument,
277 exceptionDescription);
278 }
279
280 if(fabs(utils.fTrackB.GetGlobalTime() - utils.fTrackA.GetGlobalTime()) >
281 utils.fTrackA.GetGlobalTime() * (1. - 1. / 100))
282 {
283 // DEBUG
284 G4ExceptionDescription exceptionDescription;
285 exceptionDescription
286 << "The interacting tracks are not synchronized in time" << G4endl;
287 exceptionDescription
288 << "trackB->GetGlobalTime() != fpTrackA.GetGlobalTime()" << G4endl;
289
290 exceptionDescription << "fpTrackA : trackID : "
291 << utils.fTrackA.GetTrackID()
292 << "\t Name :" << utils.fpMoleculeA->GetName()
293 << "\t fpTrackA->GetGlobalTime() = "
294 << G4BestUnit(utils.fTrackA.GetGlobalTime(), "Time")
295 << G4endl;
296
297 exceptionDescription << "trackB : trackID : " << utils.fTrackB.GetTrackID()
298 << "\t Name :" << utils.fpMoleculeB->GetName()
299 << "\t trackB->GetGlobalTime() = "
300 << G4BestUnit(utils.fTrackB.GetGlobalTime(), "Time")
301 << G4endl;
302
303 G4Exception("G4DNAIndependentReactionTimeStepper::RetrieveResults",
304 "G4DNAIndependentReactionTimeStepper004", FatalErrorInArgument,
305 exceptionDescription);
306 }
307 fReactants->push_back(const_cast<G4Track*>(&utils.fTrackB));
308}
309
310std::unique_ptr<G4ITReactionChange>
312 G4ITReactionSet* pReactionSet, const G4double& currentStepTime,
313 const G4double& /*previousStepTime*/,
314 const G4bool& /*reachedUserStepTimeLimit*/)
315{
316 if(pReactionSet == nullptr)
317 {
318 return nullptr;
319 }
320
321 G4ITReactionPerTime& reactionPerTime = pReactionSet->GetReactionsPerTime();
322 if(reactionPerTime.empty())
323 {
324 return nullptr;
325 }
326
327 for(auto reaction_i = reactionPerTime.begin();
328 reaction_i != reactionPerTime.end(); reaction_i = reactionPerTime.begin())
329 {
330 G4Track* pTrackA = (*reaction_i)->GetReactants().first;
331 if(pTrackA->GetTrackStatus() == fStopAndKill)
332 {
333 continue;
334 }
335 G4Track* pTrackB = (*reaction_i)->GetReactant(pTrackA);
336 if(pTrackB->GetTrackStatus() == fStopAndKill)
337 {
338 continue;
339 }
340
341 if(pTrackB == pTrackA)
342 {
343 G4ExceptionDescription exceptionDescription;
344 exceptionDescription << "The IT reaction process sent back a reaction "
345 "between trackA and trackB. ";
346 exceptionDescription << "The problem is trackA == trackB";
347 G4Exception("G4DNAIndependentReactionTimeStepper::FindReaction",
348 "G4DNAIndependentReactionTimeStepper02", FatalErrorInArgument,
349 exceptionDescription);
350 }
351 pReactionSet->SelectThisReaction(*reaction_i);
352 if(fpReactionProcess != nullptr &&
353 fpReactionProcess->TestReactibility(*pTrackA, *pTrackB, currentStepTime,
354 false))
355 {
356 if((fSampledPositions.find(pTrackA->GetTrackID()) ==
357 fSampledPositions.end() &&
358 (fSampledPositions.find(pTrackB->GetTrackID()) ==
359 fSampledPositions.end())))
360 {
361 G4ExceptionDescription exceptionDescription;
362 exceptionDescription
363 << "The positions of trackA and trackB have no counted ";
364 G4Exception("G4DNAIndependentReactionTimeStepper::FindReaction",
365 "G4DNAIndependentReactionTimeStepper0001",
366 FatalErrorInArgument, exceptionDescription);
367 }
368
369 pTrackA->SetPosition(fSampledPositions[pTrackA->GetTrackID()]);
370 pTrackB->SetPosition(fSampledPositions[pTrackB->GetTrackID()]);
371 auto pReactionChange =
372 fpReactionProcess->MakeReaction(*pTrackA, *pTrackB);
373 if(pReactionChange == nullptr)
374 {
375 return nullptr;
376 }
377 return pReactionChange;
378 }
379 }
380 return nullptr;
381}
382
384 G4VDNAReactionModel* pReactionModel)
385{
386 fReactionModel = pReactionModel;
387}
388
390{
391 return fReactionModel;
392}
393
395{
396 fVerbose = flag;
397}
398
399G4double G4DNAIndependentReactionTimeStepper::GetTimeToEncounter(
400 const G4Track& trackA, const G4Track& trackB)
401{
402 G4double timeToReaction =
403 dynamic_cast<G4DiffusionControlledReactionModel*>(fReactionModel)
404 ->GetTimeToEncounter(trackA, trackB);
405 return timeToReaction;
406}
407
409 G4VITReactionProcess* pReactionProcess)
410{
411 fpReactionProcess = pReactionProcess;
412}
414 G4double /*currentGlobalTime*/, G4double definedMinTimeStep)
415{
416 G4double fTSTimeStep = DBL_MAX;
417 fCheckedTracks.clear();
418
419 for(auto pTrack : *fpTrackContainer->GetMainList())
420 {
421 if(pTrack == nullptr)
422 {
423 G4ExceptionDescription exceptionDescription;
424 exceptionDescription << "No track found.";
425 G4Exception("G4DNAIndependentReactionTimeStepper::CalculateMinTimeStep",
426 "G4DNAIndependentReactionTimeStepper006",
427 FatalErrorInArgument, exceptionDescription);
428 continue;
429 }
430
431 G4TrackStatus trackStatus = pTrack->GetTrackStatus();
432 if(trackStatus == fStopAndKill || trackStatus == fStopButAlive)
433 {
434 continue;
435 }
436
437 G4double sampledMinTimeStep = CalculateStep(*pTrack, definedMinTimeStep);
438 G4TrackVectorHandle reactants = GetReactants();
439
440 if(sampledMinTimeStep < fTSTimeStep)
441 {
442 fTSTimeStep = sampledMinTimeStep;
443 fReactionSet->CleanAllReaction();
444 if(reactants)
445 {
446 fReactionSet->AddReactions(fTSTimeStep, const_cast<G4Track*>(pTrack),
447 reactants);
448
449 fSampledPositions[pTrack->GetTrackID()] = pTrack->GetPosition();
450 for(const auto& it : *fReactants)
451 {
452 auto pTrackB = it;
453 // G4cout<<"position : "<<pTrackB->GetTrackID()<<G4endl;
454 fSampledPositions[pTrackB->GetTrackID()] = pTrackB->GetPosition();
455 }
457 }
458 }
459 else if(fTSTimeStep == sampledMinTimeStep && G4bool(reactants))
460 {
461 fReactionSet->AddReactions(fTSTimeStep, const_cast<G4Track*>(pTrack),
462 reactants);
463
464 fSampledPositions[pTrack->GetTrackID()] = pTrack->GetPosition();
465 for(const auto& it : *fReactants)
466 {
467 auto pTrackB = it;
468 // G4cout<<"position : "<<pTrackB->GetTrackID()<<G4endl;
469 fSampledPositions[pTrackB->GetTrackID()] = pTrackB->GetPosition();
470 }
472 }
473 else if(reactants)
474 {
476 }
477 }
478
479 return fTSTimeStep;
480}
#define BuildChemicalMoleculeFinder()
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4shared_ptr< std::vector< G4Track * > > G4TrackVectorHandle
Definition: G4ITReaction.hh:44
std::multiset< G4ITReactionPtr, compReactionPerTime > G4ITReactionPerTime
Definition: G4ITReaction.hh:77
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
G4TrackStatus
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double CalculateStep(const G4Track &, const G4double &) override
G4double CalculateMinTimeStep(G4double, G4double) override
std::unique_ptr< G4ITReactionChange > FindReaction(G4ITReactionSet *pReactionSet, const G4double &currentStepTime=0, const G4double &previousStepTime=0, const G4bool &reachedUserStepTimeLimit=false)
void SetReactionProcess(G4VITReactionProcess *pReactionProcess)
Data * GetReactionData(Reactant *, Reactant *) const
const ReactantList * CanReactWith(Reactant *) const
static G4double GetRCutOff()
Definition: G4IRTUtils.cc:39
void SelectThisReaction(G4ITReactionPtr reaction)
void CleanAllReaction()
void AddReactions(G4double time, G4Track *trackA, G4TrackVectorHandle reactants)
G4ITReactionPerTime & GetReactionsPerTime()
G4TrackList * GetMainList(Key)
const G4String & GetName() const
Definition: G4Molecule.cc:336
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()
Definition: G4Scheduler.cc:101
G4double GetGlobalTime() const
Definition: G4Scheduler.hh:364
G4double GetStartTime() const
Definition: G4Scheduler.hh:344
G4TrackStatus GetTrackStatus() const
void SetPosition(const G4ThreeVector &aValue)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
virtual void Initialise(const G4MolecularConfiguration *, const G4Track &)
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0
virtual std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &)=0
virtual G4bool TestReactibility(const G4Track &, const G4Track &, double, bool)=0
G4TrackVectorHandle fReactants
G4TrackVectorHandle GetReactants()
static G4ThreadLocal G4double fUserMinTimeStep
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:62