Geant4  10.01.p02
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 85244 2014-10-27 08:24:13Z gcosmo $
27 //
29 #include "Randomize.hh"
30 #include "G4Track.hh"
32 #include "G4UnitsTable.hh"
33 
36 {
37  fReactionData = 0;
38 }
39 
41  G4VDNAReactionModel(__right)
42 {
43  fReactionData = 0;
44 }
45 
47 {
48  if (this == &right) return *this;
49  fReactionData = 0;
50  return *this;
51 }
52 
54 {
55  fReactionData = 0;
56 }
57 
59  const G4Track&)
60 {
62 }
63 
65 {
67 }
68 
70  const G4Molecule* mol2)
71 {
72  G4double __output = fReactionTable->GetReactionData(mol1, mol2)
74  return __output;
75 }
76 
78 {
79  G4double __output = (*fReactionData)[__i]->GetReducedReactionRadius();
80  return __output;
81 }
82 
84  const G4Track& __trackB,
85  const G4double __R,
86  G4double& __r,
87  const G4bool __alongStepReaction)
88 {
89  G4double __postStepSeparation = 0;
90  bool __do_break = false;
91  G4double __R2 = __R * __R;
92  int k = 0;
93 
94  for (; k < 3; k++)
95  {
96  __postStepSeparation += std::pow(
97  __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
98 
99  if (__postStepSeparation > __R2)
100  {
101  __do_break = true;
102  break;
103  }
104  }
105 
106  if (__do_break == false)
107  {
108  // The loop was not break
109  // => __r^2 < __R^2
110  __r = std::sqrt(__postStepSeparation);
111  return true;
112  }
113  else if (__alongStepReaction == true)
114  {
115  //G4cout << "alongStepReaction==true" << G4endl;
116  //Along step cheack and
117  // the loop has break
118 
119  // Continue loop
120  for (; k < 3; k++)
121  {
122  __postStepSeparation += std::pow(
123  __trackA.GetPosition()[k] - __trackB.GetPosition()[k], 2);
124  }
125  // Use Green approach : the Brownian bridge
126  __r = (__postStepSeparation = std::sqrt(__postStepSeparation));
127 
128  G4Molecule* __moleculeA = GetMolecule(__trackA);
129  G4Molecule* __moleculeB = GetMolecule(__trackB);
130 
131  G4double __D = __moleculeA->GetDiffusionCoefficient()
132  + __moleculeB->GetDiffusionCoefficient();
133 
134  G4ThreeVector __preStepPositionA = __trackA.GetStep()->GetPreStepPoint()
135  ->GetPosition();
136  G4ThreeVector __preStepPositionB = __trackB.GetStep()->GetPreStepPoint()
137  ->GetPosition();
138 
139  if (__preStepPositionA == __trackA.GetPosition())
140  {
141  G4ExceptionDescription exceptionDescription;
142  exceptionDescription << "The molecule : " << __moleculeA->GetName();
143  exceptionDescription << " with track ID :" << __trackA.GetTrackID();
144  exceptionDescription << " did not move since the previous step." << G4endl;
145  exceptionDescription << "Current position : "
146  << G4BestUnit(__trackA.GetPosition(), "Length")
147  << G4endl;
148  exceptionDescription << "Previous position : "
149  << G4BestUnit(__preStepPositionA, "Length") << G4endl;
150  G4Exception("G4DNASmoluchowskiReactionModel::FindReaction",
151  "G4DNASmoluchowskiReactionModel", FatalErrorInArgument,
152  exceptionDescription);
153  }
154 
155  G4double __preStepSeparation =
156  (__preStepPositionA - __preStepPositionB).mag();
157 
158  G4double __probabiltyOfEncounter = std::exp(
159  -(__preStepSeparation - __R) * (__postStepSeparation - __R) / (__D
160  * (__trackB.GetStep()->GetDeltaTime())));
161  G4double __selectedPOE = G4UniformRand();
162 
163  if (__selectedPOE <= __probabiltyOfEncounter) return true;
164  }
165 
166  return false;
167 }
G4VDNAReactionModel is an interface used by the G4DNAMolecularReaction process.
virtual G4bool FindReaction(const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)
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:465
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
virtual void Initialise(const G4Molecule *, const G4Track &)
This macro is defined in AddClone_def.
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:303
#define G4UniformRand()
Definition: Randomize.hh:93
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
virtual void InitialiseToPrint(const G4Molecule *)
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
G4DNASmoluchowskiReactionModel & operator=(const G4DNASmoluchowskiReactionModel &)
const G4DNAMolecularReactionData * GetReactionData(const G4Molecule *, const G4Molecule *) const
virtual G4double GetReactionRadius(const G4Molecule *, const G4Molecule *)
#define G4endl
Definition: G4ios.hh:61
Class Description The dynamic molecule holds all the data that change for a molecule It has a pointer...
Definition: G4Molecule.hh:93
double G4double
Definition: G4Types.hh:76