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

#include <G4DiffusionControlledReactionModel.hh>

+ Inheritance diagram for G4DiffusionControlledReactionModel:

Public Member Functions

 G4DiffusionControlledReactionModel ()
 
 ~G4DiffusionControlledReactionModel () override
 
 G4DiffusionControlledReactionModel (const G4DiffusionControlledReactionModel &)=delete
 
G4DiffusionControlledReactionModeloperator= (const G4DiffusionControlledReactionModel &)=delete
 
void Initialise (const G4MolecularConfiguration *, const G4Track &) override
 
void InitialiseToPrint (const G4MolecularConfiguration *) override
 
G4double GetReactionRadius (const G4MolecularConfiguration *, const G4MolecularConfiguration *) override
 
G4double GetReactionRadius (const G4int &) override
 
G4bool FindReaction (const G4Track &, const G4Track &, G4double, G4double &, G4bool) override
 
G4double GetTimeToEncounter (const G4Track &trackA, const G4Track &trackB)
 
- Public Member Functions inherited from G4VDNAReactionModel
 G4VDNAReactionModel ()
 
 G4VDNAReactionModel (const G4VDNAReactionModel &)=delete
 
G4VDNAReactionModeloperator= (const G4VDNAReactionModel &)=delete
 
virtual ~G4VDNAReactionModel ()
 
void SetReactionTable (const G4DNAMolecularReactionTable *)
 
const G4DNAMolecularReactionTableGetReactionTable ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VDNAReactionModel
const G4DNAMolecularReactionTablefpReactionTable {nullptr}
 

Detailed Description

Definition at line 35 of file G4DiffusionControlledReactionModel.hh.

Constructor & Destructor Documentation

◆ G4DiffusionControlledReactionModel() [1/2]

G4DiffusionControlledReactionModel::G4DiffusionControlledReactionModel ( )
default

◆ ~G4DiffusionControlledReactionModel()

G4DiffusionControlledReactionModel::~G4DiffusionControlledReactionModel ( )
overridedefault

◆ G4DiffusionControlledReactionModel() [2/2]

G4DiffusionControlledReactionModel::G4DiffusionControlledReactionModel ( const G4DiffusionControlledReactionModel & )
delete

Member Function Documentation

◆ FindReaction()

G4bool G4DiffusionControlledReactionModel::FindReaction ( const G4Track & ,
const G4Track & ,
G4double ,
G4double & ,
G4bool  )
inlineoverridevirtual

Implements G4VDNAReactionModel.

Definition at line 52 of file G4DiffusionControlledReactionModel.hh.

56 {
57 return true;
58 }

◆ GetReactionRadius() [1/2]

G4double G4DiffusionControlledReactionModel::GetReactionRadius ( const G4int & i)
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 75 of file G4DiffusionControlledReactionModel.cc.

76{
77 auto pMol1 = (*fpReactionData)[i]->GetReactant1();
78 auto pMol2 = (*fpReactionData)[i]->GetReactant2();
79 return GetReactionRadius(pMol1, pMol2);
80}
G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *) override

◆ GetReactionRadius() [2/2]

G4double G4DiffusionControlledReactionModel::GetReactionRadius ( const G4MolecularConfiguration * pMol1,
const G4MolecularConfiguration * pMol2 )
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 56 of file G4DiffusionControlledReactionModel.cc.

58{
59 auto reactionData = fpReactionTable->GetReactionData(pMol1, pMol2);
60 if(reactionData == nullptr)
61 {
62 G4ExceptionDescription exceptionDescription;
63 exceptionDescription << "No reactionData"
64 << " for : " << pMol1->GetName() << " and "
65 << pMol2->GetName();
66 G4Exception("G4DiffusionControlledReactionModel"
67 "::GetReactionRadius()",
68 "G4DiffusionControlledReactionModel00", FatalException,
69 exceptionDescription);
70 return 0.;
71 }
72 return reactionData->GetEffectiveReactionRadius();
73}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
Data * GetReactionData(Reactant *, Reactant *) const
const G4String & GetName() const
const G4DNAMolecularReactionTable * fpReactionTable

Referenced by GetReactionRadius().

◆ GetTimeToEncounter()

G4double G4DiffusionControlledReactionModel::GetTimeToEncounter ( const G4Track & trackA,
const G4Track & trackB )

Definition at line 82 of file G4DiffusionControlledReactionModel.cc.

