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

#include <G4DNASmoluchowskiReactionModel.hh>

+ Inheritance diagram for G4DNASmoluchowskiReactionModel:

Public Member Functions

 G4DNASmoluchowskiReactionModel ()
 
 ~G4DNASmoluchowskiReactionModel () override
 
 G4DNASmoluchowskiReactionModel (const G4DNASmoluchowskiReactionModel &)=delete
 
G4DNASmoluchowskiReactionModeloperator= (const G4DNASmoluchowskiReactionModel &)=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
 
- 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

G4DNASmoluchowskiReactionModel should be used for very fast reactions (high reaction rate) : the reactions between reactants occuring at encounter. When the time step is constrained this model uses brownian bridge : "Absorbing (Smoluchowski) boundary condition" Reference : On the simulation of diffusion processes close to boundaries, N. J. B. Green, Molecular Physics, 65: 6, 1399 — 1408(1988)

Definition at line 63 of file G4DNASmoluchowskiReactionModel.hh.

Constructor & Destructor Documentation

◆ G4DNASmoluchowskiReactionModel() [1/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( )
default

◆ ~G4DNASmoluchowskiReactionModel()

G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
overridedefault

◆ G4DNASmoluchowskiReactionModel() [2/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( const G4DNASmoluchowskiReactionModel & )
delete

Member Function Documentation

◆ FindReaction()

G4bool G4DNASmoluchowskiReactionModel::FindReaction ( const G4Track & __trackA,
const G4Track & __trackB,
G4double __reactionRadius,
G4double & __separationDistance,
G4bool __alongStepReaction )
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 63 of file G4DNASmoluchowskiReactionModel.cc.

68{
69 const G4double R2 = __reactionRadius * __reactionRadius;
70 G4double postStepSeparation = 0.;
71 bool do_break = false;
72 int k = 0;
73
74 for (; k < 3; ++k)
75 {
76 postStepSeparation += std::pow(
77 __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
78
79 if (postStepSeparation > R2)
80 {
81 do_break = true;
82 break;
83 }
84 }
85
86 if (!do_break)
87 {
88 // The loop was not break
89 // => r^2 < R^2
90 __separationDistance = std::sqrt(postStepSeparation);
91 return true;
92 }
93 if (__alongStepReaction)
94 {
95 //Along step check and the loop has break
96
97 // Continue loop
98 for (; k < 3; ++k)
99 {
100 postStepSeparation += std::pow(
101 __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
102 }
103 // Use Green approach : the Brownian bridge
104 __separationDistance = (postStepSeparation = std::sqrt(postStepSeparation));
105
106 auto pMoleculeA = GetMolecule(__trackA);
107 auto pMoleculeB = GetMolecule(__trackB);
108
109 G4double D = pMoleculeA->GetDiffusionCoefficient()
110 + pMoleculeB->GetDiffusionCoefficient();
111
112 const auto& preStepPositionA = __trackA.GetStep()->GetPreStepPoint()->GetPosition();
113 const auto& preStepPositionB = __trackB.GetStep()->GetPreStepPoint()->GetPosition();
114
115 G4double preStepSeparation = (preStepPositionA - preStepPositionB).mag();
116
117 //===================================
118 // Brownian bridge
119 G4double __probabiltyOfEncounter = G4Exp(-(preStepSeparation - __reactionRadius)
120 * (postStepSeparation - __reactionRadius)
121 / (D * (__trackB.GetStep()->GetDeltaTime()))
122 );
123 G4double __selectedPOE = G4UniformRand();
124
125 if (__selectedPOE <= __probabiltyOfEncounter)
126 {
127 return true;
128 }
129 //===================================
130 }
131
132 return false;
133}
G4double D(G4double temp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
const G4ThreeVector & GetPosition() const
G4double GetDeltaTime() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
const G4Step * GetStep() const

◆ GetReactionRadius() [1/2]

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4int & __i)
overridevirtual

Implements G4VDNAReactionModel.

Definition at line 57 of file G4DNASmoluchowskiReactionModel.cc.

58{
59 G4double __output = (*fpReactionData)[__i]->GetEffectiveReactionRadius();
60 return __output;
61}

◆ GetReactionRadius() [2/2]

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

Implements G4VDNAReactionModel.

Definition at line 50 of file G4DNASmoluchowskiReactionModel.cc.

52{
54 return __output;
55}
Data * GetReactionData(Reactant *, Reactant *) const
const G4DNAMolecularReactionTable * fpReactionTable

◆ Initialise()

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

Reimplemented from G4VDNAReactionModel.

Definition at line 39 of file G4DNASmoluchowskiReactionModel.cc.

41{
42 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
43}

◆ InitialiseToPrint()

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

Implements G4VDNAReactionModel.

Definition at line 45 of file G4DNASmoluchowskiReactionModel.cc.

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

◆ operator=()

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

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