Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4DNASmoluchowskiReactionModel Class Reference

#include <G4DNASmoluchowskiReactionModel.hh>

Inheritance diagram for G4DNASmoluchowskiReactionModel:
G4VDNAReactionModel

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

Definition at line 34 of file G4DNASmoluchowskiReactionModel.cc.

35 {
36  fReactionData = 0 ;
37 }
G4DNASmoluchowskiReactionModel::~G4DNASmoluchowskiReactionModel ( )
virtual

Definition at line 52 of file G4DNASmoluchowskiReactionModel.cc.

53 {
54  fReactionData = 0 ;
55 }
G4DNASmoluchowskiReactionModel::G4DNASmoluchowskiReactionModel ( const G4DNASmoluchowskiReactionModel __right)

Definition at line 39 of file G4DNASmoluchowskiReactionModel.cc.

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

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

References FatalErrorInArgument, G4BestUnit, G4endl, G4Exception(), G4UniformRand, G4Step::GetDeltaTime(), G4Molecule::GetDiffusionCoefficient(), GetMolecule(), G4Molecule::GetName(), G4StepPoint::GetPosition(), G4Track::GetPosition(), G4Step::GetPreStepPoint(), G4Track::GetStep(), and G4Track::GetTrackID().

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

Implements G4VDNAReactionModel.

Definition at line 67 of file G4DNASmoluchowskiReactionModel.cc.

References G4VDNAReactionModel::fReactionTable.

69 {
70  G4double __output = fReactionTable -> GetReactionData(mol1,mol2)->GetReducedReactionRadius();
71  return __output ;
72 }
const G4DNAMolecularReactionTable * fReactionTable
double G4double
Definition: G4Types.hh:76
G4double G4DNASmoluchowskiReactionModel::GetReactionRadius ( const G4int  __i)
virtual

Implements G4VDNAReactionModel.

Definition at line 74 of file G4DNASmoluchowskiReactionModel.cc.

75 {
76  G4double __output = (*fReactionData)[__i] -> GetReducedReactionRadius();
77  return __output ;
78 }
double G4double
Definition: G4Types.hh:76
void G4DNASmoluchowskiReactionModel::Initialise ( const G4Molecule ,
const G4Track  
)
virtual

This macro is defined in AddClone_def

Reimplemented from G4VDNAReactionModel.

Definition at line 57 of file G4DNASmoluchowskiReactionModel.cc.

References G4VDNAReactionModel::fReactionTable, and G4DNAMolecularReactionTable::GetReactionData().

58 {
59  fReactionData = fReactionTable->GetReactionData(__molecule);
60 }
const G4DNAMolecularReactionTable * fReactionTable
const G4DNAMolecularReactionData * GetReactionData(const G4Molecule *, const G4Molecule *) const
void G4DNASmoluchowskiReactionModel::InitialiseToPrint ( const G4Molecule __molecule)
virtual

Implements G4VDNAReactionModel.

Definition at line 62 of file G4DNASmoluchowskiReactionModel.cc.

References G4VDNAReactionModel::fReactionTable, and G4DNAMolecularReactionTable::GetReactionData().

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

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