Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNABornIonisationModel2.hh
이 파일의 문서화 페이지로 가기
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: G4DNABornIonisationModel2.hh 81398 2014-05-28 09:53:02Z vnivanch $
27 //
28 
29 #ifndef G4DNABornIonisationModel2_h
30 #define G4DNABornIonisationModel2_h 1
31 
32 #include "G4VEmModel.hh"
34 #include "G4ProductionCutsTable.hh"
35 
37 #include "G4Electron.hh"
38 #include "G4Proton.hh"
40 
41 #include "G4LogLogInterpolation.hh"
42 
44 #include "G4VAtomDeexcitation.hh"
45 #include "G4NistManager.hh"
46 
47 
49 {
50 
51 public:
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*>*,
67  const G4MaterialCutsCouple*,
68  const G4DynamicParticle*,
69  G4double tmin,
70  G4double maxEnergy);
71 
73  G4int /*level*/,
74  const G4ParticleDefinition*,
75  G4double /*kineticEnergy*/);
76 
77  G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
78 
79  G4double TransferedEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random) ;
80 
81  inline void SelectFasterComputation(G4bool input);
82 
83  inline void SelectStationary(G4bool input);
84 
85  inline void SelectSPScaling(G4bool input);
86 
87 protected:
88 
90 
91 private:
92 
96 
97  // Water density table
98  const std::vector<G4double>* fpMolWaterDensity;
99 
100  // Deexcitation manager to produce fluo photons and e-
102 
105 
107 
110 
111  // Cross section
114 
115  // Final state
116 
118 
119  G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
120 
121  G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
122 
124 
126  G4double e12,
127  G4double e21,
128  G4double e22,
129  G4double x11,
130  G4double x12,
131  G4double x21,
132  G4double x22,
133  G4double t1,
134  G4double t2,
135  G4double t,
136  G4double e);
137 
138  typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
139 
141  TriDimensionMap fNrjTransfData[6]; // for cumulated dcs
142 
143  std::vector<G4double> fTdummyVec;
144 
145  typedef std::map<G4double, std::vector<G4double> > VecMap;
146 
148  VecMap fProbaShellMap[6]; // for cumulated dcs
149 
150  // Partial cross section
151 
153 
154  //
155 
158 
159 };
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
164 {
165  fasterCode = input;
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
171 {
172  statCode = input;
173 }
174 
175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176 
178 {
179  spScaling = input;
180 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183 
184 #endif
TTree * t1
Definition: plottest35.C:26
G4DNACrossSectionDataSet * fTableData
const char * p
Definition: xmltok.h:285
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
G4ParticleChangeForGamma * fParticleChangeForGamma
static const G4double emax
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
std::map< G4double, std::vector< G4double > > VecMap
std::vector< G4double > fTdummyVec
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4DNABornIonisationModel2 & operator=(const G4DNABornIonisationModel2 &right)
std::map< G4double, std::map< G4double, G4double > > TriDimensionMap
G4DNABornIonisationModel2(const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel")
double energy
Definition: plottest35.C:25
virtual G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double)
G4VAtomDeexcitation * fAtomDeexcitation
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
TTree * t2
Definition: plottest35.C:36
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
int G4int
Definition: G4Types.hh:78
const G4ParticleDefinition * fParticleDef
G4DNAWaterIonisationStructure waterStructure
const std::vector< G4double > * fpMolWaterDensity
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)