Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LowPAIH2O.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// G4LowPAIxs.hh -- header file
27//
28// GEANT 4 class header file --- Copyright CERN 1995
29// CERB Geneva Switzerland
30//
31// for information related to this code, please, contact
32// CERN, CN Division, ASD Group
33//
34// Preparation of ionizing collision cross section according to Photo Absorption
35// Ionization (PAI) model for simulation of ionization energy losses in very thin
36// layers of water for low energy protons and electrons.
37//
38// Author: Vladimir.Grichine@cern.ch
39//
40// History:
41//
42
43// 2.12.24, V. Grichine: 1st version
44
45#ifndef G4LOWPAIH2O_HH
46#define G4LOWPAIH2O_HH
47
48#include "G4ios.hh"
49#include "globals.hh"
50#include "Randomize.hh"
51#include "G4PhysicsLogVector.hh"
52#include "G4DataVector.hh"
53#include "G4PhysicsTable.hh"
54#include "G4VEmModel.hh"
56
58class G4Material;
59class G4SandiaTable;
61class G4PhysicsTable;
64
66 public G4VEmModel , public G4VEmFluctuationModel
67{
68public:
69 // Constructors
70
71 explicit G4LowPAIH2O( const G4ParticleDefinition* p = nullptr,
72 const G4String& nam = "lowpaih2o");
73
74 ~G4LowPAIH2O() override;
75
76 G4LowPAIH2O & operator=(const G4LowPAIH2O &right) = delete;
77 G4LowPAIH2O(const G4LowPAIH2O&) = delete;
78
79 // methods
80
81 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
82
84 G4VEmModel* masterModel) override;
85
88 G4double kineticEnergy,
89 G4double cutEnergy,
90 G4double maxEnergy) override;
91
94 G4double kineticEnergy, G4double Z,
95 G4double A,
96 G4double cutEnergy,
97 G4double maxEnergy);
98
101 G4double kineticEnergy,
102 G4double cutEnergy,
103 G4double maxEnergy);
104
105 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
107 const G4DynamicParticle*,
108 G4double tmin,
109 G4double maxEnergy) override;
110
112 const G4DynamicParticle*,
113 const G4double tcut,
114 const G4double tmax,
115 const G4double length,
116 const G4double meanLoss) override;
119 const G4double kinEnergy,
120 const G4double cutEnergy,
121 const G4double& length,
122 G4double& eloss) override;
123
125 const G4double, const G4double, const G4double) override {return 0.;};
126
128 G4double kineticEnergy,
129 const G4Material*,
130 G4double cutEnergy = 0.0,
131 G4double maxEnergy = DBL_MAX);
132 void Initialize();
133 void InitRuthELF();
134
136 void BuildPrEnergyTable();
137 void BuildElEnergyTable();
141
143 G4double GetPrMFP( G4double Tkin);
145 G4double GetElMFP( G4double Tkin);
149
150 void SetBe2( G4double be2 ){ fBe2 = be2; };
151 G4double GetBe2(){ return fBe2; };
152 void SetOmega( G4double ww ){ fOmega = ww; };
153 G4double GetOmega(){ return fOmega;};
154 void SetBias( G4double bb ){ fBias = bb; };
155 G4double GetBias(){ return fBias;};
156
159
160 inline G4double GetSumELF( G4double energy );
161 inline G4double GetSumRuth( G4double energy );
162
163private:
164
165 // Local class members
166 G4int fTotBin{0};
167 G4int fBinTr{0};
168 G4int fBias{0};
169 G4double fCof{0.0};
170 G4double fBeta{0.0};
171 G4double fBe2{0.0};
172 G4double fTkin{0.0};
173
174 G4double fBmin{0.0};
175 G4double fBmax{0.0};
176 G4double fWmin{0.0};
177 G4double fWmax{0.0};
178 G4double fOmega{0.0};
179 G4double fElectronDensity{0.0};
180 G4double fNat{0.0};
181 G4double fNel{0.0};
182 G4double fMass{0.0};
183
184 static const G4int theBin;
185 static const G4double theEsum[1007], theELFsum[1007], theRuthSum[1007];
186
187 G4DataVector fEsum;
188 G4DataVector fELFsum, fRuthSum;
189 G4DataVector fPrWmaxVector;
190
191 G4Material* fMat{nullptr};
192
193 G4PhysicsLogVector* fBetaVector{nullptr};
194 G4PhysicsTable* fPrEnergyTable{nullptr};
195 G4PhysicsTable* fElEnergyTable{nullptr};
196 G4PhysicsLogVector* fTransferVector{nullptr};
197
198 G4ParticleDefinition* theElectron{nullptr};
199 G4ParticleDefinition* theProton{nullptr};
200 G4ParticleChangeForLoss* fParticleChange{nullptr};
201};
202
203///////////////////////////////////////////////////////////////////
204//////////////// Inline methods //////////////////////////////////
205////////////////////////////////////////////////////////////////////
206
207
208/////////////////// fast (STL) get ELF /////////////////
209
211{
212 G4double ee = energy; // /CLHEP::eV;
213 G4double elf(0.), y1(0.), y2(0.), x1(0.), x2(0.), aa(0.);
214
215 std::size_t nlow = std::lower_bound( fEsum.begin(), fEsum.end(), ee ) - fEsum.begin();
216
217 x1 = fEsum[nlow-1];
218 x2 = fEsum[nlow];
219 y1 = fELFsum[nlow-1];
220 y2 = fELFsum[nlow];
221 aa = (y2-y1)/(x2-x1);
222 elf = y1 + aa*(ee-x1);
223 return elf;
224}
225
226/////////////////////// fast (STL) get Rutherford /////////////////////
227
229{
230 G4double ee = energy; // /CLHEP::eV;
231 G4double ruth(0.), y1(0.), y2(0.), x1(0.), x2(0.), aa(0.);
232
233 std::size_t nlow = std::lower_bound( fEsum.begin(), fEsum.end(), ee ) - fEsum.begin();
234
235 x1 = fEsum[nlow-1];
236 x2 = fEsum[nlow];
237 y1 = fRuthSum[nlow-1];
238 y2 = fRuthSum[nlow];
239 aa = (y2-y1)/(x2-x1);
240 ruth = y1 + aa*(ee-x1);
241
242 return ruth;
243}
244
245///////////////////////////////////////////////
246
248 G4double Tkin,
249 const G4Material*,
250 G4double, // cutEnergy, // = 0.0,
251 G4double ) //maxEnergy ) // = DBL_MAX)
252{
253 G4double dndx(0.), mfp(DBL_MAX);
254
255 if ( pd == theProton ) dndx = GetPrdNdx(Tkin);
256 else if( pd == theElectron ) dndx = GetEldNdx(Tkin);
257 else return DBL_MAX;
258
259 if( dndx > 0.) mfp = 1./dndx;
260 else mfp = DBL_MAX;
261
262 mfp *= fBias;
263 return mfp;
264}
265
266#endif
267
268///////////////// end of G4LowPAIH2O header file //////////////////
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
void Initialize()
G4double GetPrMFP(G4double Tkin)
void BuildPrEnergyTable()
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double GetElectronTmax(G4double Tkin)
G4double GetSumELF(G4double energy)
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4double GetBias()
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void BuildPhysicsTable(const G4ParticleDefinition *pd)
G4double PrPAId2Ndxdw(G4double omega)
G4double CrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4double GetOmega()
G4double GetProtonTmax(G4double Tkin)
G4double CorrectElTransfer(G4double Tkin)
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
G4double GetElMFP(G4double Tkin)
void SetBias(G4double bb)
G4double GetPrdNdx(G4double Tkin)
void SetBe2(G4double be2)
G4double GetElTransfer(G4double Tkin)
void InitRuthELF()
void SetOmega(G4double ww)
void BuildElEnergyTable()
~G4LowPAIH2O() override
G4double ElPAId2Ndxdw(G4double omega)
G4double GetBe2()
G4LowPAIH2O(const G4ParticleDefinition *p=nullptr, const G4String &nam="lowpaih2o")
G4double GetSumRuth(G4double energy)
G4double CorrectPrTransfer(G4double Tkin)
G4double GetPrTransfer(G4double Tkin)
void CorrectionsAlongStep(const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss) override
G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double, const G4double, const G4double) override
G4LowPAIH2O(const G4LowPAIH2O &)=delete
G4double GetEldNdx(G4double Tkin)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4LowPAIH2O & operator=(const G4LowPAIH2O &right)=delete
G4VEmFluctuationModel(const G4String &nam)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
#define DBL_MAX
Definition templates.hh:62