Geant4 9.6.0
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 &)
 
virtual void Initialise (const G4Molecule *, const G4Track &)
 
virtual void InitialiseToPrint (const G4Molecule *)
 
virtual G4double GetReactionRadius (const G4Molecule *, const G4Molecule *)
 
virtual G4double GetReactionRadius (const G4int)
 
virtual G4bool FindReaction (const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)
 
- Public Member Functions inherited from G4VDNAReactionModel
 G4VDNAReactionModel ()
 
 G4VDNAReactionModel (const G4VDNAReactionModel &)
 
virtual ~G4VDNAReactionModel ()
 
virtual void Initialise (const G4Molecule *, const G4Track &)
 
virtual void InitialiseToPrint (const G4Molecule *)=0
 
virtual G4double GetReactionRadius (const G4Molecule *, const G4Molecule *)=0
 
virtual G4double GetReactionRadius (const int)=0
 
virtual G4bool FindReaction (const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)=0
 
void SetReactionTable (const G4DNAMolecularReactionTable *)
 
const G4DNAMolecularReactionTableGetReactionTable ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDNAReactionModel
G4VDNAReactionModeloperator= (const G4VDNAReactionModel &)
 
- Protected Attributes inherited from G4VDNAReactionModel
const G4DNAMolecularReactionTablefReactionTable
 

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 58 of file G4DNASmoluchowskiReactionModel.hh.

Constructor & Destructor Documentation

◆ G4DNASmoluchowskiReactionModel() [1/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( )

Definition at line 33 of file G4DNASmoluchowskiReactionModel.cc.

34{
35 fReactionData = 0 ;
36}

◆ ~G4DNASmoluchowskiReactionModel()

G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
virtual

Definition at line 51 of file G4DNASmoluchowskiReactionModel.cc.

52{
53 fReactionData = 0 ;
54}

◆ G4DNASmoluchowskiReactionModel() [2/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( const G4DNASmoluchowskiReactionModel __right)

Definition at line 38 of file G4DNASmoluchowskiReactionModel.cc.

38 :
39 G4VDNAReactionModel(__right)
40{
41 fReactionData = 0 ;
42}

Member Function Documentation

◆ FindReaction()

G4bool G4DNASmoluchowskiReactionModel::FindReaction ( const G4Track __trackA,
const G4Track __trackB,
const G4double  __R,
G4double __r,
const G4bool  __alongStepReaction 
)
virtual

Implements G4VDNAReactionModel.

Definition at line 79 of file G4DNASmoluchowskiReactionModel.cc.

84{
85 G4double __postStepSeparation = 0;
86 bool __do_break = false ;
87 G4double __R2 = __R*__R ;
88 int k = 0 ;
89
90 for(; k < 3 ; k++)
91 {
92 __postStepSeparation += std::pow(__trackA.GetPosition()[k] - __trackB.GetPosition()[k],2);
93
94 if(__postStepSeparation > __R2)
95 {
96 __do_break = true ;
97 break ;
98 }
99 }
100
101 if(__do_break == false)
102 {
103 // The loop was not break
104 // => __r^2 < __R^2
105 __r = std::sqrt(__postStepSeparation);
106 return true;
107 }
108 else if(__alongStepReaction == true)
109 {
110 //G4cout << "alongStepReaction==true" << G4endl;
111 //Along step cheack and
112 // the loop has break
113
114 // Continue loop
115 for(; k < 3 ; k++)
116 {
117 __postStepSeparation += std::pow(__trackA.GetPosition()[k] - __trackB.GetPosition()[k],2);
118 }
119 // Use Green approach : the Brownian bridge
120 __r = (__postStepSeparation = std::sqrt(__postStepSeparation) );
121
122 G4Molecule* __moleculeA = GetMolecule(__trackA);
123 G4Molecule* __moleculeB = GetMolecule(__trackB);
124
125 G4double __D = __moleculeA->GetDiffusionCoefficient() + __moleculeB->GetDiffusionCoefficient();
126
127 G4ThreeVector __preStepPositionA = __trackA.GetStep()->GetPreStepPoint() ->GetPosition();
128 G4ThreeVector __preStepPositionB = __trackB.GetStep()->GetPreStepPoint() ->GetPosition();
129
130 if(__preStepPositionA == __trackA.GetPosition())
131 {
132 G4ExceptionDescription exceptionDescription ;
133 exceptionDescription << "The molecule : " << __moleculeA->GetName();
134 exceptionDescription << " did not move since the previous step ";
135 G4Exception("G4DNASmoluchowskiReactionModel::FindReaction","G4DNASmoluchowskiReactionModel",
136 FatalErrorInArgument,exceptionDescription);
137 }
138
139 G4double __preStepSeparation = (__preStepPositionA - __preStepPositionB).mag();
140
141 G4double __probabiltyOfEncounter = std::exp(-(__preStepSeparation - __R)*(__postStepSeparation - __R)
142 / (__D* (__trackB.GetStep()->GetDeltaTime())));
143 G4double __selectedPOE = G4UniformRand();
144
145 if(__selectedPOE<=__probabiltyOfEncounter) return true;
146 }
147
148 return false ;
149}
@ FatalErrorInArgument
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
const G4String & GetName() const
Definition: G4Molecule.cc:259
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:405
const G4ThreeVector & GetPosition() const
G4double GetDeltaTime() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
const G4Step * GetStep() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ GetReactionRadius() [1/2]

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4int  __i)
virtual

Implements G4VDNAReactionModel.

Definition at line 73 of file G4DNASmoluchowskiReactionModel.cc.

74{
75 G4double __output = (*fReactionData)[__i] -> GetReducedReactionRadius();
76 return __output ;
77}

◆ GetReactionRadius() [2/2]

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4Molecule mol1,
const G4Molecule mol2 
)
virtual

Implements G4VDNAReactionModel.

Definition at line 66 of file G4DNASmoluchowskiReactionModel.cc.

68{
69 G4double __output = fReactionTable -> GetReactionData(mol1,mol2)->GetReducedReactionRadius();
70 return __output ;
71}
const G4DNAMolecularReactionTable * fReactionTable

◆ Initialise()

void G4DNASmoluchowskiReactionModel::Initialise ( const G4Molecule ,
const G4Track  
)
virtual

This macro is defined in AddClone_def

Reimplemented from G4VDNAReactionModel.

Definition at line 56 of file G4DNASmoluchowskiReactionModel.cc.

57{
58 fReactionData = fReactionTable->GetReactionData(__molecule);
59}
const G4DNAMolecularReactionData * GetReactionData(const G4Molecule *, const G4Molecule *) const

◆ InitialiseToPrint()

void G4DNASmoluchowskiReactionModel::InitialiseToPrint ( const G4Molecule __molecule)
virtual

Implements G4VDNAReactionModel.

Definition at line 61 of file G4DNASmoluchowskiReactionModel.cc.

62{
63 fReactionData = fReactionTable->GetReactionData(__molecule);
64}

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