Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 #ifndef SOURCE_PROCESSES_ELECTROMAGNETIC_DNA_PROCESSES_INCLUDE_G4DNAELECTRONHOLERECOMBINATION_HH_
34 #define SOURCE_PROCESSES_ELECTROMAGNETIC_DNA_PROCESSES_INCLUDE_G4DNAELECTRONHOLERECOMBINATION_HH_
35 
37 
39 {
40 public:
43  void Create();
44 
45  void StartTracking(G4Track*);
46 
48 
49  virtual void BuildPhysicsTable(const G4ParticleDefinition&){}
50 
52  // DoIt /////////////////
54  virtual G4VParticleChange* AtRestDoIt(const G4Track& /*track*/,
55  const G4Step& /*stepData*/);
56  // A virtual base class function that has to be overridden
57  // by any subclass. The DoIt method actually performs the
58  // physics process and determines either momentum change
59  // of the production of secondaries etc.
60  // arguments
61  // const G4Track& track:
62  // reference to the current G4Track information
63  // const G4Step& stepData:
64  // reference to the current G4Step information
65 
66  virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
67 
68 
70  // GPIL //////////////
72 
73  virtual G4double GetMeanFreePath(const G4Track& aTrack,
74  G4double previousStepSize,
76  // Calculates from the macroscopic cross section a mean
77  // free path, the value is returned in units of distance.
78 
79  virtual G4double GetMeanLifeTime(const G4Track& aTrack,
81  // Calculates the mean life-time (i.e. for decays) of the
82  // particle at rest due to the occurence of the given process,
83  // or converts the probability of interaction (i.e. for
84  // annihilation) into the life-time of the particle for the
85  // occurence of the given process.
86 
87 private:
88 
89  struct ReactionProfile
90  {
91  G4Track* fElectron;
92  G4double fDistance;
93  G4double fProbability;
94  };
95 
96  struct State : public G4ProcessStateBase<G4DNAElectronHoleRecombination>
97  {
98  std::vector<ReactionProfile> fReactants; // distanceSqr;
99  G4double fSampleProba;
100  };
101 
102  G4bool FindReactant(const G4Track& track);
103  void MakeReaction(const G4Track& track);
104 
105  const std::vector<double>* fpMoleculeDensity;
106  G4ParticleChange fParticleChange;
107  G4bool fIsInitialized;
108  std::map<int, std::pair<double, double> > fOnsagerRadiusPerMaterial;
109 
110 };
111 
112 #endif /* SOURCE_PROCESSES_ELECTROMAGNETIC_DNA_PROCESSES_INCLUDE_G4DNAELECTRONHOLERECOMBINATION_HH_ */
G4double condition(const G4ErrorSymMatrix &m)
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
#define State(X)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
bool G4bool
Definition: G4Types.hh:79
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
Definition: G4Step.hh:76
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)
double G4double
Definition: G4Types.hh:76
G4ForceCondition