Geant4  10.03.p02
 All Classes 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 97392 2016-06-02 10:10:32Z 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 
77  explicit G4MuPairProductionModel(const G4ParticleDefinition* p = nullptr,
78  const G4String& nam = "muPairProd");
79 
80  virtual ~G4MuPairProductionModel();
81 
82  virtual void Initialise(const G4ParticleDefinition*,
83  const G4DataVector&) override;
84 
85  virtual void InitialiseLocal(const G4ParticleDefinition*,
86  G4VEmModel* masterModel) override;
87 
89  const G4ParticleDefinition*,
90  G4double kineticEnergy,
92  G4double cutEnergy,
93  G4double maxEnergy) override;
94 
96  const G4ParticleDefinition*,
97  G4double kineticEnergy,
98  G4double cutEnergy) override;
99 
100  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
101  const G4MaterialCutsCouple*,
102  const G4DynamicParticle*,
103  G4double tmin,
104  G4double maxEnergy) override;
105 
106  virtual G4double MinPrimaryEnergy(const G4Material*,
107  const G4ParticleDefinition*,
108  G4double) override;
109 
110  inline void SetLowestKineticEnergy(G4double e);
111 
112  inline void SetParticle(const G4ParticleDefinition*);
113 
114 protected:
115 
117  G4double tmax);
118 
120  G4double Z,
121  G4double cut);
122 
123  virtual G4double
125  G4double pairEnergy);
126 
127  inline G4double MaxSecondaryEnergyForElement(G4double kineticEnergy,
128  G4double Z);
129 
130 private:
131 
132  void MakeSamplingTables();
133 
134  void DataCorrupted(G4int Z, G4double logTkin);
135 
136  inline G4double FindScaledEnergy(G4int Z, G4double rand, G4double logTkin,
137  G4double yymin, G4double yymax);
138 
139  // hide assignment operator
140  G4MuPairProductionModel & operator=(const G4MuPairProductionModel &right) = delete;
142 
143 protected:
144 
147 
155 
156  static const G4double xgi[8],wgi[8];
157 
158 private:
159 
160  G4ParticleDefinition* theElectron;
161  G4ParticleDefinition* thePositron;
162  G4ParticleChangeForLoss* fParticleChange;
163 
164  G4double minPairEnergy;
165  G4double lowestKinEnergy;
166 
167  G4int nzdat;
168 
169  // gamma energy bins
170  G4int nYBinPerDecade;
171  size_t nbiny;
172  size_t nbine;
173  G4double ymin;
174  G4double dy;
175  G4double emin;
176  G4double emax;
177 };
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
182 {
183  lowestKinEnergy = e;
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
187 
188 inline
190 {
191  if(!particle) {
192  particle = p;
194  }
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198 
199 inline G4double
201  G4double ZZ)
202 {
203  G4int Z = G4lrint(ZZ);
204  if(Z != currentZ) {
205  currentZ = Z;
206  z13 = nist->GetZ13(Z);
207  z23 = z13*z13;
208  lnZ = nist->GetLOGZ(Z);
209  }
210  return kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
211 }
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214 
215 inline G4double
216 G4MuPairProductionModel::FindScaledEnergy(G4int Z, G4double rand,
217  G4double logTkin,
218  G4double yymin, G4double yymax)
219 {
220  G4double res = yymin;
222  if(!pv) {
223  DataCorrupted(Z, logTkin);
224  } else {
225  G4double pmin = pv->Value(yymin, logTkin);
226  G4double pmax = pv->Value(yymax, logTkin);
227  G4double p0 = pv->Value(0.0, logTkin);
228  if(p0 <= 0.0) { DataCorrupted(Z, logTkin); }
229  else { res = pv->FindLinearX((pmin + rand*(pmax - pmin))/p0, logTkin); }
230  }
231  return res;
232 }
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235 
236 #endif
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static const G4double xgi[8]
const G4ParticleDefinition * particle
G4MuPairProductionModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="muPairProd")
const char * p
Definition: xmltok.h:285
void SetParticle(const G4ParticleDefinition *)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
int G4int
Definition: G4Types.hh:78
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
G4double GetZ13(G4double Z) const
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
double A(double temperature)
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
G4double GetLOGZ(G4int Z) const
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]
double G4double
Definition: G4Types.hh:76
G4Physics2DVector * GetElement2DData(G4int Z)
G4ElementData * fElementData
Definition: G4VEmModel.hh:423
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double FindLinearX(G4double rand, G4double y, size_t &lastidy) const
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override