Geant4 11.1.1
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 ()
 
virtual ~G4DNASmoluchowskiReactionModel ()
 
 G4DNASmoluchowskiReactionModel (const G4DNASmoluchowskiReactionModel &)=delete
 
G4DNASmoluchowskiReactionModeloperator= (const G4DNASmoluchowskiReactionModel &)=delete
 
virtual void Initialise (const G4MolecularConfiguration *, const G4Track &)
 
virtual void InitialiseToPrint (const G4MolecularConfiguration *)
 
virtual G4double GetReactionRadius (const G4MolecularConfiguration *, const G4MolecularConfiguration *)
 
virtual G4double GetReactionRadius (const G4int &)
 
virtual G4bool FindReaction (const G4Track &, const G4Track &, G4double, G4double &, G4bool)
 
- 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

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 ( )

Definition at line 35 of file G4DNASmoluchowskiReactionModel.cc.

37 , fpReactionData(nullptr)
38{
39}

◆ ~G4DNASmoluchowskiReactionModel()

G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
virtualdefault

◆ 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 
)
virtual

Implements G4VDNAReactionModel.

Definition at line 67 of file G4DNASmoluchowskiReactionModel.cc.

72{
73 const G4double R2 = __reactionRadius * __reactionRadius;
74 G4double postStepSeparation = 0.;
75 bool do_break = false;
76 int k = 0;
77
78 for (; k < 3; ++k)
79 {
80 postStepSeparation += std::pow(
81 __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
82
83 if (postStepSeparation > R2)
84 {
85 do_break = true;
86 break;
87 }
88 }
89
90 if (do_break == false)
91 {
92 // The loop was not break
93 // => r^2 < R^2
94 __separationDistance = std::sqrt(postStepSeparation);
95 return true;
96 }
97 else if (__alongStepReaction == true)
98 {
99 //Along step check and the loop has break
100
101 // Continue loop
102 for (; k < 3; ++k)
103 {
104 postStepSeparation += std::pow(
105 __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
106 }
107 // Use Green approach : the Brownian bridge
108 __separationDistance = (postStepSeparation = std::sqrt(postStepSeparation));
109
110 auto pMoleculeA = GetMolecule(__trackA);
111 auto pMoleculeB = GetMolecule(__trackB);
112
113 G4double D = pMoleculeA->GetDiffusionCoefficient()
114 + pMoleculeB->GetDiffusionCoefficient();
115
116 const auto& preStepPositionA = __trackA.GetStep()->GetPreStepPoint()->GetPosition();
117 const auto& preStepPositionB = __trackB.GetStep()->GetPreStepPoint()->GetPosition();
118
119 G4double preStepSeparation = (preStepPositionA - preStepPositionB).mag();
120
121 //===================================
122 // Brownian bridge
123 G4double __probabiltyOfEncounter = G4Exp(-(preStepSeparation - __reactionRadius)
124 * (postStepSeparation - __reactionRadius)
125 / (D * (__trackB.GetStep()->GetDeltaTime()))
126 );
127 G4double __selectedPOE = G4UniformRand();
128
129 if (__selectedPOE <= __probabiltyOfEncounter)
130 {
131 return true;
132 }
133 //===================================
134 }
135
136 return false;
137}
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)
virtual

Implements G4VDNAReactionModel.

Definition at line 61 of file G4DNASmoluchowskiReactionModel.cc.

62{
63 G4double __output = (*fpReactionData)[__i]->GetEffectiveReactionRadius();
64 return __output;
65}

◆ GetReactionRadius() [2/2]

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4MolecularConfiguration pMol1,
const G4MolecularConfiguration pMol2 
)
virtual

Implements G4VDNAReactionModel.

Definition at line 54 of file G4DNASmoluchowskiReactionModel.cc.

56{
58 return __output;
59}
Data * GetReactionData(Reactant *, Reactant *) const
const G4DNAMolecularReactionTable * fpReactionTable

◆ Initialise()

void G4DNASmoluchowskiReactionModel::Initialise ( const G4MolecularConfiguration pMolecule,
const G4Track  
)
virtual

Reimplemented from G4VDNAReactionModel.

Definition at line 43 of file G4DNASmoluchowskiReactionModel.cc.

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

◆ InitialiseToPrint()

void G4DNASmoluchowskiReactionModel::InitialiseToPrint ( const G4MolecularConfiguration pMolecule)
virtual

Implements G4VDNAReactionModel.

Definition at line 49 of file G4DNASmoluchowskiReactionModel.cc.

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

◆ operator=()

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

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