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
G4DNASmoluchowskiReactionModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
28#include "Randomize.hh"
29#include "G4Track.hh"
31#include "G4UnitsTable.hh"
32#include "G4Molecule.hh"
33#include "G4Exp.hh"
34
37 , fpReactionData(nullptr)
38{
39}
40
42
44 const G4Track&)
45{
46 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
47}
48
50{
51 fpReactionData = fpReactionTable->GetReactionData(pMolecule);
52}
53
55 const G4MolecularConfiguration* pMol2)
56{
58 return __output;
59}
60
62{
63 G4double __output = (*fpReactionData)[__i]->GetEffectiveReactionRadius();
64 return __output;
65}
66
68 const G4Track& __trackB,
69 const G4double __reactionRadius,
70 G4double& __separationDistance,
71 const G4bool __alongStepReaction)
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
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Data * GetReactionData(Reactant *, Reactant *) const
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)
virtual void Initialise(const G4MolecularConfiguration *, const G4Track &)
virtual void InitialiseToPrint(const G4MolecularConfiguration *)
virtual G4bool FindReaction(const G4Track &, const G4Track &, G4double, G4double &, G4bool)
const G4ThreeVector & GetPosition() const
G4double GetDeltaTime() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
const G4Step * GetStep() const
const G4DNAMolecularReactionTable * fpReactionTable