Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GammaConversion.cc
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: G4GammaConversion.cc 84598 2014-10-17 07:39:15Z gcosmo $
27 //
28 //
29 //------------------ G4GammaConversion physics process -------------------------
30 // by Michel Maire, 24 May 1996
31 //
32 // 11-06-96 Added SelectRandomAtom() method, M.Maire
33 // 21-06-96 SetCuts implementation, M.Maire
34 // 24-06-96 simplification in ComputeCrossSectionPerAtom, M.Maire
35 // 24-06-96 in DoIt : change the particleType stuff, M.Maire
36 // 25-06-96 modification in the generation of the teta angle, M.Maire
37 // 16-09-96 minors optimisations in DoIt. Thanks to P.Urban
38 // dynamical array PartialSumSigma
39 // 13-12-96 fast sampling of epsil below 2 MeV, L.Urban
40 // 14-01-97 crossection table + meanfreepath table.
41 // PartialSumSigma removed, M.Maire
42 // 14-01-97 in DoIt the positron is always created, even with Ekine=0,
43 // for further annihilation, M.Maire
44 // 14-03-97 new Physics scheme for geant4alpha, M.Maire
45 // 28-03-97 protection in BuildPhysicsTable, M.Maire
46 // 19-06-97 correction in ComputeCrossSectionPerAtom, L.Urban
47 // 04-06-98 in DoIt, secondary production condition:
48 // range>std::min(threshold,safety)
49 // 13-08-98 new methods SetBining() PrintInfo()
50 // 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
51 // 11-07-01 PostStepDoIt - sampling epsil: power(rndm,0.333333)
52 // 13-07-01 DoIt: suppression of production cut for the (e-,e+) (mma)
53 // 06-08-01 new methods Store/Retrieve PhysicsTable (mma)
54 // 06-08-01 BuildThePhysicsTable() called from constructor (mma)
55 // 17-09-01 migration of Materials to pure STL (mma)
56 // 20-09-01 DoIt: fminimalEnergy = 1*eV (mma)
57 // 01-10-01 come back to BuildPhysicsTable(const G4ParticleDefinition&)
58 // 11-01-02 ComputeCrossSection: correction of extrapolation below EnergyLimit
59 // 21-03-02 DoIt: correction of the e+e- angular distribution (bug 363) mma
60 // 08-11-04 Remove of Store/Retrieve tables (V.Ivantchenko)
61 // 19-04-05 Migrate to model interface and inherit
62 // from G4VEmProcess (V.Ivanchenko)
63 // 04-05-05, Make class to be default (V.Ivanchenko)
64 // 09-08-06, add SetModel(G4VEmModel*) (mma)
65 // 12-09-06, move SetModel(G4VEmModel*) in G4VEmProcess (mma)
66 // -----------------------------------------------------------------------------
67 
68 #include "G4GammaConversion.hh"
69 #include "G4PhysicalConstants.hh"
70 #include "G4SystemOfUnits.hh"
71 #include "G4BetheHeitlerModel.hh"
73 #include "G4Electron.hh"
74 #include "G4EmParameters.hh"
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77 
78 using namespace std;
79 
81  G4ProcessType type):G4VEmProcess (processName, type),
82  isInitialised(false)
83 {
87  SetBuildTableFlag(true);
89  SetLambdaBinning(220);
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
95 {}
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100 {
101  return (&p == G4Gamma::Gamma());
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107 {
108  if(!isInitialised) {
109  isInitialised = true;
111  G4double emin = std::max(param->MinKinEnergy(), 2*electron_mass_c2);
112  G4double emax = param->MaxKinEnergy();
113  G4double energyLimit = std::min(emax, 80*GeV);
114 
115  SetMinKinEnergy(emin);
116 
117  if(!EmModel(1)) { SetEmModel(new G4BetheHeitlerModel(), 1); }
118  EmModel(1)->SetLowEnergyLimit(emin);
119  EmModel(1)->SetHighEnergyLimit(energyLimit);
120  AddEmModel(1, EmModel(1));
121 
122  if(emax > energyLimit) {
123  if(!EmModel(2)) { SetEmModel(new G4PairProductionRelModel(), 2); }
124  EmModel(2)->SetLowEnergyLimit(energyLimit);
125  EmModel(2)->SetHighEnergyLimit(emax);
126  AddEmModel(2, EmModel(2));
127  }
128  }
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132 
134  const G4Material*)
135 {
136  return 2*electron_mass_c2;
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142 {}
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4GammaConversion(const G4String &processName="conv", G4ProcessType type=fElectromagnetic)
G4double MaxKinEnergy() const
void SetBuildTableFlag(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
const char * p
Definition: xmltok.h:285
virtual G4bool IsApplicable(const G4ParticleDefinition &) final
virtual void InitialiseProcess(const G4ParticleDefinition *) override
void SetStartFromNullFlag(G4bool val)
static constexpr double electron_mass_c2
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void SetLambdaBinning(G4int nbins)
void SetEmModel(G4VEmModel *, G4int index=1)
bool G4bool
Definition: G4Types.hh:79
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4double MinKinEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual void PrintInfo() override
static const G4double emax
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *) override
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetSecondaryParticle(const G4ParticleDefinition *p)
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetMinKinEnergy(G4double e)
virtual ~G4GammaConversion()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
G4ProcessType