Geant4-11
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#include "G4Threading.hh"
55
56class G4EmCorrections;
59class G4BraggIonModel;
61class G4IonICRU73Data;
62
64{
65public:
66
67 explicit G4LindhardSorensenIonModel(const G4ParticleDefinition* p = nullptr,
68 const G4String& nam = "LindhardSorensen");
69
71
72 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
73
75 const G4MaterialCutsCouple* couple) override;
76
79 G4double kineticEnergy,
80 G4double cutEnergy,
81 G4double maxEnergy);
82
85 G4double kineticEnergy,
87 G4double cutEnergy,
88 G4double maxEnergy) override;
89
92 G4double kineticEnergy,
93 G4double cutEnergy,
94 G4double maxEnergy) override;
95
98 G4double kineticEnergy,
99 G4double cutEnergy) override;
100
102 const G4Material* mat,
103 G4double kineticEnergy) override;
104
106 const G4Material* mat,
107 G4double kineticEnergy) override;
108
110 const G4DynamicParticle* dp,
111 const G4double& length,
112 G4double& eloss) override;
113
114 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
116 const G4DynamicParticle*,
117 G4double tmin,
118 G4double maxEnergy) override;
119
120 // hide assignment operator
122 (const G4LindhardSorensenIonModel &right) = delete;
124
125protected:
126
128 G4double kinEnergy) override;
129
130 inline G4double GetChargeSquareRatio() const;
131
132 inline void SetChargeSquareRatio(G4double val);
133
134private:
135
136 void SetupParameters();
137
139
142 G4double kinEnergy, G4double cutEnergy);
143
144 inline void SetParticle(const G4ParticleDefinition* p);
145
146 static const G4int MAXZION = 93;
147
150 static std::vector<G4float>* fact[MAXZION];
151
159
172
173#ifdef G4MULTITHREADED
174 static G4Mutex theLSMutex;
175#endif
176
177};
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181inline void
183{
184 if(particle != p) {
185 particle = p;
187 }
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
191
193{
194 return chargeSquare;
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198
200{
201 chargeSquare = val;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205
206#endif
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
static std::vector< G4float > * fact[MAXZION]
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4ParticleChangeForLoss * fParticleChange
G4double ComputeDEDXPerVolumeLS(const G4Material *, const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy)
static G4LindhardSorensenData * lsdata
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
void SetParticle(const G4ParticleDefinition *p)
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
const G4ParticleDefinition * particle
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