Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DiffusionControlledReactionModel.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25//
26// Author: Hoang TRAN
27
29#include "G4Track.hh"
32#include "G4Exp.hh"
33#include "G4IRTUtils.hh"
34#include "G4SystemOfUnits.hh"
35#include "G4Electron_aq.hh"
36#include "Randomize.hh"
37#include "G4Molecule.hh"
38#include "G4ErrorFunction.hh"
39
42{}
43
45 default;
46
48 const G4MolecularConfiguration* pMolecule, const G4Track&)
49{
50 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
51}
52
54 const G4MolecularConfiguration* pMolecule)
55{
56 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
57}
58
60 const G4MolecularConfiguration* pMol1, const G4MolecularConfiguration* pMol2)
61{
62 auto reactionData = fpReactionTable->GetReactionData(pMol1, pMol2);
63 if(reactionData == nullptr)
64 {
65 G4ExceptionDescription exceptionDescription;
66 exceptionDescription << "No reactionData"
67 << " for : " << pMol1->GetName() << " and "
68 << pMol2->GetName();
69 G4Exception("G4DiffusionControlledReactionModel"
70 "::GetReactionRadius()",
71 "G4DiffusionControlledReactionModel00", FatalException,
72 exceptionDescription);
73 return 0.;
74 }else
75 {
76 return reactionData->GetEffectiveReactionRadius();
77 }
78}
79
81{
82 auto pMol1 = (*fpReactionData)[i]->GetReactant1();
83 auto pMol2 = (*fpReactionData)[i]->GetReactant2();
84 return GetReactionRadius(pMol1, pMol2);
85}
86
88 const G4Track& trackA, const G4Track& trackB)
89{
90 auto pMolConfA = GetMolecule(trackA)->GetMolecularConfiguration();
91 auto pMolConfB = GetMolecule(trackB)->GetMolecularConfiguration();
92
93 G4double D =
94 pMolConfA->GetDiffusionCoefficient() + pMolConfB->GetDiffusionCoefficient();
95
96 if(D == 0)
97 {
98 G4ExceptionDescription exceptionDescription;
99 exceptionDescription << "The total diffusion coefficient for : "
100 << pMolConfA->GetName() << " and "
101 << pMolConfB->GetName() << " is null ";
102 G4Exception("G4DiffusionControlledReactionModel"
103 "::GetTimeToEncounter()",
104 "G4DiffusionControlledReactionModel03", FatalException,
105 exceptionDescription);
106 }
107
109 pMolConfA, pMolConfB);
110 G4double kobs = reactionData->GetObservedReactionRateConstant();
111 G4double distance = (trackA.GetPosition() - trackB.GetPosition()).mag();
112 G4double SmoluchowskiRadius = reactionData->GetEffectiveReactionRadius();
113
114 if(distance == 0 || distance < SmoluchowskiRadius)
115 {
116 G4ExceptionDescription exceptionDescription;
117 exceptionDescription << "distance = " << distance << " is uncorrected with "
118 << " Reff = " << SmoluchowskiRadius
119 << " for : " << pMolConfA->GetName() << " and "
120 << pMolConfB->GetName();
121 G4Exception("G4DiffusionControlledReactionModel"
122 "::GetTimeToEncounter()",
123 "G4DiffusionControlledReactionModel02", FatalException,
124 exceptionDescription);
125 }
126 else
127 {
128 G4double Winf = SmoluchowskiRadius / distance;
130 G4double X = 0;
131 G4double irt_1 = -1.0 * ps;
132 if(Winf > 0 && U < Winf)
133 {
134 G4double erfcIn = G4ErrorFunction::erfcInv(U / Winf);
135 if(erfcIn != 0)
136 {
137 G4double d =
138 (distance - SmoluchowskiRadius) / erfcIn;
139 irt_1 = (1.0 / (4 * D)) * d * d;
140 }
141 }
142
143 if(reactionData->GetReactionType() == 0) // Totally diffused contr
144 {
145 return irt_1;
146 }
147
148 if(irt_1 < 0)
149 {
150 return irt_1;
151 }
152 else
153 {
154 G4double kdif = 4 * CLHEP::pi * D * SmoluchowskiRadius * Avogadro;
155
156 if(pMolConfA == pMolConfB)
157 {
158 kdif /= 2;
159 }
160 G4double kact = G4IRTUtils::GetKact(kobs, kdif);
161 G4double sumOfk = kact + kdif;
162 if(sumOfk != 0)
163 {
164 G4double rateFactor = kact / sumOfk;
165 if(G4UniformRand() > rateFactor)
166 {
167 return -1.0 * ps;
168 }
169 G4double Y = std::abs(G4RandGauss::shoot(0.0, std::sqrt(2)));
170
171 if(Y > 0)
172 {
173 X = -(G4Log(G4UniformRand())) / Y;
174 }
175 G4double f = X * SmoluchowskiRadius * kdif / sumOfk;
176 G4double irt_2 = (f * f) / D;
177 return irt_1 + irt_2;
178 }
179 }
180 }
181 return -1.0 * ps;
182}
G4double D(G4double temp)
G4double Y(G4double density)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Data * GetReactionData(Reactant *, Reactant *) const
static G4DNAMolecularReactionTable * Instance()
G4double GetTimeToEncounter(const G4Track &trackA, const G4Track &trackB)
void InitialiseToPrint(const G4MolecularConfiguration *) override
G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *) override
void Initialise(const G4MolecularConfiguration *, const G4Track &) override
static G4double erfcInv(G4double x)
static G4double GetKact(const G4double &obs, const G4double &dif)
Definition: G4IRTUtils.hh:41
const G4String & GetName() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:530
const G4ThreeVector & GetPosition() const
const G4DNAMolecularReactionTable * fpReactionTable