Geant4-11
G4DNADingfelderChargeIncreaseModel.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
28#ifndef G4DNADingfelderChargeIncreaseModel_h
29#define G4DNADingfelderChargeIncreaseModel_h 1
30
31#include "G4VEmModel.hh"
34
35#include "G4Electron.hh"
36#include "G4Proton.hh"
38#include "G4NistManager.hh"
39
41{
42
43public:
44
46 const G4String& nam = "DNADingfelderChargeIncreaseModel");
47
49
50 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
51
53 const G4ParticleDefinition* p,
54 G4double ekin,
55 G4double emin,
57
58 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
60 const G4DynamicParticle*,
61 G4double tmin,
62 G4double maxEnergy);
63
64 inline void SelectStationary(G4bool input);
65
66protected:
67
69
70private:
71
72 // Water density table
73 const std::vector<G4double>* fpMolWaterDensity;
74
75 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
76 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
77
80
81 // Partial cross section
82
84
86
88
89 G4int numberOfPartialCrossSections[2]; // 2 is the particle type index
90
91 G4double f0[2][2];
92 G4double a0[2][2];
93 G4double a1[2][2];
94 G4double b0[2][2];
95 G4double b1[2][2];
96 G4double c0[2][2];
97 G4double d0[2][2];
98 G4double x0[2][2];
99 G4double x1[2][2];
100
101 // Final state
102
103 G4int NumberOfFinalStates(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
104
105 G4ParticleDefinition* OutgoingParticleDefinition(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
106
107 G4double WaterBindingEnergyConstant(G4ParticleDefinition * aParticleDefinition, G4int finalStateIndex);
108
110
111 G4double IncomingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
112
113 //
114
117
118};
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
122{
123 statCode = input;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
128#endif
129
static const G4double emax
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4int NumberOfFinalStates(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4double IncomingParticleBindingEnergyConstant(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4int RandomSelect(G4double energy, const G4ParticleDefinition *particle)
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
G4DNADingfelderChargeIncreaseModel(const G4DNADingfelderChargeIncreaseModel &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4ParticleDefinition * OutgoingParticleDefinition(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4double Sum(G4double energy, const G4ParticleDefinition *particle)
G4DNADingfelderChargeIncreaseModel(const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeIncreaseModel")
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)
G4double WaterBindingEnergyConstant(G4ParticleDefinition *aParticleDefinition, G4int finalStateIndex)
G4DNADingfelderChargeIncreaseModel & operator=(const G4DNADingfelderChargeIncreaseModel &right)
G4double OutgoingParticleBindingEnergyConstant(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4double energy(const ThreeVector &p, const G4double m)
string material
Definition: eplot.py:19