Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAIRT.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//
26/*
27 * G4DNAIRT.cc
28 *
29 * Created on: Jul 23, 2019
30 * Author: W. G. Shin
31 * J. Ramos-Mendez and B. Faddegon
32*/
33
34
35#include "G4DNAIRT.hh"
36#include "G4ErrorFunction.hh"
37#include "G4SystemOfUnits.hh"
39#include "Randomize.hh"
42#include "G4Molecule.hh"
43#include "G4ITReactionChange.hh"
44#include "G4ITTrackHolder.hh"
45#include "G4ITReaction.hh"
46#include "G4Scheduler.hh"
47
48using namespace std;
49
52fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable)),
53fpReactionModel(nullptr),
54fTrackHolder(G4ITTrackHolder::Instance()),
55fReactionSet(0)
56{
58 timeMax = G4Scheduler::Instance()->GetEndTime();
59
60 fXMin = 1e9*nm;
61 fYMin = 1e9*nm;
62 fZMin = 1e9*nm;
63
64 fXMax = 0e0*nm;
65 fYMax = 0e0*nm;
66 fZMax = 0e0*nm;
67
68 fNx = 0;
69 fNy = 0;
70 fNz = 0;
71
72 xiniIndex = 0, yiniIndex = 0, ziniIndex = 0;
73 xendIndex = 0, yendIndex = 0, zendIndex = 0;
74
75 fRCutOff =
76 1.45 * nm + 2 * std::sqrt(8*9.46e9*nm*nm/s * timeMax); // 95% confidence level
77
78 erfc = new G4ErrorFunction();
79}
80
81
83 : G4DNAIRT()
84{
85 fpReactionModel = pReactionModel;
86}
87
89{
90 delete erfc;
91}
92
94
95 fTrackHolder = G4ITTrackHolder::Instance();
96
97 fReactionSet = G4ITReactionSet::Instance();
98 fReactionSet->CleanAllReaction();
99 fReactionSet->SortByTime();
100
101 spaceBinned.clear();
102
103 timeMin = G4Scheduler::Instance()->GetStartTime();
104 timeMax = G4Scheduler::Instance()->GetEndTime();
105
106 xiniIndex = 0;
107 yiniIndex = 0;
108 ziniIndex = 0;
109 xendIndex = 0;
110 yendIndex = 0;
111 zendIndex = 0;
112
113 fXMin = 1e9*nm;
114 fYMin = 1e9*nm;
115 fZMin = 1e9*nm;
116
117 fXMax = 0e0*nm;
118 fYMax = 0e0*nm;
119 fZMax = 0e0*nm;
120
121 fNx = 0;
122 fNy = 0;
123 fNz = 0;
124
125 SpaceBinning(); // 1. binning the space
126 IRTSampling(); // 2. Sampling of the IRT
127
128}
129
131 auto it_begin = fTrackHolder->GetMainList()->begin();
132 while(it_begin != fTrackHolder->GetMainList()->end()){
133
134 G4ThreeVector position = it_begin->GetPosition();
135
136 if ( fXMin > position.x() ) fXMin = position.x();
137 if ( fYMin > position.y() ) fYMin = position.y();
138 if ( fZMin > position.z() ) fZMin = position.z();
139
140 if ( fXMax < position.x() ) fXMax = position.x();
141 if ( fYMax < position.y() ) fYMax = position.y();
142 if ( fZMax < position.z() ) fZMax = position.z();
143
144 ++it_begin;
145 }
146
147 fNx = G4int((fXMax-fXMin)/fRCutOff) == 0 ? 1 : G4int((fXMax-fXMin)/fRCutOff);
148 fNy = G4int((fYMax-fYMin)/fRCutOff) == 0 ? 1 : G4int((fYMax-fYMin)/fRCutOff);
149 fNz = G4int((fZMax-fZMin)/fRCutOff) == 0 ? 1 : G4int((fZMax-fZMin)/fRCutOff);
150
151}
152
154
155 auto it_begin = fTrackHolder->GetMainList()->begin();
156 while(it_begin != fTrackHolder->GetMainList()->end()){
157 G4int I = FindBin(fNx, fXMin, fXMax, it_begin->GetPosition().x());
158 G4int J = FindBin(fNy, fYMin, fYMax, it_begin->GetPosition().y());
159 G4int K = FindBin(fNz, fZMin, fZMax, it_begin->GetPosition().z());
160
161 spaceBinned[I][J][K].push_back(*it_begin);
162
163 Sampling(*it_begin);
164 ++it_begin;
165 }
166}
167
169 G4Molecule* molA = G4Molecule::GetMolecule(track);
170 const G4MolecularConfiguration* molConfA = molA->GetMolecularConfiguration();
171 if(molConfA->GetDiffusionCoefficient() == 0) return;
172
173 const vector<const G4MolecularConfiguration*>* reactivesVector =
175
176 if(reactivesVector == nullptr) return;
177
179 G4double minTime = timeMax;
180
181 xiniIndex = FindBin(fNx, fXMin, fXMax, track->GetPosition().x()-fRCutOff);
182 xendIndex = FindBin(fNx, fXMin, fXMax, track->GetPosition().x()+fRCutOff);
183 yiniIndex = FindBin(fNy, fYMin, fYMax, track->GetPosition().y()-fRCutOff);
184 yendIndex = FindBin(fNy, fYMin, fYMax, track->GetPosition().y()+fRCutOff);
185 ziniIndex = FindBin(fNz, fZMin, fZMax, track->GetPosition().z()-fRCutOff);
186 zendIndex = FindBin(fNz, fZMin, fZMax, track->GetPosition().z()+fRCutOff);
187
188 for ( int ii = xiniIndex; ii <= xendIndex; ii++ ) {
189 for ( int jj = yiniIndex; jj <= yendIndex; jj++ ) {
190 for ( int kk = ziniIndex; kk <= zendIndex; kk++ ) {
191
192 std::vector<G4Track*> spaceBin = spaceBinned[ii][jj][kk];
193 for ( int n = 0; n < (int)spaceBinned[ii][jj][kk].size(); n++ ) {
194 if(!spaceBin[n] || track == spaceBin[n]) continue;
195 if(spaceBin[n]->GetTrackStatus() == fStopButAlive) continue;
196
197 G4Molecule* molB = G4Molecule::GetMolecule(spaceBin[n]);
198 if(!molB) continue;
199
200 const G4MolecularConfiguration* molConfB = molB->GetMolecularConfiguration();
201 if(molConfB->GetDiffusionCoefficient() == 0) continue;
202
203 auto it = std::find(reactivesVector->begin(), reactivesVector->end(), molConfB);
204 if(it == reactivesVector->end()) continue;
205
206 G4ThreeVector orgPosB = spaceBin[n]->GetPosition();
207 G4double dt = track->GetGlobalTime() - spaceBin[n]->GetGlobalTime();
208 G4ThreeVector newPosB = orgPosB;
209
210 if(dt > 0){
211 G4double sigma, x, y, z;
212 G4double diffusionCoefficient = G4Molecule::GetMolecule(spaceBin[n])->GetDiffusionCoefficient();
213
214 sigma = std::sqrt(2.0 * diffusionCoefficient * dt);
215
216 x = G4RandGauss::shoot(0., 1.0)*sigma;
217 y = G4RandGauss::shoot(0., 1.0)*sigma;
218 z = G4RandGauss::shoot(0., 1.0)*sigma;
219
220 newPosB = orgPosB + G4ThreeVector(x,y,z);
221 }else if(dt < 0) continue;
222
223 G4double r0 = (newPosB - track->GetPosition()).mag();
225 molConfB,
226 r0);
227 if(irt>=0 && irt<timeMax - globalTime)
228 {
229 irt += globalTime;
230 if(irt < minTime) minTime = irt;
231#ifdef DEBUG
232 G4cout<<irt<<'\t'<<molConfA->GetName()<<" "<<track->GetTrackID()<<'\t'<<molConfB->GetName()<<" "<<spaceBin[n]->GetTrackID()<<'\n';
233#endif
234 fReactionSet->AddReaction(irt,track,spaceBin[n]);
235 }
236 }
237 spaceBin.clear();
238 }
239 }
240 }
241
242// Scavenging & first order reactions
243
244 auto fReactionDatas = fMolReactionTable->GetReactionData(molConfA);
245 G4double index = -1;
246
247 for(size_t u=0; u<fReactionDatas->size();u++){
248 if((*fReactionDatas)[u]->GetReactant2()->GetDiffusionCoefficient() == 0){
249 G4double kObs = (*fReactionDatas)[u]->GetObservedReactionRateConstant();
250 G4double time = -(std::log(1.0 - G4UniformRand())/kObs) + globalTime;
251 if( time < minTime && time >= globalTime && time < timeMax){
252 minTime = time;
253 index = (int) u;
254 }
255 }
256 }
257
258 if(index != -1){
259#ifdef DEBUG
260 G4cout<<"scavenged: "<<minTime<<'\t'<<molConfA->GetName()<<it_begin->GetTrackID()<<'\n';
261#endif
262 G4Molecule* fakeMol = new G4Molecule((*fReactionDatas)[index]->GetReactant2());
263 G4Track* fakeTrack = fakeMol->BuildTrack(globalTime,track->GetPosition());
264 fTrackHolder->Push(fakeTrack);
265 fReactionSet->AddReaction(minTime, track, fakeTrack);
266 }
267}
268
269
271 const auto pMoleculeA = molA;
272 const auto pMoleculeB = molB;
273 auto fReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
274 G4int reactionType = fReactionData->GetReactionType();
275 G4double r0 = distance;
276 if(r0 == 0) r0 += 1e-3*nm;
277 G4double irt = -1 * ps;
280 G4double rc = fReactionData->GetOnsagerRadius();
281
282 if ( reactionType == 0){
283 G4double sigma = fReactionData->GetEffectiveReactionRadius();
284
285 if( rc != 0) r0 = -rc / (1-std::exp(rc/r0));
286 if(sigma > r0) return 0; // contact reaction
287
288 G4double Winf = sigma/r0;
290
291 if ( W < Winf ) irt = (0.25/D) * std::pow( (r0-sigma)/erfc->erfcInv(r0*W/sigma), 2 );
292
293 return irt;
294 }
295 else if ( reactionType == 1 ){
296 G4double sigma = fReactionData->GetReactionRadius();
297 G4double kact = fReactionData->GetActivationRateConstant();
298 G4double kdif = fReactionData->GetDiffusionRateConstant();
299 G4double kobs = fReactionData->GetObservedReactionRateConstant();
300
301 G4double a, b, Winf;
302
303 if ( rc == 0 ) {
304 a = 1/sigma * kact / kobs;
305 b = (r0 - sigma) / 2;
306 } else {
307 G4double v = kact/Avogadro/(4*CLHEP::pi*pow(sigma,2) * exp(-rc / sigma));
308 G4double alpha = v+rc*D/(pow(sigma,2)*(1-exp(-rc/sigma)));
309 a = 4*pow(sigma,2)*alpha/(D*pow(rc,2))*pow(sinh(rc/(2*sigma)),2);
310 b = rc/4*(cosh(rc/(2*r0))/sinh(rc/(2*r0))-cosh(rc/(2*sigma))/sinh(rc/(2*sigma)));
311 r0 = -rc/(1-std::exp(rc/r0));
312 sigma = fReactionData->GetEffectiveReactionRadius();
313 }
314
315 if(sigma > r0){
316 if(fReactionData->GetProbability() > G4UniformRand()) return 0;
317 else return irt;
318 }
319 Winf = sigma / r0 * kobs / kdif;
320
321 if(Winf > G4UniformRand()) irt = SamplePDC(a,b)/D;
322 return irt;
323 }
324
325 return -1 * ps;
326}
327
329
330 G4int bin = -1;
331 if ( value <= xmin )
332 bin = 0; //1;
333 else if ( value >= xmax) //!(xmax < value) ) //value >= xmax )
334 bin = n-1; //n;
335 else
336 bin = G4int( n * ( value - xmin )/( xmax - xmin ) ); //bin = 1 + G4int( n * ( value - xmin )/( xmax - xmin ) );
337
338 if ( bin < 0 ) bin = 0;
339 if ( bin >= n ) bin = n-1;
340
341 return bin;
342}
343
345
346 G4double p = 2.0 * std::sqrt(2.0*b/a);
347 G4double q = 2.0 / std::sqrt(2.0*b/a);
348 G4double M = max(1.0/(a*a),3.0*b/a);
349
350 G4double X, U, lambdax;
351
352 G4int ntrials = 0;
353 while(1) {
354
355 // Generate X
356 U = G4UniformRand();
357 if ( U < p/(p + q * M) ) X = pow(U * (p + q * M) / 2, 2);
358 else X = pow(2/((1-U)*(p+q*M)/M),2);
359
360 U = G4UniformRand();
361
362 lambdax = std::exp(-b*b/X) * ( 1.0 - a * std::sqrt(CLHEP::pi * X) * erfc->erfcx(b/std::sqrt(X) + a*std::sqrt(X)));
363
364 if ((X <= 2.0*b/a && U <= lambdax) ||
365 (X >= 2.0*b/a && U*M/X <= lambdax)) break;
366
367 ntrials++;
368
369 if ( ntrials > 10000 ){
370 G4cout<<"Totally rejected"<<'\n';
371 return -1.0;
372 }
373 }
374 return X;
375}
376
377std::unique_ptr<G4ITReactionChange> G4DNAIRT::MakeReaction(const G4Track& trackA,
378 const G4Track& trackB)
379{
380
381 std::unique_ptr<G4ITReactionChange> pChanges(new G4ITReactionChange());
382 pChanges->Initialize(trackA, trackB);
383
384 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
385 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
386 const auto pReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
387
389 G4double effectiveReactionRadius = pReactionData->GetEffectiveReactionRadius();
390
391 const G4double D1 = pMoleculeA->GetDiffusionCoefficient();
392 const G4double D2 = pMoleculeB->GetDiffusionCoefficient();
393
394 G4ThreeVector r1 = trackA.GetPosition();
395 G4ThreeVector r2 = trackB.GetPosition();
396
397 if(r1 == r2) r2 += G4ThreeVector(0,0,1e-3*nm);
398
399 G4ThreeVector S1 = r1 - r2;
400
401 G4double r0 = S1.mag();
402
403 S1.setMag(effectiveReactionRadius);
404
405 G4double dt = globalTime - trackA.GetGlobalTime();
406
407 if(dt != 0){
408 G4double s12 = 2.0 * D1 * dt;
409 G4double s22 = 2.0 * D2 * dt;
410 if(s12 == 0) r2 = r1;
411 else if(s22 == 0) r1 = r2;
412 else{
413 G4double alpha = effectiveReactionRadius * r0 / (2*(D1 + D2)*dt);
414 G4ThreeVector S2 = (r1 + (s12 / s22)*r2) + G4ThreeVector(G4RandGauss::shoot(0, s12 + s22 * s22 / s12),
415 G4RandGauss::shoot(0, s12 + s22 * s22 / s12),
416 G4RandGauss::shoot(0, s12 + s22 * s22 / s12));
417
418 S1.setPhi(rad * G4UniformRand() * 2.0 * CLHEP::pi);
419 S1.setTheta(rad * std::acos(1.0 + 1./alpha * std::log(1.0 - G4UniformRand() * (1 - std::exp(-2.0 * alpha)))));
420
421 r1 = (D1 * S1 + D2 * S2) / (D1 + D2);
422 r2 = D2 * (S2 - S1) / (D1 + D2);
423 }
424 }
425
426 auto pTrackA = const_cast<G4Track*>(pChanges->GetTrackA());
427 auto pTrackB = const_cast<G4Track*>(pChanges->GetTrackB());
428
429 pTrackA->SetPosition(r1);
430 pTrackB->SetPosition(r2);
431
432 pTrackA->SetGlobalTime(globalTime);
433 pTrackB->SetGlobalTime(globalTime);
434
435 pTrackA->SetTrackStatus(fStopButAlive);
436 pTrackB->SetTrackStatus(fStopButAlive);
437
438 const G4int nbProducts = pReactionData->GetNbProducts();
439
440 if(nbProducts){
441
442 const G4double sqrD1 = D1 == 0. ? 0. : std::sqrt(D1);
443 const G4double sqrD2 = D2 == 0. ? 0. : std::sqrt(D2);
444 const G4double inv_numerator = 1./(sqrD1 + sqrD2);
445 const G4ThreeVector reactionSite = sqrD2 * inv_numerator * trackA.GetPosition()
446 + sqrD1 * inv_numerator * trackB.GetPosition();
447
448 std::vector<G4ThreeVector> position;
449
450 if(nbProducts == 1){
451 position.push_back(reactionSite);
452 }else if(nbProducts == 2){
453 position.push_back(trackA.GetPosition());
454 position.push_back(trackB.GetPosition());
455 }else if (nbProducts == 3){
456 position.push_back(reactionSite);
457 position.push_back(trackA.GetPosition());
458 position.push_back(trackB.GetPosition());
459 }
460
461 for(G4int u = 0; u < nbProducts; u++){
462
463 auto product = new G4Molecule(pReactionData->GetProduct(u));
464 auto productTrack = product->BuildTrack(globalTime,
465 position[u]);
466
467 productTrack->SetTrackStatus(fAlive);
468 fTrackHolder->Push(productTrack);
469
470 pChanges->AddSecondary(productTrack);
471
472 G4int I = FindBin(fNx, fXMin, fXMax, position[u].x());
473 G4int J = FindBin(fNy, fYMin, fYMax, position[u].y());
474 G4int K = FindBin(fNz, fZMin, fZMax, position[u].z());
475
476 spaceBinned[I][J][K].push_back(productTrack);
477
478 Sampling(productTrack);
479 }
480 }
481
482 fTrackHolder->MergeSecondariesWithMainList();
483 pChanges->KillParents(true);
484 return pChanges;
485}
486
487
488std::vector<std::unique_ptr<G4ITReactionChange>> G4DNAIRT::FindReaction(
489 G4ITReactionSet* pReactionSet,
490 const double /*currentStepTime*/,
491 const double fGlobalTime,
492 const bool /*reachedUserStepTimeLimit*/)
493{
494 std::vector<std::unique_ptr<G4ITReactionChange>> fReactionInfo;
495 fReactionInfo.clear();
496
497 if (pReactionSet == nullptr)
498 {
499 return fReactionInfo;
500 }
501
502 auto fReactionsetInTime = pReactionSet->GetReactionsPerTime();
503 assert(fReactionsetInTime.begin() != fReactionsetInTime.end());
504
505 auto it_begin = fReactionsetInTime.begin();
506 while(it_begin != fReactionsetInTime.end())
507 {
508 G4double irt = it_begin->get()->GetTime();
509
510 if(fGlobalTime < irt) break;
511
512 pReactionSet->SelectThisReaction(*it_begin);
513
514 G4Track* pTrackA = it_begin->get()->GetReactants().first;
515 G4Track* pTrackB = it_begin->get()->GetReactants().second;
516 auto pReactionChange = MakeReaction(*pTrackA, *pTrackB);
517
518 if(pReactionChange){
519 fReactionInfo.push_back(std::move(pReactionChange));
520 }
521
522 fReactionsetInTime = pReactionSet->GetReactionsPerTime();
523 it_begin = fReactionsetInTime.begin();
524 }
525
526 return fReactionInfo;
527}
528
530 const G4Track& /*trackB*/,
531 double /*currentStepTime*/,
532 bool /*userStepTimeLimit*/) /*const*/
533{
534 return true;
535}
536
538{
539 fpReactionModel = model;
540}
double D(double temp)
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
ReturnType & reference_cast(OriginalType &source)
CLHEP::Hep3Vector G4ThreeVector
@ fAlive
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
void setTheta(double)
double mag() const
void setMag(double)
Definition: ThreeVector.cc:20
void setPhi(double)
void IRTSampling()
Definition: G4DNAIRT.cc:153
void SetReactionModel(G4VDNAReactionModel *)
Definition: G4DNAIRT.cc:537
G4double SamplePDC(G4double, G4double)
Definition: G4DNAIRT.cc:344
const G4DNAMolecularReactionTable *& fMolReactionTable
Definition: G4DNAIRT.hh:93
void Initialize() override
Definition: G4DNAIRT.cc:93
G4int FindBin(G4int, G4double, G4double, G4double)
Definition: G4DNAIRT.cc:328
G4double GetIndependentReactionTime(const G4MolecularConfiguration *, const G4MolecularConfiguration *, G4double)
Definition: G4DNAIRT.cc:270
G4DNAIRT()
Definition: G4DNAIRT.cc:50
void SpaceBinning()
Definition: G4DNAIRT.cc:130
std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &) override
Definition: G4DNAIRT.cc:377
std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction(G4ITReactionSet *, const double, const double, const bool) override
Definition: G4DNAIRT.cc:488
~G4DNAIRT() override
Definition: G4DNAIRT.cc:88
G4bool TestReactibility(const G4Track &, const G4Track &, double, bool) override
Definition: G4DNAIRT.cc:529
G4VDNAReactionModel * fpReactionModel
Definition: G4DNAIRT.hh:94
void Sampling(G4Track *)
Definition: G4DNAIRT.cc:168
Data * GetReactionData(Reactant *, Reactant *) const
const ReactantList * CanReactWith(Reactant *) const
static G4double erfcx(G4double x)
static G4double erfcInv(G4double x)
iterator begin()
iterator end()
static G4ITReactionSet * Instance()
void SelectThisReaction(G4ITReactionPtr reaction)
void CleanAllReaction()
void AddReaction(double time, G4Track *trackA, G4Track *trackB)
G4ITReactionPerTime & GetReactionsPerTime()
G4TrackList * GetMainList(Key)
virtual void Push(G4Track *)
void MergeSecondariesWithMainList()
static G4ITTrackHolder * Instance()
const G4String & GetName() const
static G4Molecule * GetMolecule(const G4Track *)
Definition: G4Molecule.cc:90
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
G4Track * BuildTrack(G4double globalTime, const G4ThreeVector &Position)
Definition: G4Molecule.cc:373
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:516
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4double GetEndTime() const
Definition: G4Scheduler.hh:335
G4double GetGlobalTime() const
Definition: G4Scheduler.hh:350
G4double GetStartTime() const
Definition: G4Scheduler.hh:330
void SetPosition(const G4ThreeVector &aValue)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const