Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel.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// $Id$
27//
28
29#ifndef G4DNABornIonisationModel_h
30#define G4DNABornIonisationModel_h 1
31
32#include "G4VEmModel.hh"
35
37#include "G4Electron.hh"
38#include "G4Proton.hh"
40
42
45#include "G4NistManager.hh"
46
47
49{
50
51public:
52
54 const G4String& nam = "DNABornIonisationModel");
55
57
58 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector()));
59
60 virtual G4double CrossSectionPerVolume( const G4Material* material,
61 const G4ParticleDefinition* p,
62 G4double ekin,
63 G4double emin,
64 G4double emax);
65
66 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
68 const G4DynamicParticle*,
69 G4double tmin,
70 G4double maxEnergy);
71
72 double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
73
74protected:
75
77
78private:
79
80 // Water density table
81 const std::vector<G4double>* fpMolWaterDensity;
82
83 //deexcitation manager to produce fluo photns and e-
84 G4VAtomDeexcitation* fAtomDeexcitation;
85
86 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
87 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
88
89 G4bool isInitialised;
90 G4int verboseLevel;
91
92 // Cross section
93
94 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
95 MapFile tableFile;
96
97 typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
98 MapData tableData;
99
100 // Final state
101
102 G4DNAWaterIonisationStructure waterStructure;
103
104 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
105
106 void RandomizeEjectedElectronDirection(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4double
107 outgoingParticleEnergy, G4double & cosTheta, G4double & phi );
108
109 G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
110
111 G4double QuadInterpolator( G4double e11,
112 G4double e12,
113 G4double e21,
114 G4double e22,
115 G4double x11,
116 G4double x12,
117 G4double x21,
118 G4double x22,
119 G4double t1,
120 G4double t2,
121 G4double t,
122 G4double e);
123
124 typedef std::map<double, std::map<double, double> > TriDimensionMap;
125 TriDimensionMap eDiffCrossSectionData[6];
126 TriDimensionMap pDiffCrossSectionData[6];
127 std::vector<double> eTdummyVec;
128 std::vector<double> pTdummyVec;
129
130 typedef std::map<double, std::vector<double> > VecMap;
131 VecMap eVecm;
132 VecMap pVecm;
133
134 // Partial cross section
135
136 G4int RandomSelect(G4double energy,const G4String& particle );
137
138 //
139
140 G4DNABornIonisationModel & operator=(const G4DNABornIonisationModel &right);
142
143};
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146
147#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4ParticleChangeForGamma * fParticleChangeForGamma
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)