Geant4 11.1.1
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 ()
 
virtual void Initialise (const G4MolecularConfiguration *, const G4Track &)
 
virtual void InitialiseToPrint (const G4MolecularConfiguration *)=0
 
virtual G4double GetReactionRadius (const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0
 
virtual G4double GetReactionRadius (const G4int &)=0
 
virtual G4bool FindReaction (const G4Track &, const G4Track &, G4double, G4double &, G4bool)=0
 
void SetReactionTable (const G4DNAMolecularReactionTable *)
 
const G4DNAMolecularReactionTableGetReactionTable ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VDNAReactionModel
const G4DNAMolecularReactionTablefpReactionTable
 

Detailed Description

Definition at line 35 of file G4DiffusionControlledReactionModel.hh.

Constructor & Destructor Documentation

◆ G4DiffusionControlledReactionModel() [1/2]

G4DiffusionControlledReactionModel::G4DiffusionControlledReactionModel ( )

◆ ~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 80 of file G4DiffusionControlledReactionModel.cc.

81{
82 auto pMol1 = (*fpReactionData)[i]->GetReactant1();
83 auto pMol2 = (*fpReactionData)[i]->GetReactant2();
84 return GetReactionRadius(pMol1, pMol2);
85}
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 59 of file G4DiffusionControlledReactionModel.cc.

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}
@ 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
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 87 of file G4DiffusionControlledReactionModel.cc.

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)
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
Definition: G4Molecule.cc:530
const G4ThreeVector & GetPosition() const

◆ Initialise()

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

Reimplemented from G4VDNAReactionModel.

Definition at line 47 of file G4DiffusionControlledReactionModel.cc.

49{
50 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
51}

◆ InitialiseToPrint()

void G4DiffusionControlledReactionModel::InitialiseToPrint ( const G4MolecularConfiguration pMolecule)
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 53 of file G4DiffusionControlledReactionModel.cc.

55{
56 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
57}

◆ operator=()

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

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