Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuIonisation.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: G4MuIonisation.cc 68035 2013-03-13 14:12:34Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4MuIonisation
34 //
35 // Author: Laszlo Urban
36 //
37 // Creation date: 30.09.1997
38 //
39 // Modifications:
40 //
41 // 08-04-98 remove 'tracking cut' of the ionizing particle (mma)
42 // 26-10-98 new stuff from R.Kokoulin + cleanup , L.Urban
43 // 10-02-00 modifications , new e.m. structure, L.Urban
44 // 23-03-01 R.Kokoulin's correction is commented out, L.Urban
45 // 29-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
46 // 10-08-01 new methods Store/Retrieve PhysicsTable (mma)
47 // 28-08-01 new function ComputeRestrictedMeandEdx() + 'cleanup' (mma)
48 // 17-09-01 migration of Materials to pure STL (mma)
49 // 26-09-01 completion of RetrievePhysicsTable (mma)
50 // 29-10-01 all static functions no more inlined (mma)
51 // 07-11-01 correction(Tmax+xsection computation) L.Urban
52 // 08-11-01 particleMass becomes a local variable (mma)
53 // 10-05-02 V.Ivanchenko update to new design
54 // 04-12-02 V.Ivanchenko the low energy limit for Kokoulin model to 10 GeV
55 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
56 // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
57 // 13-02-03 SubCutoff regime is assigned to a region (V.Ivanchenko)
58 // 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko)
59 // 03-06-03 Add SetIntegral method to choose fluctuation model (V.Ivanchenko)
60 // 03-06-03 Fix initialisation problem for STD ionisation (V.Ivanchenko)
61 // 04-08-03 Set integral=false to be default (V.Ivanchenko)
62 // 08-08-03 STD substitute standard (V.Ivanchenko)
63 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
64 // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.Ivanchenko)
65 // 27-05-04 Set integral to be a default regime (V.Ivanchenko)
66 // 17-08-04 Utilise mu+ tables for mu- (V.Ivanchenko)
67 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
68 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
69 // 12-08-05 SetStepLimits(0.2, 0.1*mm) (mma)
70 // 02-09-05 SetStepLimits(0.2, 1*mm) (V.Ivantchenko)
71 // 12-08-05 SetStepLimits(0.2, 0.1*mm) + integral off (V.Ivantchenko)
72 // 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
73 //
74 // -------------------------------------------------------------------
75 //
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
79 #include "G4MuIonisation.hh"
80 #include "G4PhysicalConstants.hh"
81 #include "G4SystemOfUnits.hh"
82 #include "G4Electron.hh"
83 #include "G4MuonPlus.hh"
84 #include "G4MuonMinus.hh"
85 #include "G4BraggModel.hh"
86 #include "G4BetheBlochModel.hh"
87 #include "G4MuBetheBlochModel.hh"
89 #include "G4IonFluctuations.hh"
90 #include "G4BohrFluctuations.hh"
91 #include "G4UnitsTable.hh"
92 #include "G4ICRU73QOModel.hh"
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
96 using namespace std;
97 
99  : G4VEnergyLossProcess(name),
100  theParticle(0),
101  theBaseParticle(0),
102  isInitialised(false)
103 {
104  mass = ratio = 0;
105  // SetStepFunction(0.2, 1*mm);
106  //SetIntegral(true);
107  //SetVerboseLevel(1);
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113 
115 {}
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118 
120 {
121  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV);
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125 
127  const G4Material*,
128  G4double cut)
129 {
130  G4double x = 0.5*cut/electron_mass_c2;
131  G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
132  return mass*(gam - 1.0);
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136 
138  const G4ParticleDefinition* bpart)
139 {
140  if(!isInitialised) {
141 
142  theParticle = part;
143  theBaseParticle = bpart;
144 
145  mass = theParticle->GetPDGMass();
146  G4double q = theParticle->GetPDGCharge();
147  G4double elow = 0.2*MeV;
148 
149  // Bragg peak model
150  if (!EmModel(1)) {
151  if(q > 0.0) { SetEmModel(new G4BraggModel(),1); }
152  else { SetEmModel(new G4ICRU73QOModel(),1); }
153  }
155  EmModel(1)->SetHighEnergyLimit(elow);
156  AddEmModel(1, EmModel(1), new G4IonFluctuations());
157 
158  // high energy fluctuation model
160 
161  // moderate energy model
162  if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
163  EmModel(2)->SetLowEnergyLimit(elow);
164  EmModel(2)->SetHighEnergyLimit(1.0*GeV);
165  AddEmModel(2, EmModel(2), FluctModel());
166 
167  // high energy model
168  if (!EmModel(3)) { SetEmModel(new G4MuBetheBlochModel(),3); }
169  EmModel(3)->SetLowEnergyLimit(1.0*GeV);
171  AddEmModel(3, EmModel(3), FluctModel());
172 
173  ratio = electron_mass_c2/mass;
174  isInitialised = true;
175  }
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
181 {}
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184 
185 
186 
187 
virtual ~G4MuIonisation()
const char * p
Definition: xmltok.h:285
const XML_Char * name
void SetFluctModel(G4VEmFluctuationModel *)
G4MuIonisation(const G4String &name="muIoni")
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4VEmFluctuationModel * FluctModel()
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut)
bool G4bool
Definition: G4Types.hh:79
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
float electron_mass_c2
Definition: hepunit.py:274
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4double GetPDGMass() const
void SetEmModel(G4VEmModel *, G4int index=1)
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)
static G4Electron * Electron()
Definition: G4Electron.cc:94
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
G4double GetPDGCharge() const
virtual void PrintInfo()
G4VEmModel * EmModel(G4int index=1) const