Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuPairProductionModel.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//
29// GEANT4 Class header file
30//
31//
32// File name: G4MuPairProductionModel
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 18.05.2002
37//
38// Modifications:
39//
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 27-01-03 Make models region aware (V.Ivanchenko)
42// 13-02-03 Add name (V.Ivanchenko)
43// 10-02-04 Update parameterisation using R.Kokoulin model (V.Ivanchenko)
44// 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
45// 13-02-06 Add ComputeCrossSectionPerAtom (mma)
46// 12-05-06 Add parameter to SelectRandomAtom (A.Bogdanov)
47// 11-10-07 Add ignoreCut flag (V.Ivanchenko)
48// 28-02-08 Reorganized protected methods and members (V.Ivanchenko)
49//
50// Class Description:
51//
52// Implementation of e+e- pair production by muons
53// A.G. Bogdanov et al., IEEE Trans. Nuc. Sci., Vol.53, No.2, 2006
54//
55// -------------------------------------------------------------------
56//
57
58#ifndef G4MuPairProductionModel_h
59#define G4MuPairProductionModel_h 1
60
61#include "G4VEmModel.hh"
62#include "G4NistManager.hh"
63#include "G4ElementData.hh"
64#include "G4Physics2DVector.hh"
65#include <vector>
66
67class G4Element;
70
72{
73public:
74
75 explicit G4MuPairProductionModel(const G4ParticleDefinition* p = nullptr,
76 const G4String& nam = "muPairProd");
77
79
80 virtual void Initialise(const G4ParticleDefinition*,
81 const G4DataVector&) override;
82
83 virtual void InitialiseLocal(const G4ParticleDefinition*,
84 G4VEmModel* masterModel) override;
85
88 G4double kineticEnergy,
90 G4double cutEnergy,
91 G4double maxEnergy) override;
92
95 G4double kineticEnergy,
96 G4double cutEnergy) override;
97
98 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
100 const G4DynamicParticle*,
101 G4double tmin,
102 G4double maxEnergy) override;
103
104 virtual G4double MinPrimaryEnergy(const G4Material*,
106 G4double) override;
107
108 inline void SetLowestKineticEnergy(G4double e);
109
110 inline void SetParticle(const G4ParticleDefinition*);
111
112 // hide assignment operator and copy constructor
114 (const G4MuPairProductionModel &right) = delete;
116
117protected:
118
120 G4double tmax);
121
123 G4double Z,
124 G4double cut);
125
126 virtual G4double
128 G4double pairEnergy);
129
131 G4double Z);
132
133private:
134
135 void MakeSamplingTables();
136
137 void StoreTables() const;
138
139 G4bool RetrieveTables();
140
141 void DataCorrupted(G4int Z, G4double logTkin) const;
142
143 inline G4double FindScaledEnergy(G4int Z, G4double rand, G4double logTkin,
144 G4double yymin, G4double yymax);
145
146protected:
147
150
158
162
165
166 // gamma energy bins
168 size_t nbiny;
169 size_t nbine;
174
176};
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
181{
182 lowestKinEnergy = e;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186
187inline
189{
190 if(!particle) {
191 particle = p;
193 }
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197
198inline G4double
200 G4double ZZ)
201{
202 G4int Z = G4lrint(ZZ);
203 if(Z != currentZ) {
204 currentZ = Z;
205 z13 = nist->GetZ13(Z);
206 z23 = z13*z13;
207 lnZ = nist->GetLOGZ(Z);
208 }
209 return kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
210}
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213
214inline G4double
215G4MuPairProductionModel::FindScaledEnergy(G4int Z, G4double rand,
216 G4double logTkin,
217 G4double yymin, G4double yymax)
218{
219 G4double res = yymin;
221 if(!pv) {
222 DataCorrupted(Z, logTkin);
223 } else {
224 G4double pmin = pv->Value(yymin, logTkin);
225 G4double pmax = pv->Value(yymax, logTkin);
226 G4double p0 = pv->Value(0.0, logTkin);
227 if(p0 <= 0.0) { DataCorrupted(Z, logTkin); }
228 else { res = pv->FindLinearX((pmin + rand*(pmax - pmin))/p0, logTkin); }
229 }
230 return res;
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234
235#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4Physics2DVector * GetElement2DData(G4int Z)
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
void SetParticle(const G4ParticleDefinition *)
G4ParticleDefinition * thePositron
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4MuPairProductionModel(const G4MuPairProductionModel &)=delete
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
const G4ParticleDefinition * particle
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
G4ParticleChangeForLoss * fParticleChange
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
G4ParticleDefinition * theElectron
G4double GetZ13(G4double Z) const
G4double GetLOGZ(G4int Z) const
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
G4double FindLinearX(G4double rand, G4double y, std::size_t &lastidy) const
G4ElementData * fElementData
Definition: G4VEmModel.hh:438
int G4lrint(double ad)
Definition: templates.hh:134