ATLAS Offline Software
Loading...
Searching...
No Matches
G4mplAtlasIonisationModel.cxx
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: G4mplAtlasIonisationModel.cxx 729653 2016-03-14 15:55:47Z jchapman $
27// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4mplAtlasIonisationModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 06.09.2005
39//
40// Modifications:
41//
42//
43// -------------------------------------------------------------------
44//
45
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50// class header
51#include "G4mplAtlasIonisationModel.hh"
52// Geant4 headers
53#include "Randomize.hh"
54#include "G4LossTableManager.hh"
55#include "G4ParticleChangeForLoss.hh"
56#include "G4Version.hh"
57// CLHEP headers
58#include "CLHEP/Units/SystemOfUnits.h"
59#include "CLHEP/Units/PhysicalConstants.h"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63using namespace std;
64
65G4mplAtlasIonisationModel::G4mplAtlasIonisationModel(G4double mCharge,
66 const G4String& nam)
67 : G4VEmModel(nam),G4VEmFluctuationModel(nam),
68 monopole(0),
69 fParticleChange(0),
70 mass(100000.),
71 magCharge(mCharge),
72 twoln10(2.0*log(10.0)),
73 beta2low(0.0001),
74 beta2lim(0.01),
75 bg2lim(beta2lim*(1.0 + beta2lim))
76{
77 std::cout <<"!!! G4mplAtlasIonisationModel constructor"<<std::endl;
78
79
80 nmpl = G4int(abs(magCharge)/68.0);
81 if(nmpl > 6) nmpl = 6;
82 else if(nmpl < 1) nmpl = 1;
83 G4double x = 45.0*CLHEP::GeV*G4double(nmpl)/CLHEP::cm;
84 factlow = x*x;
85 chargeSquare = magCharge*magCharge;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
90G4mplAtlasIonisationModel::~G4mplAtlasIonisationModel()
91{
92 std::cout <<"!!! G4mplAtlasIonisationModel destructor"<<std::endl;
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
97void G4mplAtlasIonisationModel::Initialise(const G4ParticleDefinition* p,
98 const G4DataVector&)
99{
100
101 std::cout <<"!!! G4mplAtlasIonisationModel::Initialise"<<std::endl;
102
103 monopole = p;
104 mass = monopole->GetPDGMass();
105
106 if(pParticleChange)
107 fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
108 else
109 fParticleChange = new G4ParticleChangeForLoss();
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
114G4double G4mplAtlasIonisationModel::ComputeDEDXPerVolume(const G4Material* material,
115 const G4ParticleDefinition*,
116 G4double kineticEnergy,
117 G4double)
118{
119 G4double tau = kineticEnergy/mass;
120 G4double gam = tau + 1.0;
121 G4double bg2 = tau * (tau+2.0);
122 G4double beta2 = bg2/(gam*gam);
123
124 G4double dedx0 = factlow*abs(beta2);
125
126 if(beta2 > beta2low) {
127
128 G4double b2 = beta2;
129 if(beta2 < beta2lim) {
130 beta2= beta2lim;
131 bg2 = bg2lim;
132 }
133
134 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
135 G4double cden = material->GetIonisation()->GetCdensity();
136 G4double mden = material->GetIonisation()->GetMdensity();
137 G4double aden = material->GetIonisation()->GetAdensity();
138 G4double x0den = material->GetIonisation()->GetX0density();
139 G4double x1den = material->GetIonisation()->GetX1density();
140
141 // std::cout <<"SB: G4mplAtlasIonisationModel"<<std::endl;
142
143 G4double eDensity = material->GetElectronDensity();
144
145 G4double dedx = 2.0*log(2.0*CLHEP::electron_mass_c2*bg2/eexc) - 1.0;
146
147 G4double k = 0.406;
148 if(nmpl > 1) k = 0.346;
149 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
150
151 dedx += k - B[nmpl];
152
153 // density correction
154 G4double x = log(bg2)/twoln10;
155 if ( x >= x0den ) {
156 dedx -= twoln10*x - cden ;
157 if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
158 }
159
160 // now compute the total ionization loss
161
162 if (dedx < 0.0) dedx = 0.0 ;
163
164 dedx *= CLHEP::twopi_mc2_rcl2*chargeSquare*eDensity;
165
166 // extrapolate between two formula
167 if(beta2 < beta2lim) {
168 x = log(dedx0) + log(dedx/dedx0)*log(b2/beta2low)/log(beta2lim/beta2low);
169 dedx = exp(x);
170 }
171 dedx0 = dedx;
172 }
173
174 return dedx0;
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
179G4double G4mplAtlasIonisationModel::SampleFluctuations(
180#if G4VERSION_NUMBER < 1100
181 const G4MaterialCutsCouple* material,
182 const G4DynamicParticle* dp,
183 G4double tmax,
184 G4double length,
185 G4double meanLoss)
186#else
187 const G4MaterialCutsCouple* material,
188 const G4DynamicParticle* dp,
189 const G4double tcut,
190 const G4double tmax,
191 const G4double length,
192 const G4double meanLoss)
193
194#endif
195{
196#if G4VERSION_NUMBER < 1100
197 G4double siga = Dispersion(material->GetMaterial(),dp,tmax,length);
198#else
199 G4double siga = Dispersion(material->GetMaterial(),dp,tcut,tmax,length);
200#endif
201 G4double loss = meanLoss;
202 siga = sqrt(siga);
203 G4double twomeanLoss = meanLoss + meanLoss;
204
205 if(twomeanLoss < siga) {
206 G4double x;
207 do {
208 loss = twomeanLoss*G4UniformRand();
209 x = (loss - meanLoss)/siga;
210 } while (1.0 - 0.5*x*x < G4UniformRand());
211 } else {
212 do {
213 loss = G4RandGauss::shoot(meanLoss,siga);
214 } while (0.0 > loss || loss > twomeanLoss);
215 }
216 //G4cout << "G4mplAtlasIonisationModel::SampleFluctuations: loss= " << loss
217 //<< " siga= " << siga << G4endl;
218 return loss;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double length(const pvec &v)
#define x
constexpr int pow(int base, int exp) noexcept
STL namespace.