Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4MuPairProductionModel
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 18.05.2002
38//
39// Modifications:
40//
41// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42// 27-01-03 Make models region aware (V.Ivanchenko)
43// 13-02-03 Add name (V.Ivanchenko)
44// 10-02-04 Update parameterisation using R.Kokoulin model (V.Ivanchenko)
45// 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
46// 13-02-06 Add ComputeCrossSectionPerAtom (mma)
47// 12-05-06 Add parameter to SelectRandomAtom (A.Bogdanov)
48// 11-10-07 Add ignoreCut flag (V.Ivanchenko)
49// 28-02-08 Reorganized protected methods and members (V.Ivanchenko)
50
51//
52// Class Description:
53//
54// Implementation of e+e- pair production by muons
55//
56
57// -------------------------------------------------------------------
58//
59
60#ifndef G4MuPairProductionModel_h
61#define G4MuPairProductionModel_h 1
62
63#include "G4VEmModel.hh"
64#include "G4NistManager.hh"
65#include <vector>
66
67class G4Element;
70
72{
73public:
74
76 const G4String& nam = "muPairProd");
77
79
80 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
81
84 G4double kineticEnergy,
85 G4double Z, G4double A,
86 G4double cutEnergy,
87 G4double maxEnergy);
88
91 G4double kineticEnergy,
92 G4double cutEnergy);
93
94 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
96 const G4DynamicParticle*,
97 G4double tmin,
98 G4double maxEnergy);
99
101 const G4MaterialCutsCouple*);
102
103 inline void SetLowestKineticEnergy(G4double e);
104
105 inline void SetParticle(const G4ParticleDefinition*);
106
107protected:
108
110 G4double tmax);
111
113 G4double Z,
114 G4double cut);
115
117 G4double Z,
118 G4double pairEnergy);
119
121 G4double kineticEnergy);
122
123 inline void SetCurrentElement(G4double Z);
124
125private:
126
127 const G4Element* SelectRandomAtom(G4double kinEnergy,
128 G4double dt,
129 G4int it,
130 const G4MaterialCutsCouple* couple,
131 G4double tmin);
132
133 void MakeSamplingTables();
134
135 inline G4double InterpolatedIntegralCrossSection(
136 G4double dt, G4double dz, G4int iz,
137 G4int it, G4int iy, G4double z);
138
139 // hide assignment operator
140 G4MuPairProductionModel & operator=(const G4MuPairProductionModel &right);
142
143protected:
144
147
155
156 static G4double xgi[8],wgi[8];
157
158private:
159
160 G4ParticleDefinition* theElectron;
161 G4ParticleDefinition* thePositron;
162 G4ParticleChangeForLoss* fParticleChange;
163
164 G4double minPairEnergy;
165 G4double lowestKinEnergy;
166
167 // tables for sampling
168 G4int nzdat;
169 G4int ntdat;
170 G4int nbiny;
171 size_t nmaxElements;
172 static G4double zdat[5], adat[5], tdat[8];
173 G4double ya[1001], proba[5][8][1001];
174
175 G4double ymin;
176 G4double ymax;
177 G4double dy;
178
179 G4bool samplingTablesAreFilled;
180 std::vector<G4double> partialSum;
181};
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
186{
187 lowestKinEnergy = e;
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191
192inline
194{
195 if(!particle) {
196 particle = p;
198 }
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
204{
205 if(Z != currentZ) {
206 currentZ = Z;
207 G4int iz = G4int(Z);
208 z13 = nist->GetZ13(iz);
209 z23 = z13*z13;
210 lnZ = nist->GetLOGZ(iz);
211 }
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
216inline G4double G4MuPairProductionModel::InterpolatedIntegralCrossSection(
217 G4double dt, G4double dz,
218 G4int iz, G4int it, G4int iy, G4double z)
219{
220 G4double fac = 1./(zdat[iz] *(zdat[iz] +1.));
221 G4double fac1 = 1./(zdat[iz-1]*(zdat[iz-1]+1.));
222 G4double f0 = fac1*proba[iz-1][it-1][iy] +
223 (fac*proba[iz][it-1][iy]-fac1*proba[iz-1][it-1][iy])*dz;
224 G4double f1 = fac1*proba[iz-1][it ][iy] +
225 (fac*proba[iz][it ][iy]-fac1*proba[iz-1][it ][iy])*dz;
226 return (f0 + (f1-f0)*dt)*z*(z+1.);
227}
228
229//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230
231#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
void SetParticle(const G4ParticleDefinition *)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
const G4ParticleDefinition * particle
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4double GetZ13(G4double Z)
G4double GetLOGZ(G4int Z)