Geant4  10.02.p01
G4mplIonisationWithDeltaModel.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: G4mplIonisationWithDeltaModel.cc 91869 2015-08-07 15:21:02Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4mplIonisationWithDeltaModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 06.09.2005
38 //
39 // Modifications:
40 // 12.08.2007 Changing low energy approximation and extrapolation.
41 // Small bug fixing and refactoring (M. Vladymyrov)
42 // 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko)
43 //
44 //
45 // -------------------------------------------------------------------
46 // References
47 // [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles,
48 // S.P. Ahlen, Rev. Mod. Phys 52(1980), p121
49 // [2] K.A. Milton arXiv:hep-ex/0602040
50 // [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347
51 
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 
57 #include "Randomize.hh"
58 #include "G4PhysicalConstants.hh"
59 #include "G4SystemOfUnits.hh"
61 #include "G4Electron.hh"
62 #include "G4DynamicParticle.hh"
63 #include "G4ProductionCutsTable.hh"
64 #include "G4MaterialCutsCouple.hh"
65 #include "G4Log.hh"
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68 
69 using namespace std;
70 
71 std::vector<G4double>* G4mplIonisationWithDeltaModel::dedx0 = 0;
72 
74  const G4String& nam)
76  magCharge(mCharge),
77  twoln10(log(100.0)),
78  betalow(0.01),
79  betalim(0.1),
80  beta2lim(betalim*betalim),
81  bg2lim(beta2lim*(1.0 + beta2lim))
82 {
83  nmpl = G4lrint(std::fabs(magCharge) * 2 * fine_structure_const);
84  if(nmpl > 6) { nmpl = 6; }
85  else if(nmpl < 1) { nmpl = 1; }
86  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
88  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
89  fParticleChange = 0;
91  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
92  << magCharge/eplus << G4endl;
93  monopole = 0;
94  mass = 0.0;
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100 {
101  if(IsMaster()) { delete dedx0; }
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107 {
108  monopole = p;
109  mass = monopole->GetPDGMass();
110  G4double emin =
111  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
112  G4double emax =
113  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
114  SetLowEnergyLimit(emin);
115  SetHighEnergyLimit(emax);
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119 
120 void
122  const G4DataVector&)
123 {
124  if(!monopole) { SetParticle(p); }
126  if(IsMaster()) {
127  if(!dedx0) { dedx0 = new std::vector<G4double>; }
128  G4ProductionCutsTable* theCoupleTable=
130  G4int numOfCouples = theCoupleTable->GetTableSize();
131  G4int n = dedx0->size();
132  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
133 
134  // initialise vector
135  for(G4int i=0; i<numOfCouples; ++i) {
136 
137  const G4Material* material =
138  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
139  G4double eDensity = material->GetElectronDensity();
140  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
141  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
142  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
143  }
144  }
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
149 G4double
151  const G4ParticleDefinition* p,
152  G4double kineticEnergy,
153  G4double maxEnergy)
154 {
155  if(!monopole) { SetParticle(p); }
156  G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
157  G4double cutEnergy = std::min(tmax, maxEnergy);
158  cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
159  G4double tau = kineticEnergy / mass;
160  G4double gam = tau + 1.0;
161  G4double bg2 = tau * (tau + 2.0);
162  G4double beta2 = bg2 / (gam * gam);
163  G4double beta = sqrt(beta2);
164 
165  // low-energy asymptotic formula
166  //G4double dedx = dedxlim*beta*material->GetDensity();
167  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
168 
169  // above asymptotic
170  if(beta > betalow) {
171 
172  // high energy
173  if(beta >= betalim) {
174  dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
175 
176  } else {
177 
178  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
179  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
180  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
181 
182  // extrapolation between two formula
183  G4double kapa2 = beta - betalow;
184  G4double kapa1 = betalim - beta;
185  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
186  }
187  }
188  return dedx;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
193 G4double
195  G4double bg2,
196  G4double cutEnergy)
197 {
198  G4double eDensity = material->GetElectronDensity();
199  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
200 
201  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
202  G4double dedx =
203  0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
204 
205  // Kazama et al. cross-section correction
206  G4double k = 0.406;
207  if(nmpl > 1) { k = 0.346; }
208 
209  // Bloch correction
210  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
211 
212  dedx += 0.5 * k - B[nmpl];
213 
214  // density effect correction
215  G4double x = G4Log(bg2)/twoln10;
216  dedx -= material->GetIonisation()->DensityCorrection(x);
217 
218  // now compute the total ionization loss
219  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
220 
221  if (dedx < 0.0) { dedx = 0; }
222  return dedx;
223 }
224 
225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
226 
227 G4double
229  const G4ParticleDefinition* p,
230  G4double kineticEnergy,
231  G4double cut,
232  G4double maxKinEnergy)
233 {
234  if(!monopole) { SetParticle(p); }
235  G4double cross = 0.0;
236  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
237  G4double maxEnergy = std::min(tmax,maxKinEnergy);
238  G4double cutEnergy = std::max(LowEnergyLimit(), cut);
239  if(cutEnergy < maxEnergy) {
240  cross = (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl;
241  }
242  return cross;
243 }
244 
245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
246 
247 G4double
249  const G4ParticleDefinition* p,
250  G4double kineticEnergy,
251  G4double Z, G4double,
252  G4double cutEnergy,
253  G4double maxEnergy)
254 {
255  G4double cross =
256  Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
257  return cross;
258 }
259 
260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
261 
262 void
263 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
264  const G4MaterialCutsCouple*,
265  const G4DynamicParticle* dp,
266  G4double minKinEnergy,
267  G4double maxEnergy)
268 {
269  G4double kineticEnergy = dp->GetKineticEnergy();
270  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
271 
272  G4double maxKinEnergy = std::min(maxEnergy,tmax);
273  if(minKinEnergy >= maxKinEnergy) { return; }
274 
275  //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
276  // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
277  // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
278 
279  G4double totEnergy = kineticEnergy + mass;
280  G4double etot2 = totEnergy*totEnergy;
281  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
282 
283  // sampling without nuclear size effect
284  G4double q = G4UniformRand();
285  G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
286  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
287 
288  // delta-electron is produced
289  G4double totMomentum = totEnergy*sqrt(beta2);
290  G4double deltaMomentum =
291  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
292  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
293  (deltaMomentum * totMomentum);
294  if(cost > 1.0) { cost = 1.0; }
295 
296  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
297 
298  G4double phi = twopi * G4UniformRand() ;
299 
300  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
301  G4ThreeVector direction = dp->GetMomentumDirection();
302  deltaDirection.rotateUz(direction);
303 
304  // create G4DynamicParticle object for delta ray
305  G4DynamicParticle* delta =
306  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
307 
308  vdp->push_back(delta);
309 
310  // Change kinematics of primary particle
311  kineticEnergy -= deltaKinEnergy;
312  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
313  finalP = finalP.unit();
314 
317 }
318 
319 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
320 
322  const G4MaterialCutsCouple* couple,
323  const G4DynamicParticle* dp,
324  G4double tmax,
325  G4double length,
326  G4double meanLoss)
327 {
328  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
329  G4double loss = meanLoss;
330  siga = sqrt(siga);
331  G4double twomeanLoss = meanLoss + meanLoss;
332 
333  if(twomeanLoss < siga) {
334  G4double x;
335  do {
336  loss = twomeanLoss*G4UniformRand();
337  x = (loss - meanLoss)/siga;
338  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
339  } while (1.0 - 0.5*x*x < G4UniformRand());
340  } else {
341  do {
342  loss = G4RandGauss::shoot(meanLoss,siga);
343  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
344  } while (0.0 > loss || loss > twomeanLoss);
345  }
346  return loss;
347 }
348 
349 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350 
351 G4double
353  const G4DynamicParticle* dp,
354  G4double tmax,
355  G4double length)
356 {
357  G4double siga = 0.0;
358  G4double tau = dp->GetKineticEnergy()/mass;
359  if(tau > 0.0) {
360  G4double electronDensity = material->GetElectronDensity();
361  G4double gam = tau + 1.0;
362  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
363  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
364  * electronDensity * chargeSquare;
365  }
366  return siga;
367 }
368 
369 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
370 
371 G4double
373  G4double kinEnergy)
374 {
375  G4double tau = kinEnergy/mass;
376  return 2.0*electron_mass_c2*tau*(tau + 2.);
377 }
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
ThreeVector shoot(const G4int Ap, const G4int Af)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
static const double cm2
Definition: G4SIunits.hh:119
G4mplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)
void SetParticle(const G4ParticleDefinition *p)
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4ParticleDefinition * GetDefinition() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double B(double temperature)
int G4int
Definition: G4Types.hh:78
static const G4double bg2lim
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4double ComputeDEDXAhlen(const G4Material *material, G4double bg2, G4double cut)
static std::vector< G4double > * dedx0
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:217
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const double twopi
Definition: G4SIunits.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static const G4double beta2lim
static const double GeV
Definition: G4SIunits.hh:214
const G4int n
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static const G4double emax
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
static G4ProductionCutsTable * GetProductionCutsTable()
G4double DensityCorrection(G4double x)
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static const double pi
Definition: G4SIunits.hh:74
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const double g
Definition: G4SIunits.hh:180
const G4double x[NPOINTSGL]
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const G4double twoln10
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static const double eplus
Definition: G4SIunits.hh:196
const G4Material * GetMaterial() const