Geant4-11
G4DNATotallyDiffusionControlled.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// 20/2/2019
27// Author : HoangTRAN
28
30#include "G4IRTUtils.hh"
33#include "G4OctreeFinder.hh"
34#include "G4SystemOfUnits.hh"
37#include "Randomize.hh"
38#include "G4Molecule.hh"
39#include "G4Electron_aq.hh"
40#include "G4Hydrogen.hh"
41#include "G4ErrorFunction.hh"
44{}
46
48 const G4Track& trackB)
49{
50 auto pMolConfA = GetMolecule(trackA)->GetMolecularConfiguration();
51 auto pMolConfB = GetMolecule(trackB)->GetMolecularConfiguration();
52
53 G4double D = GetDiffusionCoefficient(pMolConfA, pMolConfB);
54 auto reactionData = G4DNAMolecularReactionTable::Instance()
55 ->GetReactionData(pMolConfA, pMolConfB);
56 G4double kobs = reactionData->GetObservedReactionRateConstant();
57 G4double distance = (trackA.GetPosition() - trackB.GetPosition()).mag();
58 G4double Reff = kobs / ( 4 * CLHEP::pi * D * Avogadro );
59
60 if( distance < Reff )
61 {
62 G4ExceptionDescription exceptionDescription;
63 exceptionDescription << "distance = "<< distance
64 << " is uncorrected with "
65 <<" Reff = "<<Reff
66 <<" for : "<<pMolConfA->GetName()
67 <<" and "<<pMolConfB->GetName();
68 G4Exception("G4DNATotallyDiffusionControlled"
69 "::GetTimeToEncounter()", "G4DNATotallyDiffusionControlled02",
70 FatalException, exceptionDescription);
71 }
72
73 G4double Winf = Reff / distance;
75 G4double irt = -1.0 * ps;
76
77 if ( U < Winf )
78 {
79 G4double d = ( distance - Reff ) /
80 G4ErrorFunction::erfcInv( U / Winf );
81 irt = ( 1.0 / ( 4 * D ) ) * d * d;
82 }
83 return irt;
84}
87 const G4MolecularConfiguration* pMolB)
88{
91 {
92 G4bool spinA;
93 G4bool spinB;
94 spinA = G4UniformRand() < 0.5;
95 if(spinA &&
98 {
99 spinB = G4UniformRand() < 0.5;
100 if( !spinB )
101 {
102 return true;
103 }
104 }
105 return false;
106 }
107 return true;
108}
109
112 const G4MolecularConfiguration* mB)
113{
114 G4double D;
115 if(mA == mB)
116 {
117 D = (mA->GetDiffusionCoefficient());
118 }
119 else
120 {
121 D = (mA->GetDiffusionCoefficient() +
123 }
124 return D;
125}
126
127
128
129
G4double D(G4double temp)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
static constexpr double ps
Definition: G4SIunits.hh:157
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
Data * GetReactionData(Reactant *, Reactant *) const
static G4DNAMolecularReactionTable * Instance()
G4double GetDiffusionCoefficient(const G4MolecularConfiguration *pMA, const G4MolecularConfiguration *pMB)
G4bool GeminateRecombinationProbability(const G4MolecularConfiguration *pConfMolA, const G4MolecularConfiguration *pConfMolB) override
G4double GetTimeToEncounter(const G4Track &trackA, const G4Track &trackB) override
static G4Electron_aq * Definition()
static G4double erfcInv(G4double x)
static G4Hydrogen * Definition()
Definition: G4Hydrogen.cc:45
const G4MoleculeDefinition * GetDefinition() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
const G4ThreeVector & GetPosition() const
static constexpr double pi
Definition: SystemOfUnits.h:55
float Avogadro
Definition: hepunit.py:252