Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30//
31// File name: G4MuIonisation
32//
33// Author: Laszlo Urban
34//
35// Creation date: 30.09.1997
36//
37// Modifications:
38//
39// 08-04-98 remove 'tracking cut' of the ionizing particle (mma)
40// 26-10-98 new stuff from R.Kokoulin + cleanup , L.Urban
41// 10-02-00 modifications , new e.m. structure, L.Urban
42// 23-03-01 R.Kokoulin's correction is commented out, L.Urban
43// 29-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
44// 10-08-01 new methods Store/Retrieve PhysicsTable (mma)
45// 28-08-01 new function ComputeRestrictedMeandEdx() + 'cleanup' (mma)
46// 17-09-01 migration of Materials to pure STL (mma)
47// 26-09-01 completion of RetrievePhysicsTable (mma)
48// 29-10-01 all static functions no more inlined (mma)
49// 07-11-01 correction(Tmax+xsection computation) L.Urban
50// 08-11-01 particleMass becomes a local variable (mma)
51// 10-05-02 V.Ivanchenko update to new design
52// 04-12-02 V.Ivanchenko the low energy limit for Kokoulin model to 10 GeV
53// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
54// 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
55// 13-02-03 SubCutoff regime is assigned to a region (V.Ivanchenko)
56// 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko)
57// 03-06-03 Add SetIntegral method to choose fluctuation model (V.Ivanchenko)
58// 03-06-03 Fix initialisation problem for STD ionisation (V.Ivanchenko)
59// 04-08-03 Set integral=false to be default (V.Ivanchenko)
60// 08-08-03 STD substitute standard (V.Ivanchenko)
61// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
62// 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I.)
63// 27-05-04 Set integral to be a default regime (V.Ivanchenko)
64// 17-08-04 Utilise mu+ tables for mu- (V.Ivanchenko)
65// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
66// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
67// 12-08-05 SetStepLimits(0.2, 0.1*mm) (mma)
68// 02-09-05 SetStepLimits(0.2, 1*mm) (V.Ivantchenko)
69// 12-08-05 SetStepLimits(0.2, 0.1*mm) + integral off (V.Ivantchenko)
70// 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
71//
72// -------------------------------------------------------------------
73//
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77#include "G4MuIonisation.hh"
79#include "G4SystemOfUnits.hh"
80#include "G4Electron.hh"
81#include "G4BraggModel.hh"
82#include "G4BetheBlochModel.hh"
84#include "G4EmStandUtil.hh"
85#include "G4ICRU73QOModel.hh"
86#include "G4EmParameters.hh"
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
92{
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
100{
101 return (p.GetPDGCharge() != 0.0);
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107 const G4Material*,
108 G4double cut)
109{
110 G4double x = 0.5*cut/CLHEP::electron_mass_c2;
111 G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
112 return mass*(gam - 1.0);
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
117void
119 const G4ParticleDefinition* bpart)
120{
121 if(!isInitialised) {
122
123 theParticle = part;
124 theBaseParticle = bpart;
125
126 mass = theParticle->GetPDGMass();
127 G4double q = theParticle->GetPDGCharge();
128
130 G4double elow = 0.2*CLHEP::MeV;
131 G4double emax = param->MaxKinEnergy();
132
133 // Bragg peak model
134 if (nullptr == EmModel(0)) {
135 if(q > 0.0) { SetEmModel(new G4BraggModel()); }
136 else { SetEmModel(new G4ICRU73QOModel()); }
137 }
139 EmModel(0)->SetHighEnergyLimit(elow);
140
141 // high energy fluctuation model
142 if (nullptr == FluctModel()) {
144 }
145 // low-energy fluctuation model
147 AddEmModel(1, EmModel(0), f);
148
149 // high energy model
150 if (nullptr == EmModel(1)) { SetEmModel(new G4MuBetheBlochModel()); }
151 EmModel(1)->SetLowEnergyLimit(elow);
152 EmModel(1)->SetHighEnergyLimit(emax);
153 AddEmModel(1, EmModel(1), FluctModel());
154
155 ratio = CLHEP::electron_mass_c2/mass;
156 isInitialised = true;
157 }
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
162void G4MuIonisation::ProcessDescription(std::ostream& out) const
163{
164 out << " Muon ionisation";
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fIonisation
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
void ProcessDescription(std::ostream &) const override
G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) override
void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
G4MuIonisation(const G4String &name="muIoni")
G4bool IsApplicable(const G4ParticleDefinition &p) override
G4double GetPDGCharge() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void ProcessDescription(std::ostream &outFile) const override
G4VEmModel * EmModel(std::size_t index=0) const
void SetEmModel(G4VEmModel *, G4int index=0)
G4VEmFluctuationModel * FluctModel() const
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410