84{
85 auto pMolConfA = GetMolecule(trackA)->GetMolecularConfiguration();
86 auto pMolConfB = GetMolecule(trackB)->GetMolecularConfiguration();
87
88 G4double D =
89 pMolConfA->GetDiffusionCoefficient() + pMolConfB->GetDiffusionCoefficient();
90
91 if(D == 0)
92 {
93 G4ExceptionDescription exceptionDescription;
94 exceptionDescription << "The total diffusion coefficient for : "
95 << pMolConfA->GetName() << " and "
96 << pMolConfB->GetName() << " is null ";
97 G4Exception("G4DiffusionControlledReactionModel"
98 "::GetTimeToEncounter()",
99 "G4DiffusionControlledReactionModel03", FatalException,
100 exceptionDescription);
101 }
102
104 pMolConfA, pMolConfB);
105 G4double kobs = reactionData->GetObservedReactionRateConstant();
106 G4double distance = (trackA.GetPosition() - trackB.GetPosition()).mag();
107 G4double SmoluchowskiRadius = reactionData->GetEffectiveReactionRadius();
108
109 if(distance == 0 || distance < SmoluchowskiRadius)
110 {
111 G4ExceptionDescription exceptionDescription;
112 exceptionDescription << "distance = " << distance << " is uncorrected with "
113 << " Reff = " << SmoluchowskiRadius
114 << " for : " << pMolConfA->GetName() << " and "
115 << pMolConfB->GetName();
116 G4Exception("G4DiffusionControlledReactionModel"
117 "::GetTimeToEncounter()",
118 "G4DiffusionControlledReactionModel02", FatalException,
119 exceptionDescription);
120 }
121 else
122 {
123 G4double Winf = SmoluchowskiRadius / distance;
125 G4double X = 0;
126 G4double irt_1 = -1.0 * ps;
127 if(Winf > 0 && U < Winf)
128 {
129 G4double erfcIn = G4ErrorFunction::erfcInv(U / Winf);
130 if(erfcIn != 0)
131 {
132 G4double d =
133 (distance - SmoluchowskiRadius) / erfcIn;
134 irt_1 = (1.0 / (4 * D)) * d * d;
135 }
136 }
137
138 if(reactionData->GetReactionType() == 0) // Totally diffused contr
139 {
140 return irt_1;
141 }
142
143 if(irt_1 < 0)
144 {
145 return irt_1;
146 }
147
148 G4double kdif = 4 * CLHEP::pi * D * SmoluchowskiRadius * Avogadro;
149
150 if(pMolConfA == pMolConfB)
151 {
152 kdif /= 2;
153 }
154 G4double kact = G4IRTUtils::GetKact(kobs, kdif);
155 G4double sumOfk = kact + kdif;
156 if(sumOfk != 0)
157 {
158 G4double rateFactor = kact / sumOfk;
159 if(G4UniformRand() > rateFactor)
160 {
161 return -1.0 * ps;
162 }
163 G4double Y = std::abs(G4RandGauss::shoot(0.0, std::sqrt(2)));
164
165 if(Y > 0)
166 {
167 X = -(G4Log(G4UniformRand())) / Y;
168 }
169 G4double f = X * SmoluchowskiRadius * kdif / sumOfk;
170 G4double irt_2 = (f * f) / D;
171 return irt_1 + irt_2;
172 }
173 }
174 return -1.0 * ps;
175}
G4double D(G4double temp)
G4double Y(G4double density)
G4double G4Log(G4double x)
Definition G4Log.hh:227
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
static G4DNAMolecularReactionTable * Instance()
static G4double erfcInv(G4double x)
static G4double GetKact(const G4double &obs, const G4double &dif)
Definition G4IRTUtils.hh:41
const G4MolecularConfiguration * GetMolecularConfiguration() const
const G4ThreeVector & GetPosition() const

◆ Initialise()

void G4DiffusionControlledReactionModel::Initialise ( const G4MolecularConfiguration * pMolecule,
const G4Track &  )
overridevirtual

Reimplemented from G4VDNAReactionModel.

Definition at line 44 of file G4DiffusionControlledReactionModel.cc.

46{
47 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
48}

◆ InitialiseToPrint()

void G4DiffusionControlledReactionModel::InitialiseToPrint ( const G4MolecularConfiguration * pMolecule)
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 50 of file G4DiffusionControlledReactionModel.cc.

52{
53 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
54}

◆ operator=()

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

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