Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4MuPairProductionModel.hh 74544 2013-10-14 12:40:29Z gcosmo $
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 "G4ElementData.hh"
66 #include "G4Physics2DVector.hh"
67 #include <vector>
68 
69 class G4Element;
72 
74 {
75 public:
76 
78  const G4String& nam = "muPairProd");
79 
80  virtual ~G4MuPairProductionModel();
81 
82  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83 
84  virtual void InitialiseLocal(const G4ParticleDefinition*,
85  G4VEmModel* masterModel);
86 
88  const G4ParticleDefinition*,
89  G4double kineticEnergy,
90  G4double Z, G4double A,
91  G4double cutEnergy,
92  G4double maxEnergy);
93 
95  const G4ParticleDefinition*,
96  G4double kineticEnergy,
97  G4double cutEnergy);
98 
99  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
100  const G4MaterialCutsCouple*,
101  const G4DynamicParticle*,
102  G4double tmin,
103  G4double maxEnergy);
104 
105  virtual G4double MinPrimaryEnergy(const G4Material*,
106  const G4ParticleDefinition*,
107  G4double);
108 
109  inline void SetLowestKineticEnergy(G4double e);
110 
111  inline void SetParticle(const G4ParticleDefinition*);
112 
113 protected:
114 
116  G4double tmax);
117 
119  G4double Z,
120  G4double cut);
121 
123  G4double Z,
124  G4double pairEnergy);
125 
126  inline G4double MaxSecondaryEnergyForElement(G4double kineticEnergy,
127  G4double Z);
128 
129 private:
130 
131  void MakeSamplingTables();
132 
133  void DataCorrupted(G4int Z, G4double logTkin);
134 
135  inline G4double FindScaledEnergy(G4int Z, G4double rand, G4double logTkin,
136  G4double yymin, G4double yymax);
137 
138  // hide assignment operator
141 
142 protected:
143 
146 
154 
155  static const G4double xgi[8],wgi[8];
156 
157 private:
158 
159  G4ParticleDefinition* theElectron;
160  G4ParticleDefinition* thePositron;
161  G4ParticleChangeForLoss* fParticleChange;
162 
163  G4double minPairEnergy;
164  G4double lowestKinEnergy;
165 
166  G4int nzdat;
167 
168  // gamma energy bins
169  G4int nYBinPerDecade;
170  size_t nbiny;
171  size_t nbine;
172  G4double ymin;
173  G4double dy;
174  G4double emin;
175  G4double emax;
176 
177  static const G4int zdat[5];
178  static const G4double adat[5];
179 };
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182 
184 {
185  lowestKinEnergy = e;
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
190 inline
192 {
193  if(!particle) {
194  particle = p;
196  }
197 }
198 
199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200 
201 inline G4double
203  G4double ZZ)
204 {
205  G4int Z = G4lrint(ZZ);
206  if(Z != currentZ) {
207  currentZ = Z;
208  z13 = nist->GetZ13(Z);
209  z23 = z13*z13;
210  lnZ = nist->GetLOGZ(Z);
211  }
212  return kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 
217 inline G4double
218 G4MuPairProductionModel::FindScaledEnergy(G4int Z, G4double rand,
219  G4double logTkin,
220  G4double yymin, G4double yymax)
221 {
222  G4double res = yymin;
224  if(!pv) {
225  DataCorrupted(Z, logTkin);
226  } else {
227  G4double pmin = pv->Value(yymin, logTkin);
228  G4double pmax = pv->Value(yymax, logTkin);
229  G4double p0 = pv->Value(0.0, logTkin);
230  if(p0 <= 0.0) { DataCorrupted(Z, logTkin); }
231  else { res = pv->FindLinearX((pmin + rand*(pmax - pmin))/p0, logTkin); }
232  }
233  return res;
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
238 #endif
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
static const G4double xgi[8]
G4MuPairProductionModel(const G4ParticleDefinition *p=0, const G4String &nam="muPairProd")
const G4ParticleDefinition * particle
const char * p
Definition: xmltok.h:285
G4double GetZ13(G4double Z)
void SetParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
static const G4double wgi[8]
G4double GetLOGZ(G4int Z)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
G4Physics2DVector * GetElement2DData(G4int Z)
G4ElementData * fElementData
Definition: G4VEmModel.hh:397
G4double FindLinearX(G4double rand, G4double y, size_t &lastidy) const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)