49static double onsager_constant = e_squared / (4. *
pi * epsilon0 * k_Boltzmann);
55double Y(
double density)
57 return 1. / (1. + 0.0012 / (density * density));
60double A(
double temperature)
62 double temp_inverse = 1 / temperature;
64 + 642.0 * temp_inverse
65 - 1.167e5 * temp_inverse * temp_inverse
66 + 9.190e6 * temp_inverse * temp_inverse * temp_inverse;
69double B(
double temperature)
71 double temp_inverse = 1 / temperature;
73 + 275.4 * temp_inverse
74 + 0.3245e5 * temp_inverse * temp_inverse;
79 double temp_inverse = 1 / temp;
82 - 11.41 * temp_inverse
83 - 35260.0 * temp_inverse * temp_inverse;
88 return A(temp) -
B(temp) - 3;
96double epsilon(
double density,
double temperature)
98 return 1 +
G4Exp(std::log(10.) *
100 (
C(temperature) + (
S(temperature) - 1) * std::log(density) / std::log(10.))
101 +
D(temperature) + std::log(density) / std::log(10.)));
128 fIsInitialized =
false;
130 fpMoleculeDensity =
nullptr;
145 return &fParticleChange;
159void G4DNAElectronHoleRecombination::MakeReaction(
const G4Track& track)
163 double random = pState->fSampleProba;
164 std::vector<ReactantInfo>& reactants = pState->fReactants;
166 G4Track* pSelectedReactant =
nullptr;
168 for (
const auto& reactantInfo : reactants)
170 if (reactantInfo.fElectron->GetTrackStatus() !=
fAlive)
174 if (reactantInfo.fProbability > random)
176 pSelectedReactant = reactantInfo.fElectron;
181 if (pSelectedReactant)
186 RemoveAMoleculeAtTime(
GetMolecule(track)->GetMolecularConfiguration(),
195 AddAMoleculeAtTime(
GetMolecule(track)->GetMolecularConfiguration(),
215G4bool G4DNAElectronHoleRecombination::FindReactant(
const G4Track& track)
223 const auto pDensityTable =
227 double density = (*pDensityTable)[track.
GetMaterial()->
GetIndex()] / (g / (1e-2 * m * 1e-2 * m * 1e-2 * m));
228 double eps =
epsilon(density, temperature);
230 double onsagerRadius = onsager_constant * 1. / (temperature * eps);
236 e_aq.GetMoleculeID(),
237 10. * onsagerRadius);
239 if (results == 0 || results->GetSize() == 0)
247 std::vector<ReactantInfo>& reactants = pState->fReactants;
250 reactants.resize(results->GetSize());
252 for (
size_t i = 0; !results->End(); results->Next(), ++i)
254 reactants[i].fElectron = results->GetItem<
G4IT>()->GetTrack();
255 reactants[i].fDistance = std::sqrt(results->GetDistanceSqr());
257 if (reactants[i].fDistance != 0)
259 reactants[i].fProbability = 1. -
G4Exp(-onsagerRadius / reactants[i].fDistance);
263 reactants[i].fProbability = 1.;
267 return reactants.empty() ? false : reactants[0].fProbability > pState->fSampleProba;
272void G4DNAElectronHoleRecombination::BuildDissociationChannels()
276 auto pH2ODef = pMoleculeTable->GetMoleculeDefinition(
"H2O",
false);
278 if (pH2ODef ==
nullptr)
285 assert(pH2Ovib !=
nullptr);
287 const auto pH2 = pMoleculeTable->GetConfiguration(
"H2",
false);
288 const auto pOH = pMoleculeTable->GetConfiguration(
"OH",
false);
289 const auto pH = pMoleculeTable->GetConfiguration(
"H",
false);
291 double probaRemaining = 1.;
295 auto pDissocationChannel2 =
299 pDissocationChannel2->AddProduct(pH2);
303 pDissocationChannel2->AddProduct(pOH);
304 pDissocationChannel2->AddProduct(pOH);
307 double channelProbability = 0.15;
308 pDissocationChannel2->SetProbability(channelProbability);
309 probaRemaining -= channelProbability;
311 B1A1_DissociationDecay);
312 pH2ODef->AddDecayChannel(pH2Ovib, pDissocationChannel2);
317 auto pDissociationChannel3 =
321 pDissociationChannel3->AddProduct(pOH);
325 pDissociationChannel3->AddProduct(pH);
327 double channelProbability = 0.55;
328 pDissociationChannel3->SetProbability(channelProbability);
329 probaRemaining -= channelProbability;
331 A1B1_DissociationDecay);
332 pH2ODef->AddDecayChannel(pH2Ovib, pDissociationChannel3);
335 auto pDissociationChannel1 =
337 pDissociationChannel1->SetProbability(probaRemaining);
338 pH2ODef->AddDecayChannel(pH2Ovib, pDissociationChannel1);
352 BuildDissociationChannels();
364 if (FindReactant(track))
377 if (FindReactant(track))
double B(double temperature)
double epsilon(double density, double temperature)
double A(double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4Molecule * GetMolecule(const G4Track &track)
G4DNAElectronHoleRecombination()
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void StartTracking(G4Track *)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual ~G4DNAElectronHoleRecombination()
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
static G4DNAMolecularMaterial * Instance()
const std::vector< G4double > * GetDensityTableFor(const G4Material *) const
Retrieve a table of volumetric mass densities (mass per unit volume) in the G4 unit system for chosen...
static G4Electron_aq * Definition()
static G4H2O * DefinitionIfExists()
static G4H2O * Definition()
G4KDTreeResultHandle FindNearestInRange(const T *point, int key, G4double)
static G4ITFinder * Instance()
G4double GetTemperature() const
G4MolecularConfiguration * NewConfiguration(const G4String &excitedStateLabel)
static G4MoleculeTable * Instance()
void ChangeConfigurationToLabel(const G4String &label)
virtual void Initialize(const G4Track &)
void SetTrackStatus(const G4TrackStatus aTrackStatus)
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4Material * GetMaterial() const
virtual void StartTracking(G4Track *)
void SetInstantiateProcessState(G4bool flag)
G4shared_ptr< G4ProcessState > fpState
virtual void ClearInteractionTimeLeft()
virtual void ClearNumberOfInteractionLengthLeft()
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
static G4VMoleculeCounter * Instance()
void ProposeTrackStatus(G4TrackStatus status)
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
virtual void StartTracking(G4Track *)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange