source: trunk/source/processes/electromagnetic/lowenergy/include/G4DNARuddIonisationModel.hh @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 5.1 KB
Line 
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: G4DNARuddIonisationModel.hh,v 1.4 2010/01/07 18:10:19 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29
30#ifndef G4DNARuddIonisationModel_h
31#define G4DNARuddIonisationModel_h 1
32
33#include "G4VEmModel.hh"
34#include "G4ParticleChangeForGamma.hh"
35#include "G4ProductionCutsTable.hh"
36
37#include "G4DNAGenericIonsManager.hh"
38//#include "G4DNAGenericMoleculeManager.hh"
39#include "G4DNACrossSectionDataSet.hh"
40#include "G4Electron.hh"
41#include "G4Proton.hh"
42#include "G4LogLogInterpolation.hh"
43
44#include "G4WaterIonisationStructure.hh"
45
46class G4DNARuddIonisationModel : public G4VEmModel
47{
48
49public:
50
51  G4DNARuddIonisationModel(const G4ParticleDefinition* p = 0, 
52                           const G4String& nam = "DNARuddIonisationModel");
53
54  virtual ~G4DNARuddIonisationModel();
55
56  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
57
58  virtual G4double CrossSectionPerVolume(  const G4Material* material,
59                                           const G4ParticleDefinition* p,
60                                           G4double ekin,
61                                           G4double emin,
62                                           G4double emax);
63
64  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
65                                 const G4MaterialCutsCouple*,
66                                 const G4DynamicParticle*,
67                                 G4double tmin,
68                                 G4double maxEnergy);
69
70protected:
71
72  G4ParticleChangeForGamma* fParticleChangeForGamma;
73
74private:
75
76  std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
77  std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
78
79  G4double lowEnergyLimitForZ1; 
80  G4double lowEnergyLimitForZ2; 
81  G4double lowEnergyLimitOfModelForZ1; 
82  G4double lowEnergyLimitOfModelForZ2; 
83  G4double killBelowEnergyForZ1; 
84  G4double killBelowEnergyForZ2; 
85
86  G4bool isInitialised;
87  G4int verboseLevel;
88 
89  // Cross section
90
91  typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
92  MapFile tableFile;
93
94  typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
95  MapData tableData;
96 
97  // Final state
98 
99  G4WaterIonisationStructure waterStructure;
100
101  G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition, 
102                                          G4double incomingParticleEnergy, 
103                                          G4int shell);
104
105  void RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition, 
106                                         G4double incomingParticleEnergy, 
107                                         G4double outgoingParticleEnergy, 
108                                         G4double & cosTheta, 
109                                         G4double & phi);
110   
111  G4double  DifferentialCrossSection(G4ParticleDefinition* particleDefinition, 
112                                   G4double k, 
113                                   G4double energyTransfer, 
114                                   G4int shell);
115
116  G4double CorrectionFactor(G4ParticleDefinition* particleDefinition, G4double k);
117
118  G4double S_1s(G4double t, 
119                G4double energyTransferred, 
120                G4double slaterEffectiveChg, 
121                G4double shellNumber);
122
123  G4double S_2s(G4double t, 
124                G4double energyTransferred, 
125                G4double slaterEffectiveChg, 
126                G4double shellNumber);
127
128
129  G4double S_2p(G4double t, 
130                G4double energyTransferred, 
131                G4double slaterEffectiveChg, 
132                G4double shellNumber);
133
134  G4double R(G4double t, 
135             G4double energyTransferred, 
136             G4double slaterEffectiveChg, 
137             G4double shellNumber) ;
138
139  G4double slaterEffectiveCharge[3];
140  G4double sCoefficient[3];
141 
142  // Partial cross section
143 
144  G4double PartialCrossSection(const G4Track& track);
145                       
146  G4double Sum(G4double energy, const G4String& particle);
147
148  G4int RandomSelect(G4double energy,const G4String& particle );
149 
150  //
151   
152  G4DNARuddIonisationModel & operator=(const  G4DNARuddIonisationModel &right);
153  G4DNARuddIonisationModel(const  G4DNARuddIonisationModel&);
154
155};
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158
159#endif
Note: See TracBrowser for help on using the repository browser.