Geant4  10.02.p03
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 ()
 

Private Member Functions

G4DNASmoluchowskiReactionModeloperator= (const G4DNASmoluchowskiReactionModel &)
 

Private Attributes

const std::vector< const G4DNAMolecularReactionData * > * fReactionData
 

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() [1/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( )

Definition at line 36 of file G4DNASmoluchowskiReactionModel.cc.

36  :
38 {
39  fReactionData = 0;
40 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData

◆ ~G4DNASmoluchowskiReactionModel()

G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
virtual

Definition at line 55 of file G4DNASmoluchowskiReactionModel.cc.

56 {
57  fReactionData = 0;
58 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData

◆ G4DNASmoluchowskiReactionModel() [2/2]

G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( const G4DNASmoluchowskiReactionModel __right)

Definition at line 42 of file G4DNASmoluchowskiReactionModel.cc.

42  :
43  G4VDNAReactionModel(__right)
44 {
45  fReactionData = 0;
46 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData

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 88 of file G4DNASmoluchowskiReactionModel.cc.

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

◆ GetReactionRadius() [1/2]

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

Implements G4VDNAReactionModel.

Definition at line 74 of file G4DNASmoluchowskiReactionModel.cc.

76 {
77  G4double __output = fReactionTable->GetReactionData(__mol1, __mol2)
79  return __output;
80 }
const G4DNAMolecularReactionTable * fReactionTable
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetReactionRadius() [2/2]

G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4int  __i)
virtual

Implements G4VDNAReactionModel.

Definition at line 82 of file G4DNASmoluchowskiReactionModel.cc.

83 {
84  G4double __output = (*fReactionData)[__i]->GetEffectiveReactionRadius();
85  return __output;
86 }
double G4double
Definition: G4Types.hh:76

◆ Initialise()

void G4DNASmoluchowskiReactionModel::Initialise ( G4MolecularConfiguration __molecule,
const G4Track &   
)
virtual

Reimplemented from G4VDNAReactionModel.

Definition at line 60 of file G4DNASmoluchowskiReactionModel.cc.

62 {
64 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData
const G4DNAMolecularReactionTable * fReactionTable
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
Here is the call graph for this function:

◆ InitialiseToPrint()

void G4DNASmoluchowskiReactionModel::InitialiseToPrint ( G4MolecularConfiguration __molecule)
virtual

Implements G4VDNAReactionModel.

Definition at line 68 of file G4DNASmoluchowskiReactionModel.cc.

69 {
71 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData
const G4DNAMolecularReactionTable * fReactionTable
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4DNASmoluchowskiReactionModel & G4DNASmoluchowskiReactionModel::operator= ( const G4DNASmoluchowskiReactionModel right)
private

Definition at line 48 of file G4DNASmoluchowskiReactionModel.cc.

49 {
50  if (this == &right) return *this;
51  fReactionData = 0;
52  return *this;
53 }
const std::vector< const G4DNAMolecularReactionData * > * fReactionData

Member Data Documentation

◆ fReactionData

const std::vector<const G4DNAMolecularReactionData*>* G4DNASmoluchowskiReactionModel::fReactionData
private

Definition at line 90 of file G4DNASmoluchowskiReactionModel.hh.


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