Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAElectronHoleRecombination.hh
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.hh
28 *
29 * Created on: Jun 17, 2015
30 * Author: mkaramit
31 */
32
33#pragma once
35
37{
38public:
41 void Create();
42
44
46
48
49 ////////////////////////////
50 // DoIt /////////////////
51 ///////////////////////////
52 virtual G4VParticleChange* AtRestDoIt(const G4Track& /*track*/,
53 const G4Step& /*stepData*/);
54 // A virtual base class function that has to be overridden
55 // by any subclass. The DoIt method actually performs the
56 // physics process and determines either momentum change
57 // of the production of secondaries etc.
58 // arguments
59 // const G4Track& track:
60 // reference to the current G4Track information
61 // const G4Step& stepData:
62 // reference to the current G4Step information
63
64 virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
65
66
67 //////////////////////////
68 // GPIL //////////////
69 /////////////////////////
70
71 virtual G4double GetMeanFreePath(const G4Track& aTrack,
72 G4double previousStepSize,
74 // Calculates from the macroscopic cross section a mean
75 // free path, the value is returned in units of distance.
76
77 virtual G4double GetMeanLifeTime(const G4Track& aTrack,
79 // Calculates the mean life-time (i.e. for decays) of the
80 // particle at rest due to the occurrence of the given process,
81 // or converts the probability of interaction (i.e. for
82 // annihilation) into the life-time of the particle for the
83 // occurence of the given process.
84
85private:
86 struct ReactantInfo
87 {
88 G4Track* fElectron;
89 G4double fDistance;
90 G4double fProbability;
91 };
92
93 struct State : public G4ProcessStateBase<G4DNAElectronHoleRecombination>
94 {
95 std::vector<ReactantInfo> fReactants; // distanceSqr;
96 G4double fSampleProba;
97 };
98
99 G4bool FindReactant(const G4Track& track);
100 void MakeReaction(const G4Track& track);
101 void BuildDissociationChannels();
102
103 const std::vector<double>* fpMoleculeDensity;
104 G4ParticleChange fParticleChange;
105 G4bool fIsInitialized;
106 std::map<int, std::pair<double, double> > fOnsagerRadiusPerMaterial;
107};
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
#define State(X)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
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 &)
Definition: G4Step.hh:62
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.