Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LindhardSorensenIonModel.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// GEANT4 Class header file
29//
30//
31// File name: G4LindhardSorensenIonModel
32//
33// Author: Alexander Bagulya & Vladimir Ivanchenko
34//
35// Creation date: 16.04.2018
36//
37//
38// Class Description:
39//
40// Implementation of ion ionisation energy loss and delta-electron
41// production by heavy charged particles according to
42// J. Lindhard & A.H. Sorensen, Phys. Rev. A 53 (1996) 2443-2455
43
44// -------------------------------------------------------------------
45//
46
47#ifndef G4LindhardSorensenIonModel_h
48#define G4LindhardSorensenIonModel_h 1
49
50#include <vector>
51
52#include "G4VEmModel.hh"
53#include "G4NistManager.hh"
54
55class G4EmCorrections;
58class G4BraggModel;
60class G4IonICRU73Data;
61
63{
64public:
65
66 explicit G4LindhardSorensenIonModel(const G4ParticleDefinition* p = nullptr,
67 const G4String& nam = "LindhardSorensen");
68
70
71 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
72
74 const G4MaterialCutsCouple* couple) override;
75
78 G4double kineticEnergy,
79 G4double cutEnergy,
80 G4double maxEnergy);
81
84 G4double kineticEnergy,
86 G4double cutEnergy,
87 G4double maxEnergy) override;
88
91 G4double kineticEnergy,
92 G4double cutEnergy,
93 G4double maxEnergy) override;
94
97 G4double kineticEnergy,
98 G4double cutEnergy) override;
99
101 const G4Material* mat,
102 G4double kineticEnergy) override;
103
105 const G4Material* mat,
106 G4double kineticEnergy) override;
107
109 const G4DynamicParticle* dp,
110 const G4double& length,
111 G4double& eloss) override;
112
113 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
115 const G4DynamicParticle*,
116 G4double tmin,
117 G4double maxEnergy) override;
118
119 // hide assignment operator
121 (const G4LindhardSorensenIonModel &right) = delete;
123
124protected:
125
127 G4double kinEnergy) override;
128
129 inline void SetChargeSquareRatio(G4double val);
130
131private:
132
133 void SetupParameters();
134
135 void InitialiseLS();
136
137 G4double ComputeDEDXPerVolumeLS(const G4Material*,
139 G4double kinEnergy, G4double cutEnergy);
140
141 inline void SetParticle(const G4ParticleDefinition* p);
142
143 // static const G4int MAXZION = 93;
144
145 static G4IonICRU73Data* fIonData;
146 static G4LindhardSorensenData* lsdata;
147
148 const G4ParticleDefinition* particle = nullptr;
149 G4ParticleDefinition* theElectron;
150 G4EmCorrections* corr;
151 G4ParticleChangeForLoss* fParticleChange = nullptr;
152 G4NistManager* nist;
153 G4BraggModel* fBraggModel;
154 G4BetheBlochModel* fBBModel;
155
156 G4int Zin = 1;
157 G4double mass = 0.0;
158 G4double tlimit = DBL_MAX;
159 G4double spin = 0.0;
160 G4double magMoment2 = 0.0;
161 G4double chargeSquare = 1.0;
162 G4double charge = 1.0;
163 G4double eRatio = 0.0;
164 G4double pRatio = 1.0;
165 G4double formfact = 0.0;
166 G4double twoln10;
167 G4double fElimit;
168 G4bool isFirst = false;
169};
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
173inline void
174G4LindhardSorensenIonModel::SetParticle(const G4ParticleDefinition* p)
175{
176 if(particle != p) {
177 particle = p;
178 SetupParameters();
179 }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
185{
186 chargeSquare = val;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190
191#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double GetChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LindhardSorensenIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LindhardSorensen")
void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, const G4double &length, G4double &eloss) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4LindhardSorensenIonModel(const G4LindhardSorensenIonModel &)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
#define DBL_MAX
Definition templates.hh:62