Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNASmoluchowskiReactionModel Class Reference

#include <G4DNASmoluchowskiReactionModel.hh>

Inheritance diagram for G4DNASmoluchowskiReactionModel:
Collaboration diagram for G4DNASmoluchowskiReactionModel:

Public Member Functions

 G4DNASmoluchowskiReactionModel ()
 
virtual ~G4DNASmoluchowskiReactionModel ()
 
 G4DNASmoluchowskiReactionModel (const G4DNASmoluchowskiReactionModel &)
 
virtual void Initialise (G4MolecularConfiguration *, const G4Track &)
 
virtual void InitialiseToPrint (G4MolecularConfiguration *)
 
virtual G4double GetReactionRadius (G4MolecularConfiguration *, G4MolecularConfiguration *)
 
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 ()
 
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 67 of file G4DNASmoluchowskiReactionModel.hh.

Constructor & Destructor Documentation

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( )

Definition at line 37 of file G4DNASmoluchowskiReactionModel.cc.

37  :
39 {
40  fReactionData = 0;
41 }
G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
virtual

Definition at line 56 of file G4DNASmoluchowskiReactionModel.cc.

57 {
58  fReactionData = 0;
59 }
G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( const G4DNASmoluchowskiReactionModel __right)

Definition at line 43 of file G4DNASmoluchowskiReactionModel.cc.

43  :
44  G4VDNAReactionModel(__right)
45 {
46  fReactionData = 0;
47 }

Member Function Documentation

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

Implements G4VDNAReactionModel.

Definition at line 89 of file G4DNASmoluchowskiReactionModel.cc.

94 {
95  G4double postStepSeparation = 0;
96  bool do_break = false;
97  G4double R2 = __R * __R;
98  int k = 0;
99 
100  for (; k < 3; k++)
101  {
102  postStepSeparation += std::pow(
103  __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
104 
105  if (postStepSeparation > R2)
106  {
107  do_break = true;
108  break;
109  }
110  }
111 
112  if (do_break == false)
113  {
114  // The loop was not break
115  // => __r^2 < __R^2
116  __r = std::sqrt(postStepSeparation);
117  return true;
118  }
119  else if (__alongStepReaction == true)
120  {
121  //G4cout << "alongStepReaction==true" << G4endl;
122  //Along step cheack and
123  // the loop has break
124 
125  // Continue loop
126  for (; k < 3; k++)
127  {
128  postStepSeparation += std::pow(
129  __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
130  }
131  // Use Green approach : the Brownian bridge
132  __r = (postStepSeparation = std::sqrt(postStepSeparation));
133 
134  G4Molecule* __moleculeA = GetMolecule(__trackA);
135  G4Molecule* __moleculeB = GetMolecule(__trackB);
136 
137  G4double __D = __moleculeA->GetDiffusionCoefficient()
138  + __moleculeB->GetDiffusionCoefficient();
139 
140  G4ThreeVector __preStepPositionA = __trackA.GetStep()->GetPreStepPoint()
141  ->GetPosition();
142  G4ThreeVector __preStepPositionB = __trackB.GetStep()->GetPreStepPoint()
143  ->GetPosition();
144 
145  if (__preStepPositionA == __trackA.GetPosition())
146  {
147  G4ExceptionDescription exceptionDescription;
148  exceptionDescription << "The molecule : " << __moleculeA->GetName();
149  exceptionDescription << " with track ID :" << __trackA.GetTrackID();
150  exceptionDescription << " did not move since the previous step." << G4endl;
151  exceptionDescription << "Current position : "
152  << G4BestUnit(__trackA.GetPosition(), "Length")
153  << G4endl;
154  exceptionDescription << "Previous position : "
155  << G4BestUnit(__preStepPositionA, "Length") << G4endl;
156  G4Exception("G4DNASmoluchowskiReactionModel::FindReaction",
157  "G4DNASmoluchowskiReactionModel", FatalErrorInArgument,
158  exceptionDescription);
159  }
160 
161  G4double __preStepSeparation =
162  (__preStepPositionA - __preStepPositionB).mag();
163 
164  //===================================
165  // Brownian bridge
166 
167 
168 // if(G4Scheduler::Instance()->GetTimeStep() != __trackB.GetStep()->GetDeltaTime())
169 // {
170 // G4cout << G4Scheduler::Instance()->GetTimeStep() << G4endl;
171 // G4cout << __trackB.GetStep()->GetDeltaTime() << G4endl;
172 // assert(G4Scheduler::Instance()->GetTimeStep() == __trackB.GetStep()->GetDeltaTime());
173 // }
174 
175  G4double __probabiltyOfEncounter = G4Exp(
176  -(__preStepSeparation - __R) * (postStepSeparation - __R) / (__D
177  * (__trackB.GetStep()->GetDeltaTime())));
178  G4double __selectedPOE = G4UniformRand();
179 
180  if (__selectedPOE <= __probabiltyOfEncounter) return true;
181  //===================================
182  }
183 
184  return false;
185 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4ThreeVector & GetPosition() const
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:560
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4Step * GetStep() const
G4StepPoint * GetPreStepPoint() const
const G4String & GetName() const
Definition: G4Molecule.cc:356
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetPosition() const
G4double GetDeltaTime() const
G4int GetTrackID() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( G4MolecularConfiguration __mol1,
G4MolecularConfiguration __mol2 
)
virtual

Implements G4VDNAReactionModel.

Definition at line 75 of file G4DNASmoluchowskiReactionModel.cc.

77 {
78  G4double __output = fReactionTable->GetReactionData(__mol1, __mol2)
80  return __output;
81 }
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
const G4DNAMolecularReactionTable * fReactionTable
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4int  __i)
virtual

Implements G4VDNAReactionModel.

Definition at line 83 of file G4DNASmoluchowskiReactionModel.cc.

84 {
85  G4double __output = (*fReactionData)[__i]->GetEffectiveReactionRadius();
86  return __output;
87 }
double G4double
Definition: G4Types.hh:76
void G4DNASmoluchowskiReactionModel::Initialise ( G4MolecularConfiguration ,
const G4Track  
)
virtual

This macro is defined in AddClone_def

Reimplemented from G4VDNAReactionModel.

Definition at line 61 of file G4DNASmoluchowskiReactionModel.cc.

63 {
64  fReactionData = fReactionTable->GetReactionData(__molecule);
65 }
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
const G4DNAMolecularReactionTable * fReactionTable

Here is the call graph for this function:

void G4DNASmoluchowskiReactionModel::InitialiseToPrint ( G4MolecularConfiguration __molecule)
virtual

Implements G4VDNAReactionModel.

Definition at line 69 of file G4DNASmoluchowskiReactionModel.cc.

70 {
71  fReactionData = fReactionTable->GetReactionData(__molecule);
72 }
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
const G4DNAMolecularReactionTable * fReactionTable

Here is the call graph for this function:


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