Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAElectronHoleRecombination.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/*
27 * G4DNAElectronHoleRecombination.cc
28 *
29 * Created on: Jun 17, 2015
30 * Author: mkaramit
31 */
32
34#include <G4MoleculeFinder.hh>
36#include "G4Electron_aq.hh"
37#include "G4MoleculeTable.hh"
39#include "G4H2.hh"
40#include "G4H2O.hh"
42#include "G4MoleculeTable.hh"
45#include "G4SystemOfUnits.hh"
46#include "G4VMoleculeCounter.hh"
47#include "G4Exp.hh"
48
49static double onsager_constant = e_squared / (4. * pi * epsilon0 * k_Boltzmann);
50
51//------------------------------------------------------------------------------
52
53// Parameterisation of dielectric constant vs temperature and density
54
55double Y(double density)
56{
57 return 1. / (1. + 0.0012 / (density * density));
58}
59
60double A(double temperature)
61{
62 double temp_inverse = 1 / temperature;
63 return 0.7017
64 + 642.0 * temp_inverse
65 - 1.167e5 * temp_inverse * temp_inverse
66 + 9.190e6 * temp_inverse * temp_inverse * temp_inverse;
67}
68
69double B(double temperature)
70{
71 double temp_inverse = 1 / temperature;
72 return -2.71
73 + 275.4 * temp_inverse
74 + 0.3245e5 * temp_inverse * temp_inverse;
75}
76
77double S(double temp)
78{
79 double temp_inverse = 1 / temp;
80
81 return 1.667
82 - 11.41 * temp_inverse
83 - 35260.0 * temp_inverse * temp_inverse;
84}
85
86double C(double temp)
87{
88 return A(temp) - B(temp) - 3;
89}
90
91double D(double temp)
92{
93 return B(temp) + 3;
94}
95
96double epsilon(double density, double temperature)
97{
98 return 1 + G4Exp(std::log(10.) *
99 (Y(density) *
100 (C(temperature) + (S(temperature) - 1) * std::log(density) / std::log(10.))
101 + D(temperature) + std::log(density) / std::log(10.)));
102}
103
104//------------------------------------------------------------------------------
105
107 : G4VITRestDiscreteProcess("G4DNAElectronHoleRecombination",
109{
110 Create();
111}
112
114
116{
117 pParticleChange = &fParticleChange;
118 enableAtRestDoIt = true;
119 enableAlongStepDoIt = false;
120 enablePostStepDoIt = true;
121
123
125 // ie G4DNAElectronHoleRecombination uses a state class
126 // inheriting from G4ProcessState
127
128 fIsInitialized = false;
129 fProposesTimeStep = true;
130 fpMoleculeDensity = nullptr;
131
132 verboseLevel = 0;
133}
134
135//______________________________________________________________________________
136
139 const G4Step& /*stepData*/)
140{
141 fParticleChange.Initialize(track);
144 MakeReaction(track);
145 return &fParticleChange;
146}
147
148//______________________________________________________________________________
149
151{
153 G4VITProcess::fpState.reset(new State());
155}
156
157//______________________________________________________________________________
158
159void G4DNAElectronHoleRecombination::MakeReaction(const G4Track& track)
160{
161 fParticleChange.Initialize(track);
162 auto pState = fpState->GetState<State>();
163 double random = pState->fSampleProba;
164 std::vector<ReactantInfo>& reactants = pState->fReactants;
165
166 G4Track* pSelectedReactant = nullptr;
167
168 for (const auto& reactantInfo : reactants)
169 {
170 if (reactantInfo.fElectron->GetTrackStatus() != fAlive)
171 {
172 continue;
173 }
174 if (reactantInfo.fProbability > random)
175 {
176 pSelectedReactant = reactantInfo.fElectron;
177 }
178 break;
179 }
180
181 if (pSelectedReactant)
182 {
184 {
186 RemoveAMoleculeAtTime(GetMolecule(track)->GetMolecularConfiguration(),
187 track.GetGlobalTime(),
188 &(track.GetPosition()));
189 }
190 GetMolecule(track)->ChangeConfigurationToLabel("H2Ovib");
191
193 {
195 AddAMoleculeAtTime(GetMolecule(track)->GetMolecularConfiguration(),
196 track.GetGlobalTime(),
197 &(track.GetPosition()));
198 }
199
200 // fParticleChange.ProposeTrackStatus(fStopAndKill);
201 fParticleChange.ProposeTrackStatus(fStopButAlive);
202
203 pSelectedReactant->SetTrackStatus(fStopAndKill);
204 // G4TrackList::Pop(pSelectedReactant);
205 // G4ITTrackHolder::Instance()->PushToKill(pSelectedReactant);
206 }
207 else
208 {
209 fParticleChange.ProposeTrackStatus(fStopButAlive);
210 }
211}
212
213//______________________________________________________________________________
214
215G4bool G4DNAElectronHoleRecombination::FindReactant(const G4Track& track)
216{
217 if (GetMolecule(track)->GetCharge() <= 0)
218 {
219 // cannot react with electrons
220 return false;
221 }
222
223 const auto pDensityTable =
225
226 double temperature = track.GetMaterial()->GetTemperature();
227 double density = (*pDensityTable)[track.GetMaterial()->GetIndex()] / (g / (1e-2 * m * 1e-2 * m * 1e-2 * m));
228 double eps = epsilon(density, temperature);
229
230 double onsagerRadius = onsager_constant * 1. / (temperature * eps);
231
233
236 e_aq.GetMoleculeID(),
237 10. * onsagerRadius);
238
239 if (results == 0 || results->GetSize() == 0)
240 {
241 return false;
242 }
243
244 results->Sort();
245
246 auto pState = fpState->GetState<State>();
247 std::vector<ReactantInfo>& reactants = pState->fReactants;
248 pState->fSampleProba = G4UniformRand();
249
250 reactants.resize(results->GetSize());
251
252 for (size_t i = 0; !results->End(); results->Next(), ++i)
253 {
254 reactants[i].fElectron = results->GetItem<G4IT>()->GetTrack();
255 reactants[i].fDistance = std::sqrt(results->GetDistanceSqr());
256
257 if (reactants[i].fDistance != 0)
258 {
259 reactants[i].fProbability = 1. - G4Exp(-onsagerRadius / reactants[i].fDistance);
260 }
261 else
262 {
263 reactants[i].fProbability = 1.;
264 }
265 }
266
267 return reactants.empty() ? false : reactants[0].fProbability > pState->fSampleProba;
268}
269
270//______________________________________________________________________________
271
272void G4DNAElectronHoleRecombination::BuildDissociationChannels()
273{
274 auto pMoleculeTable = G4MoleculeTable::Instance();
275
276 auto pH2ODef = pMoleculeTable->GetMoleculeDefinition("H2O", false);
277
278 if (pH2ODef == nullptr) // if this condition does not hold => process cannot be applied
279 {
280 return;
281 }
282
283 const auto pH2Ovib = G4H2O::Definition()->NewConfiguration("H2Ovib");
284
285 assert(pH2Ovib != nullptr);
286
287 const auto pH2 = pMoleculeTable->GetConfiguration("H2", false);
288 const auto pOH = pMoleculeTable->GetConfiguration("OH", false);
289 const auto pH = pMoleculeTable->GetConfiguration("H", false);
290
291 double probaRemaining = 1.;
292
293 if (pOH || pH2)
294 {
295 auto pDissocationChannel2 =
296 new G4MolecularDissociationChannel("H2Ovib_DissociativeDecay1");
297 if (pH2)
298 {
299 pDissocationChannel2->AddProduct(pH2);
300 }
301 if (pOH)
302 {
303 pDissocationChannel2->AddProduct(pOH);
304 pDissocationChannel2->AddProduct(pOH);
305 }
306
307 double channelProbability = 0.15;
308 pDissocationChannel2->SetProbability(channelProbability);
309 probaRemaining -= channelProbability;
310 pDissocationChannel2->SetDisplacementType(G4DNAWaterDissociationDisplacer::
311 B1A1_DissociationDecay);
312 pH2ODef->AddDecayChannel(pH2Ovib, pDissocationChannel2);
313 }
314
315 if (pOH || pH)
316 {
317 auto pDissociationChannel3 =
318 new G4MolecularDissociationChannel("H2Ovib_DissociativeDecay2");
319 if (pOH)
320 {
321 pDissociationChannel3->AddProduct(pOH);
322 }
323 if (pH)
324 {
325 pDissociationChannel3->AddProduct(pH);
326 }
327 double channelProbability = 0.55;
328 pDissociationChannel3->SetProbability(channelProbability);
329 probaRemaining -= channelProbability;
330 pDissociationChannel3->SetDisplacementType(G4DNAWaterDissociationDisplacer::
331 A1B1_DissociationDecay);
332 pH2ODef->AddDecayChannel(pH2Ovib, pDissociationChannel3);
333 }
334
335 auto pDissociationChannel1 =
336 new G4MolecularDissociationChannel("H2Ovib_NonDissociative");
337 pDissociationChannel1->SetProbability(probaRemaining);
338 pH2ODef->AddDecayChannel(pH2Ovib, pDissociationChannel1);
339}
340
341G4bool
343IsApplicable(const G4ParticleDefinition& particle)
344{
345 if (&particle != G4H2O::DefinitionIfExists())
346 {
347 return false;
348 }
349
351 {
352 BuildDissociationChannels();
353 }
354
355 return true;
356}
357
358//______________________________________________________________________________
359
361 G4double,
363{
364 if (FindReactant(track))
365 {
366 return 0.;
367 }
368
369 return DBL_MAX;
370}
371
372//______________________________________________________________________________
373
376{
377 if (FindReactant(track))
378 {
379 return 0.;
380 }
381
382 return DBL_MAX;
383}
384
386 const G4Step& step)
387{
388 return AtRestDoIt(track, step);
389}
double B(double temperature)
double epsilon(double density, double temperature)
double S(double temp)
double C(double temp)
double Y(double density)
double A(double temperature)
double D(double temp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4ForceCondition
#define State(X)
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
@ fElectromagnetic
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
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()
Definition: G4H2O.cc:80
static G4H2O * Definition()
Definition: G4H2O.cc:42
G4KDTreeResultHandle FindNearestInRange(const T *point, int key, G4double)
static G4ITFinder * Instance()
Definition: G4IT.hh:88
G4double GetTemperature() const
Definition: G4Material.hh:180
size_t GetIndex() const
Definition: G4Material.hh:258
G4MolecularConfiguration * NewConfiguration(const G4String &excitedStateLabel)
static G4MoleculeTable * Instance()
void ChangeConfigurationToLabel(const G4String &label)
Definition: G4Molecule.cc:546
virtual void Initialize(const G4Track &)
Definition: G4Step.hh:62
void SetTrackStatus(const G4TrackStatus aTrackStatus)
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4Material * GetMaterial() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:85
void SetInstantiateProcessState(G4bool flag)
G4shared_ptr< G4ProcessState > fpState
virtual void ClearInteractionTimeLeft()
virtual void ClearNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
static G4VMoleculeCounter * Instance()
void ProposeTrackStatus(G4TrackStatus status)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
G4int verboseLevel
Definition: G4VProcess.hh:356
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4double pi
G4bool IsMasterThread()
Definition: G4Threading.cc:124
#define DBL_MAX
Definition: templates.hh:62