Geant4  10.02.p01
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 // $Id: G4DNASmoluchowskiReactionModel.cc 94218 2015-11-09 08:24:48Z gcosmo $
27 //
29 #include "Randomize.hh"
30 #include "G4Track.hh"
32 #include "G4UnitsTable.hh"
33 #include "G4Molecule.hh"
34 //#include "G4Scheduler.hh"
35 
38 {
39  fReactionData = 0;
40 }
41 
43  G4VDNAReactionModel(__right)
44 {
45  fReactionData = 0;
46 }
47 
49 {
50  if (this == &right) return *this;
51  fReactionData = 0;
52  return *this;
53 }
54 
56 {
57  fReactionData = 0;
58 }
59 
61  const G4Track&)
62 {
64 }
65 
66 void
69 {
71 }
72 
76 {
77  G4double __output = fReactionTable->GetReactionData(__mol1, __mol2)
79  return __output;
80 }
81 
83 {
84  G4double __output = (*fReactionData)[__i]->GetEffectiveReactionRadius();
85  return __output;
86 }
87 
89  const G4Track& __trackB,
90  const G4double __R,
91  G4double& __r,
92  const G4bool __alongStepReaction)
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 }
The pointer G4MolecularConfiguration will be shared by all the molecules having the same molecule def...
G4VDNAReactionModel is an interface used by the G4DNAMolecularReaction process.
virtual G4bool FindReaction(const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
G4DNASmoluchowskiReactionModel should be used for very fast reactions (high reaction rate) : the reac...
const G4ThreeVector & GetPosition() const
G4double GetDiffusionCoefficient() const
Returns the diffusion coefficient D.
Definition: G4Molecule.cc:549
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
const std::vector< const G4DNAMolecularReactionData * > * fReactionData
const G4DNAMolecularReactionTable * fReactionTable
G4StepPoint * GetPreStepPoint() const
const G4String & GetName() const
Returns the name of the molecule.
Definition: G4Molecule.cc:345
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetReactionRadius(G4MolecularConfiguration *, G4MolecularConfiguration *)
G4double GetDeltaTime() const
G4int GetTrackID() const
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
G4DNASmoluchowskiReactionModel & operator=(const G4DNASmoluchowskiReactionModel &)
#define G4endl
Definition: G4ios.hh:61
virtual void Initialise(G4MolecularConfiguration *, const G4Track &)
This macro is defined in AddClone_def.
virtual void InitialiseToPrint(G4MolecularConfiguration *)
Class Description The dynamic molecule holds all the data that change for a molecule It has a pointer...
Definition: G4Molecule.hh:94
double G4double
Definition: G4Types.hh:76