Geant4-11
G4DNAPartiallyDiffusionControlled.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
28#include "G4IRTUtils.hh"
31#include "G4SystemOfUnits.hh"
34#include "Randomize.hh"
35#include "G4Molecule.hh"
36#include "G4ITReactionChange.hh"
37#include "G4VReactionType.hh"
38#include "G4Electron_aq.hh"
39#include "G4ErrorFunction.hh"
42{}
43
45
49{
50 auto reactionData = G4DNAMolecularReactionTable::Instance()
51 ->GetReactionData(mA, mB);
52
55
56 const G4double Rs = 0.3 * nm;
57 G4double kobs = reactionData->GetObservedReactionRateConstant() / Avogadro;
58 if(mA->GetCharge() * mB->GetCharge() == 0)
59 {
60 G4double kdif = 4 * CLHEP::pi * D * R * Avogadro;
61 G4double kact = G4IRTUtils::GetKact(kobs, kdif);
62 return G4UniformRand() < Rs / ( Rs + ( kdif / kact ) * ( R + Rs ));
63 }
64 else
65 {
66 G4double rc = 0.71 * nm * mA->GetCharge() *
67 mB->GetCharge();
69 G4double kdif = 4 * CLHEP::pi * D * sigmaEff;
70 G4double kact = G4IRTUtils::GetKact(kobs, kdif);
71 G4double a = std::exp( -rc / R );
72 G4double b = std::exp( -rc / ( R + Rs ) );
73 G4double Preact = ( a - b ) / ( a - b - ( kdif / kact ) * ( 1 - a ) );
74
75 return G4UniformRand() < Preact;
76 }
77}
78
82{
83 G4double D;
84
85 if(mA == mB)
86 {
87 D = (mA->GetDiffusionCoefficient());
88 }
89 else
90 {
91 D = (mA->GetDiffusionCoefficient() +
93 }
94 return D;
95}
96
98 const G4Track& trackB)
99{
100 auto pMolConfA = GetMolecule(trackA)->GetMolecularConfiguration();
101 auto pMolConfB = GetMolecule(trackB)->GetMolecularConfiguration();
102
103 G4double D = GetDiffusionCoefficient(pMolConfA, pMolConfB);
104 auto reactionData = G4DNAMolecularReactionTable::Instance()
105 ->GetReactionData(pMolConfA, pMolConfB);
106 G4double Reff;
107 G4double kobs = reactionData->GetObservedReactionRateConstant();
108 G4double distance = (trackA.GetPosition() - trackB.GetPosition()).mag();
109 G4double SmoluchowskiRadius;
110 G4double RVal = pMolConfA->GetVanDerVaalsRadius() + pMolConfB->GetVanDerVaalsRadius();
111
112 if((pMolConfA->GetCharge() != 0) &&
113 (pMolConfB->GetCharge() != 0))
114 {
115 G4double rc = 0.71 * nm * pMolConfA->GetCharge() *
116 pMolConfB->GetCharge();
117 distance = G4IRTUtils::EffectiveDistance( rc, distance );
118 Reff = G4IRTUtils::EffectiveDistance( rc, RVal );
119 SmoluchowskiRadius = Reff;
120 }
121 else
122 {
123 SmoluchowskiRadius = RVal;
124 }
125
126 G4double Winf = SmoluchowskiRadius / distance;
127 G4double U1 = G4UniformRand();
128 G4double U2 = G4UniformRand();
130 G4double X = 0;
131 G4double irt_1 = -1.0 * ps;
132 G4double irt_2;
133
134 G4double kdif = 4 * CLHEP::pi * D * SmoluchowskiRadius * Avogadro;
135 G4double kact = G4IRTUtils::GetKact(kobs, kdif);
136
137 if ( U < Winf )
138 {
139 G4double d = ( distance - SmoluchowskiRadius ) /
140 G4ErrorFunction::erfcInv( U / Winf );
141 irt_1 = ( 1.0 / ( 4 * D ) ) * d * d;
142 }
143
144 if( irt_1 < 0)
145 {
146 return irt_1;
147 }
148 else
149 {
150 G4double rateFactor = kact / ( kact + kdif );
151 if( U1 > rateFactor )
152 {
153 return -1.0 * ps;
154 }
155 G4double Y = std::abs(G4RandGauss::shoot(0.0,std::sqrt(2)));
156
157 if( Y > 0)
158 {
159 X = - ( G4Log( U2 ) ) / Y;
160 }
161
162 G4double f = X * SmoluchowskiRadius * kdif / ( kact + kdif );
163 irt_2 = ( f * f ) / D ;
164 }
165
166 return irt_1 + irt_2;
167}
168
G4double D(G4double temp)
G4double Y(G4double density)
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
static constexpr double nm
Definition: G4SIunits.hh:92
static constexpr double ps
Definition: G4SIunits.hh:157
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
Data * GetReactionData(Reactant *, Reactant *) const
static G4DNAMolecularReactionTable * Instance()
G4double GetTimeToEncounter(const G4Track &trackA, const G4Track &trackB) override
G4double GetDiffusionCoefficient(const G4MolecularConfiguration *, const G4MolecularConfiguration *)
G4bool GeminateRecombinationProbability(const G4MolecularConfiguration *, const G4MolecularConfiguration *) override
static G4double erfcInv(G4double x)
static G4double EffectiveDistance(const G4double &rc, const G4double &r0)
Definition: G4IRTUtils.cc:32
static G4double GetKact(const G4double &obs, const G4double &dif)
Definition: G4IRTUtils.hh:41
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
const G4ThreeVector & GetPosition() const
static constexpr double pi
Definition: SystemOfUnits.h:55
ThreeVector shoot(const G4int Ap, const G4int Af)
float Avogadro
Definition: hepunit.py:252