Geant4-11
G4DiffusionControlledReactionModel.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25//
26// Author: Hoang TRAN
27
29#include "G4Track.hh"
32#include "G4Exp.hh"
33#include "G4IRTUtils.hh"
34#include "G4SystemOfUnits.hh"
36#include "G4Electron_aq.hh"
39 , fpReactionData(nullptr)
40 , fReactionTypeManager(nullptr)
41{
42}
43
45
47 const G4Track&)
48{
50}
51
53{
55}
56
58 const G4MolecularConfiguration* pMol2)
59{
60 auto reactionData = fpReactionTable->GetReactionData(pMol1, pMol2);
61 if(reactionData == nullptr)
62 {
63 G4ExceptionDescription exceptionDescription;
64 exceptionDescription << "No reactionData"
65 <<" for : "<<pMol1->GetName()
66 <<" and "<<pMol2->GetName();
67 G4Exception("G4DiffusionControlledReactionModel"
68 "::GetReactionRadius()", "G4DiffusionControlledReactionModel00",
69 FatalException, exceptionDescription);
70 }
71 G4double kobs = reactionData->GetObservedReactionRateConstant();
72 G4double D;
73 if(pMol1 == pMol2)
74 {
75 D = (pMol1->GetDiffusionCoefficient());
76 }
77 else
78 {
79 D = (pMol1->GetDiffusionCoefficient() +
80 pMol2->GetDiffusionCoefficient());//
81 }
82
83 if ( D == 0 )
84 {
85 G4ExceptionDescription exceptionDescription;
86 exceptionDescription << "D = "<< D
87 << " is uncorrected"
88 <<" for : "<<pMol1->GetName()
89 <<" and "<<pMol2->GetName();
90 G4Exception("G4DiffusionControlledReactionModel"
91 "::GetReactionRadius()", "G4DiffusionControlledReactionModel01",
92 FatalException, exceptionDescription);
93 }
94
95 G4double Reff = kobs / ( 4 * CLHEP::pi * D * Avogadro );
96 return Reff;
97}
98
100{
101 auto pMol1 = (*fpReactionData)[i]->GetReactant1();
102 auto pMol2 = (*fpReactionData)[i]->GetReactant2();
103 return GetReactionRadius(pMol1,pMol2);
104}
105
107{
109}
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
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
Data * GetReactionData(Reactant *, Reactant *) const
void InitialiseToPrint(const G4MolecularConfiguration *) override
const std::vector< const G4DNAMolecularReactionData * > * fpReactionData
G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *) override
void Initialise(const G4MolecularConfiguration *, const G4Track &) override
void SetReactionTypeManager(G4VReactionTypeManager *typeManager)
const G4String & GetName() const
const G4DNAMolecularReactionTable * fpReactionTable
static constexpr double pi
Definition: SystemOfUnits.h:55
float Avogadro
Definition: hepunit.py:252