Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 65695 2012-11-27 11:39:12Z gunter $
27 //
29 #include "Randomize.hh"
30 #include "G4Track.hh"
32 
34 {
35  fReactionData = 0 ;
36 }
37 
39  G4VDNAReactionModel(__right)
40 {
41  fReactionData = 0 ;
42 }
43 
44 G4DNASmoluchowskiReactionModel& G4DNASmoluchowskiReactionModel::operator=(const G4DNASmoluchowskiReactionModel& right)
45 {
46  if(this == &right) return *this;
47  fReactionData = 0;
48  return *this;
49 }
50 
52 {
53  fReactionData = 0 ;
54 }
55 
57 {
58  fReactionData = fReactionTable->GetReactionData(__molecule);
59 }
60 
62 {
63  fReactionData = fReactionTable->GetReactionData(__molecule);
64 }
65 
67  const G4Molecule* mol2)
68 {
69  G4double __output = fReactionTable -> GetReactionData(mol1,mol2)->GetReducedReactionRadius();
70  return __output ;
71 }
72 
74 {
75  G4double __output = (*fReactionData)[__i] -> GetReducedReactionRadius();
76  return __output ;
77 }
78 
80  const G4Track& __trackB,
81  const G4double __R,
82  G4double& __r,
83  const G4bool __alongStepReaction)
84 {
85  G4double __postStepSeparation = 0;
86  bool __do_break = false ;
87  G4double __R2 = __R*__R ;
88  int k = 0 ;
89 
90  for(; k < 3 ; k++)
91  {
92  __postStepSeparation += std::pow(__trackA.GetPosition()[k] - __trackB.GetPosition()[k],2);
93 
94  if(__postStepSeparation > __R2)
95  {
96  __do_break = true ;
97  break ;
98  }
99  }
100 
101  if(__do_break == false)
102  {
103  // The loop was not break
104  // => __r^2 < __R^2
105  __r = std::sqrt(__postStepSeparation);
106  return true;
107  }
108  else if(__alongStepReaction == true)
109  {
110  //G4cout << "alongStepReaction==true" << G4endl;
111  //Along step cheack and
112  // the loop has break
113 
114  // Continue loop
115  for(; k < 3 ; k++)
116  {
117  __postStepSeparation += std::pow(__trackA.GetPosition()[k] - __trackB.GetPosition()[k],2);
118  }
119  // Use Green approach : the Brownian bridge
120  __r = (__postStepSeparation = std::sqrt(__postStepSeparation) );
121 
122  G4Molecule* __moleculeA = GetMolecule(__trackA);
123  G4Molecule* __moleculeB = GetMolecule(__trackB);
124 
125  G4double __D = __moleculeA->GetDiffusionCoefficient() + __moleculeB->GetDiffusionCoefficient();
126 
127  G4ThreeVector __preStepPositionA = __trackA.GetStep()->GetPreStepPoint() ->GetPosition();
128  G4ThreeVector __preStepPositionB = __trackB.GetStep()->GetPreStepPoint() ->GetPosition();
129 
130  if(__preStepPositionA == __trackA.GetPosition())
131  {
132  G4ExceptionDescription exceptionDescription ;
133  exceptionDescription << "The molecule : " << __moleculeA->GetName();
134  exceptionDescription << " did not move since the previous step ";
135  G4Exception("G4DNASmoluchowskiReactionModel::FindReaction","G4DNASmoluchowskiReactionModel",
136  FatalErrorInArgument,exceptionDescription);
137  }
138 
139  G4double __preStepSeparation = (__preStepPositionA - __preStepPositionB).mag();
140 
141  G4double __probabiltyOfEncounter = std::exp(-(__preStepSeparation - __R)*(__postStepSeparation - __R)
142  / (__D* (__trackB.GetStep()->GetDeltaTime())));
143  G4double __selectedPOE = G4UniformRand();
144 
145  if(__selectedPOE<=__probabiltyOfEncounter) return true;
146  }
147 
148  return false ;
149 }