Geant4-11
G4BraggModel.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//
29// GEANT4 Class file
30//
31//
32// File name: G4BraggModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 03.01.2002
37//
38// Modifications:
39//
40// 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
41// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42// 27-01-03 Make models region aware (V.Ivanchenko)
43// 13-02-03 Add name (V.Ivanchenko)
44// 04-06-03 Fix compilation warnings (V.Ivanchenko)
45// 12-09-04 Add lowestKinEnergy and change order of if in DEDX method (VI)
46// 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47// 16-06-05 Fix problem of chemical formula (V.Ivantchenko)
48// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49// 25-04-06 Add stopping data from PSTAR (V.Ivanchenko)
50// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
51// CorrectionsAlongStep needed for ions(V.Ivanchenko)
52
53// Class Description:
54//
55// Implementation of energy loss and delta-electron production by
56// slow charged heavy particles
57
58// -------------------------------------------------------------------
59//
60
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65#include "G4BraggModel.hh"
67#include "G4SystemOfUnits.hh"
68#include "Randomize.hh"
69#include "G4Electron.hh"
71#include "G4LossTableManager.hh"
72#include "G4EmCorrections.hh"
73#include "G4EmParameters.hh"
74#include "G4DeltaAngle.hh"
76#include "G4NistManager.hh"
77#include "G4Log.hh"
78#include "G4Exp.hh"
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
83
85 : G4VEmModel(nam),
86 protonMassAMU(1.007276)
87{
89
93 expStopPower125 = 0.0;
94
96 if(nullptr != p) { SetParticle(p); }
97 else { SetParticle(theElectron); }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
103{
104 if(IsMaster()) {
105 delete fPSTAR;
106 fPSTAR = nullptr;
107 }
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111
113 const G4DataVector&)
114{
115 if(p != particle) { SetParticle(p); }
116
117 // always false before the run
118 SetDeexcitationFlag(false);
119
120 if(IsMaster()) {
121 if(nullptr == fPSTAR) { fPSTAR = new G4PSTARStopping(); }
123 if(G4EmParameters::Instance()->UseICRU90Data()) {
124 if(!fICRU90) {
126 } else if(particle->GetPDGMass() < CLHEP::GeV) { fICRU90->Initialise(); }
127 }
128 }
129
130 if(nullptr == fParticleChange) {
131
134 }
136 if(particle->GetParticleType() == "nucleus" &&
137 pname != "deuteron" && pname != "triton" &&
138 pname != "alpha+" && pname != "helium" &&
139 pname != "hydrogen") { isIon = true; }
140
142 }
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
148 const G4Material* mat,
149 G4double kineticEnergy)
150{
151 // this method is called only for ions
152 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
154 return q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
160 const G4MaterialCutsCouple* couple)
161{
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166
168 const G4Material* mat,
169 G4double kineticEnergy)
170{
171 // this method is called only for ions, so no check if it is an ion
172 return corr->GetParticleCharge(p,mat,kineticEnergy);
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176
178 const G4ParticleDefinition* p,
179 G4double kineticEnergy,
180 G4double cut,
181 G4double maxKinEnergy)
182{
183 G4double cross = 0.0;
184 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
185 const G4double maxEnergy = std::min(tmax, maxKinEnergy);
186 const G4double cutEnergy = std::max(cut, lowestKinEnergy*massRate);
187 if(cutEnergy < maxEnergy) {
188
189 const G4double energy = kineticEnergy + mass;
190 const G4double energy2 = energy*energy;
191 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
192 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
193 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
194
195 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
196
197 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
198 }
199 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
200 // << " tmax= " << tmax << " cross= " << cross << G4endl;
201 return cross;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
208 G4double kineticEnergy,
210 G4double cutEnergy,
211 G4double maxEnergy)
212{
213 return
214 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
218
220 const G4ParticleDefinition* p,
221 G4double kineticEnergy,
222 G4double cutEnergy,
223 G4double maxEnergy)
224{
225 return material->GetElectronDensity()
226 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
227}
228
229//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230
232 const G4ParticleDefinition* p,
233 G4double kineticEnergy,
234 G4double cut)
235{
236 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
237 const G4double tkin = kineticEnergy/massRate;
238 const G4double cutEnergy = std::max(cut, lowestKinEnergy*massRate);
239 G4double dedx = 0.0;
240
241 if(tkin < lowestKinEnergy) {
242 dedx = DEDX(material, lowestKinEnergy)*std::sqrt(tkin/lowestKinEnergy);
243 } else {
244 dedx = DEDX(material, tkin);
245
246 if (cutEnergy < tmax) {
247 const G4double tau = kineticEnergy/mass;
248 const G4double x = cutEnergy/tmax;
249
250 dedx += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x) *
251 CLHEP::twopi_mc2_rcl2 * material->GetElectronDensity();
252 }
253 }
254 dedx = std::max(dedx, 0.0) * chargeSquare;
255
256 //G4cout << "E(MeV)= " << tkin/MeV << " dedx= " << dedx
257 // << " " << material->GetName() << G4endl;
258 return dedx;
259}
260
261//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262
263void G4BraggModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
264 const G4MaterialCutsCouple* couple,
265 const G4DynamicParticle* dp,
266 G4double minEnergy,
267 G4double maxEnergy)
268{
269 const G4double tmax = MaxSecondaryKinEnergy(dp);
270 const G4double xmax = std::min(tmax, maxEnergy);
271 const G4double xmin = std::max(lowestKinEnergy*massRate, minEnergy);
272 if(xmin >= xmax) { return; }
273
274 G4double kineticEnergy = dp->GetKineticEnergy();
275 const G4double energy = kineticEnergy + mass;
276 const G4double energy2 = energy*energy;
277 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
278 const G4double grej = 1.0;
279 G4double deltaKinEnergy, f;
280
281 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
282 G4double rndm[2];
283
284 // sampling follows ...
285 do {
286 rndmEngineMod->flatArray(2, rndm);
287 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
288
289 f = 1.0 - beta2*deltaKinEnergy/tmax;
290
291 if(f > grej) {
292 G4cout << "G4BraggModel::SampleSecondary Warning! "
293 << "Majorant " << grej << " < "
294 << f << " for e= " << deltaKinEnergy
295 << G4endl;
296 }
297
298 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
299 } while( grej*rndm[1] >= f );
300
301 G4ThreeVector deltaDirection;
302
304 const G4Material* mat = couple->GetMaterial();
306
307 deltaDirection =
308 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
309
310 } else {
311
312 G4double deltaMomentum =
313 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
314 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
315 (deltaMomentum * dp->GetTotalMomentum());
316 if(cost > 1.0) { cost = 1.0; }
317 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
318
319 G4double phi = twopi*rndmEngineMod->flat();
320
321 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
322 deltaDirection.rotateUz(dp->GetMomentumDirection());
323 }
324
325 // create G4DynamicParticle object for delta ray
326 G4DynamicParticle* delta =
327 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
328
329 // Change kinematics of primary particle
330 kineticEnergy -= deltaKinEnergy;
331 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
332 finalP = finalP.unit();
333
336
337 vdp->push_back(delta);
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
341
343 G4double kinEnergy)
344{
345 if(pd != particle) { SetParticle(pd); }
346 G4double tau = kinEnergy/mass;
347 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
348 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
349 return tmax;
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
353
355{
356 const G4String& chFormula = mat->GetChemicalFormula();
357 if(chFormula.empty()) { return; }
358
359 // ICRU Report N49, 1993. Power's model for H
360 static const size_t numberOfMolecula = 11;
361 static const G4String molName[numberOfMolecula] = {
362 "Al_2O_3", "CO_2", "CH_4",
363 "(C_2H_4)_N-Polyethylene", "(C_2H_4)_N-Polypropylene", "(C_8H_8)_N",
364 "C_3H_8", "SiO_2", "H_2O",
365 "H_2O-Gas", "Graphite" } ;
366
367 // Search for the material in the table
368 for (size_t i=0; i<numberOfMolecula; ++i) {
369 if (chFormula == molName[i]) {
370 iMolecula = i;
371 return;
372 }
373 }
374 return;
375}
376
377//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
378
380 G4double kineticEnergy)
381{
382 G4double ionloss = 0.0 ;
383
384 if (iMolecula >= 0) {
385
386 // The data and the fit from:
387 // ICRU Report N49, 1993. Ziegler's model for protons.
388 // Proton kinetic energy for parametrisation (keV/amu)
389
390 G4double T = kineticEnergy/(keV*protonMassAMU) ;
391
392 static const G4float a[11][5] = {
393 {1.187E+1f, 1.343E+1f, 1.069E+4f, 7.723E+2f, 2.153E-2f},
394 {7.802E+0f, 8.814E+0f, 8.303E+3f, 7.446E+2f, 7.966E-3f},
395 {7.294E+0f, 8.284E+0f, 5.010E+3f, 4.544E+2f, 8.153E-3f},
396 {8.646E+0f, 9.800E+0f, 7.066E+3f, 4.581E+2f, 9.383E-3f},
397 {1.286E+1f, 1.462E+1f, 5.625E+3f, 2.621E+3f, 3.512E-2f},
398 {3.229E+1f, 3.696E+1f, 8.918E+3f, 3.244E+3f, 1.273E-1f},
399 {1.604E+1f, 1.825E+1f, 6.967E+3f, 2.307E+3f, 3.775E-2f},
400 {8.049E+0f, 9.099E+0f, 9.257E+3f, 3.846E+2f, 1.007E-2f},
401 {4.015E+0f, 4.542E+0f, 3.955E+3f, 4.847E+2f, 7.904E-3f},
402 {4.571E+0f, 5.173E+0f, 4.346E+3f, 4.779E+2f, 8.572E-3f},
403 {2.631E+0f, 2.601E+0f, 1.701E+3f, 1.279E+3f, 1.638E-2f} };
404
405 static const G4float atomicWeight[11] = {
406 101.96128f, 44.0098f, 16.0426f, 28.0536f, 42.0804f,
407 104.1512f, 44.665f, 60.0843f, 18.0152f, 18.0152f, 12.0f};
408
409 if ( T < 10.0 ) {
410 ionloss = ((G4double)(a[iMolecula][0])) * std::sqrt(T) ;
411
412 } else if ( T < 10000.0 ) {
413 G4double x1 = (G4double)(a[iMolecula][1]);
414 G4double x2 = (G4double)(a[iMolecula][2]);
415 G4double x3 = (G4double)(a[iMolecula][3]);
416 G4double x4 = (G4double)(a[iMolecula][4]);
417 G4double slow = x1 * G4Exp(G4Log(T)* 0.45);
418 G4double shigh = G4Log( 1.0 + x3/T + x4*T ) * x2/T;
419 ionloss = slow*shigh / (slow + shigh) ;
420 }
421
422 ionloss = std::max(ionloss, 0.0);
423 if ( 10 == iMolecula ) {
424 static const G4double invLog10 = 1.0/G4Log(10.);
425
426 if (T < 100.0) {
427 ionloss *= (1.0+0.023+0.0066*G4Log(T)*invLog10);
428 }
429 else if (T < 700.0) {
430 ionloss *=(1.0+0.089-0.0248*G4Log(T-99.)*invLog10);
431 }
432 else if (T < 10000.0) {
433 ionloss *=(1.0+0.089-0.0248*G4Log(700.-99.)*invLog10);
434 }
435 }
436 ionloss /= (G4double)atomicWeight[iMolecula];
437
438 // pure material (normally not the case for this function)
439 } else if(1 == (material->GetNumberOfElements())) {
440 G4double z = material->GetZ() ;
441 ionloss = ElectronicStoppingPower( z, kineticEnergy ) ;
442 }
443
444 return ionloss;
445}
446
447//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
448
450 G4double kineticEnergy) const
451{
452 G4double ionloss ;
453 G4int i = std::min(std::max(G4lrint(z)-1,0),91); // index of atom
454
455 // The data and the fit from:
456 // ICRU Report 49, 1993. Ziegler's type of parametrisations.
457 // Proton kinetic energy for parametrisation (keV/amu)
458
459 G4double T = kineticEnergy/(keV*protonMassAMU) ;
460
461 static const G4float a[92][5] = {
462 {1.254E+0f, 1.440E+0f, 2.426E+2f, 1.200E+4f, 1.159E-1f},
463 {1.229E+0f, 1.397E+0f, 4.845E+2f, 5.873E+3f, 5.225E-2f},
464 {1.411E+0f, 1.600E+0f, 7.256E+2f, 3.013E+3f, 4.578E-2f},
465 {2.248E+0f, 2.590E+0f, 9.660E+2f, 1.538E+2f, 3.475E-2f},
466 {2.474E+0f, 2.815E+0f, 1.206E+3f, 1.060E+3f, 2.855E-2f},
467 {2.631E+0f, 2.601E+0f, 1.701E+3f, 1.279E+3f, 1.638E-2f},
468 {2.954E+0f, 3.350E+0f, 1.683E+3f, 1.900E+3f, 2.513E-2f},
469 {2.652E+0f, 3.000E+0f, 1.920E+3f, 2.000E+3f, 2.230E-2f},
470 {2.085E+0f, 2.352E+0f, 2.157E+3f, 2.634E+3f, 1.816E-2f},
471 {1.951E+0f, 2.199E+0f, 2.393E+3f, 2.699E+3f, 1.568E-2f},
472 // Z= 11-20
473 {2.542E+0f, 2.869E+0f, 2.628E+3f, 1.854E+3f, 1.472E-2f},
474 {3.791E+0f, 4.293E+0f, 2.862E+3f, 1.009E+3f, 1.397E-2f},
475 {4.154E+0f, 4.739E+0f, 2.766E+3f, 1.645E+2f, 2.023E-2f},
476 {4.914E+0f, 5.598E+0f, 3.193E+3f, 2.327E+2f, 1.419E-2f},
477 {3.232E+0f, 3.647E+0f, 3.561E+3f, 1.560E+3f, 1.267E-2f},
478 {3.447E+0f, 3.891E+0f, 3.792E+3f, 1.219E+3f, 1.211E-2f},
479 {5.301E+0f, 6.008E+0f, 3.969E+3f, 6.451E+2f, 1.183E-2f},
480 {5.731E+0f, 6.500E+0f, 4.253E+3f, 5.300E+2f, 1.123E-2f},
481 {5.152E+0f, 5.833E+0f, 4.482E+3f, 5.457E+2f, 1.129E-2f},
482 {5.521E+0f, 6.252E+0f, 4.710E+3f, 5.533E+2f, 1.112E-2f},
483 // Z= 21-30
484 {5.201E+0f, 5.884E+0f, 4.938E+3f, 5.609E+2f, 9.995E-3f},
485 {4.858E+0f, 5.489E+0f, 5.260E+3f, 6.511E+2f, 8.930E-3f},
486 {4.479E+0f, 5.055E+0f, 5.391E+3f, 9.523E+2f, 9.117E-3f},
487 {3.983E+0f, 4.489E+0f, 5.616E+3f, 1.336E+3f, 8.413E-3f},
488 {3.469E+0f, 3.907E+0f, 5.725E+3f, 1.461E+3f, 8.829E-3f},
489 {3.519E+0f, 3.963E+0f, 6.065E+3f, 1.243E+3f, 7.782E-3f},
490 {3.140E+0f, 3.535E+0f, 6.288E+3f, 1.372E+3f, 7.361E-3f},
491 {3.553E+0f, 4.004E+0f, 6.205E+3f, 5.551E+2f, 8.763E-3f},
492 {3.696E+0f, 4.194E+0f, 4.649E+3f, 8.113E+1f, 2.242E-2f},
493 {4.210E+0f, 4.750E+0f, 6.953E+3f, 2.952E+2f, 6.809E-3f},
494 // Z= 31-40
495 {5.041E+0f, 5.697E+0f, 7.173E+3f, 2.026E+2f, 6.725E-3f},
496 {5.554E+0f, 6.300E+0f, 6.496E+3f, 1.100E+2f, 9.689E-3f},
497 {5.323E+0f, 6.012E+0f, 7.611E+3f, 2.925E+2f, 6.447E-3f},
498 {5.874E+0f, 6.656E+0f, 7.395E+3f, 1.175E+2f, 7.684E-3f},
499 {6.658E+0f, 7.536E+0f, 7.694E+3f, 2.223E+2f, 6.509E-3f},
500 {6.413E+0f, 7.240E+0f, 1.185E+4f, 1.537E+2f, 2.880E-3f},
501 {5.694E+0f, 6.429E+0f, 8.478E+3f, 2.929E+2f, 6.087E-3f},
502 {6.339E+0f, 7.159E+0f, 8.693E+3f, 3.303E+2f, 6.003E-3f},
503 {6.407E+0f, 7.234E+0f, 8.907E+3f, 3.678E+2f, 5.889E-3f},
504 {6.734E+0f, 7.603E+0f, 9.120E+3f, 4.052E+2f, 5.765E-3f},
505 // Z= 41-50
506 {6.901E+0f, 7.791E+0f, 9.333E+3f, 4.427E+2f, 5.587E-3f},
507 {6.424E+0f, 7.248E+0f, 9.545E+3f, 4.802E+2f, 5.376E-3f},
508 {6.799E+0f, 7.671E+0f, 9.756E+3f, 5.176E+2f, 5.315E-3f},
509 {6.109E+0f, 6.887E+0f, 9.966E+3f, 5.551E+2f, 5.151E-3f},
510 {5.924E+0f, 6.677E+0f, 1.018E+4f, 5.925E+2f, 4.919E-3f},
511 {5.238E+0f, 5.900E+0f, 1.038E+4f, 6.300E+2f, 4.758E-3f},
512 // {5.623f, 6.354f, 7160.0f, 337.6f, 0.013940f}, // Ag Ziegler77
513 {5.345E+0f, 6.038E+0f, 6.790E+3f, 3.978E+2f, 1.676E-2f}, // Ag ICRU49
514 {5.814E+0f, 6.554E+0f, 1.080E+4f, 3.555E+2f, 4.626E-3f},
515 {6.229E+0f, 7.024E+0f, 1.101E+4f, 3.709E+2f, 4.540E-3f},
516 {6.409E+0f, 7.227E+0f, 1.121E+4f, 3.864E+2f, 4.474E-3f},
517 // Z= 51-60
518 {7.500E+0f, 8.480E+0f, 8.608E+3f, 3.480E+2f, 9.074E-3f},
519 {6.979E+0f, 7.871E+0f, 1.162E+4f, 3.924E+2f, 4.402E-3f},
520 {7.725E+0f, 8.716E+0f, 1.183E+4f, 3.948E+2f, 4.376E-3f},
521 {8.337E+0f, 9.425E+0f, 1.051E+4f, 2.696E+2f, 6.206E-3f},
522 {7.287E+0f, 8.218E+0f, 1.223E+4f, 3.997E+2f, 4.447E-3f},
523 {7.899E+0f, 8.911E+0f, 1.243E+4f, 4.021E+2f, 4.511E-3f},
524 {8.041E+0f, 9.071E+0f, 1.263E+4f, 4.045E+2f, 4.540E-3f},
525 {7.488E+0f, 8.444E+0f, 1.283E+4f, 4.069E+2f, 4.420E-3f},
526 {7.291E+0f, 8.219E+0f, 1.303E+4f, 4.093E+2f, 4.298E-3f},
527 {7.098E+0f, 8.000E+0f, 1.323E+4f, 4.118E+2f, 4.182E-3f},
528 // Z= 61-70
529 {6.909E+0f, 7.786E+0f, 1.343E+4f, 4.142E+2f, 4.058E-3f},
530 {6.728E+0f, 7.580E+0f, 1.362E+4f, 4.166E+2f, 3.976E-3f},
531 {6.551E+0f, 7.380E+0f, 1.382E+4f, 4.190E+2f, 3.877E-3f},
532 {6.739E+0f, 7.592E+0f, 1.402E+4f, 4.214E+2f, 3.863E-3f},
533 {6.212E+0f, 6.996E+0f, 1.421E+4f, 4.239E+2f, 3.725E-3f},
534 {5.517E+0f, 6.210E+0f, 1.440E+4f, 4.263E+2f, 3.632E-3f},
535 {5.220E+0f, 5.874E+0f, 1.460E+4f, 4.287E+2f, 3.498E-3f},
536 {5.071E+0f, 5.706E+0f, 1.479E+4f, 4.330E+2f, 3.405E-3f},
537 {4.926E+0f, 5.542E+0f, 1.498E+4f, 4.335E+2f, 3.342E-3f},
538 {4.788E+0f, 5.386E+0f, 1.517E+4f, 4.359E+2f, 3.292E-3f},
539 // Z= 71-80
540 {4.893E+0f, 5.505E+0f, 1.536E+4f, 4.384E+2f, 3.243E-3f},
541 {5.028E+0f, 5.657E+0f, 1.555E+4f, 4.408E+2f, 3.195E-3f},
542 {4.738E+0f, 5.329E+0f, 1.574E+4f, 4.432E+2f, 3.186E-3f},
543 {4.587E+0f, 5.160E+0f, 1.541E+4f, 4.153E+2f, 3.406E-3f},
544 {5.201E+0f, 5.851E+0f, 1.612E+4f, 4.416E+2f, 3.122E-3f},
545 {5.071E+0f, 5.704E+0f, 1.630E+4f, 4.409E+2f, 3.082E-3f},
546 {4.946E+0f, 5.563E+0f, 1.649E+4f, 4.401E+2f, 2.965E-3f},
547 {4.477E+0f, 5.034E+0f, 1.667E+4f, 4.393E+2f, 2.871E-3f},
548 // {4.856f, 5.460f, 18320.0f, 438.5f, 0.002542f}, //Ziegler77
549 {4.844E+0f, 5.458E+0f, 7.852E+3f, 9.758E+2f, 2.077E-2f}, //ICRU49
550 {4.307E+0f, 4.843E+0f, 1.704E+4f, 4.878E+2f, 2.882E-3f},
551 // Z= 81-90
552 {4.723E+0f, 5.311E+0f, 1.722E+4f, 5.370E+2f, 2.913E-3f},
553 {5.319E+0f, 5.982E+0f, 1.740E+4f, 5.863E+2f, 2.871E-3f},
554 {5.956E+0f, 6.700E+0f, 1.780E+4f, 6.770E+2f, 2.660E-3f},
555 {6.158E+0f, 6.928E+0f, 1.777E+4f, 5.863E+2f, 2.812E-3f},
556 {6.203E+0f, 6.979E+0f, 1.795E+4f, 5.863E+2f, 2.776E-3f},
557 {6.181E+0f, 6.954E+0f, 1.812E+4f, 5.863E+2f, 2.748E-3f},
558 {6.949E+0f, 7.820E+0f, 1.830E+4f, 5.863E+2f, 2.737E-3f},
559 {7.506E+0f, 8.448E+0f, 1.848E+4f, 5.863E+2f, 2.727E-3f},
560 {7.648E+0f, 8.609E+0f, 1.866E+4f, 5.863E+2f, 2.697E-3f},
561 {7.711E+0f, 8.679E+0f, 1.883E+4f, 5.863E+2f, 2.641E-3f},
562 // Z= 91-92
563 {7.407E+0f, 8.336E+0f, 1.901E+4f, 5.863E+2f, 2.603E-3f},
564 {7.290E+0f, 8.204E+0f, 1.918E+4f, 5.863E+2f, 2.673E-3f}
565 };
566
567 G4double fac = 1.0 ;
568
569 // Carbon specific case for E < 40 keV
570 if ( T < 40.0 && 5 == i) {
571 fac = std::sqrt(T*0.025);
572 T = 40.0;
573
574 // Free electron gas model
575 } else if ( T < 10.0 ) {
576 fac = std::sqrt(T*0.1) ;
577 T = 10.0;
578 }
579
580 // Main parametrisation
581 G4double x1 = (G4double)(a[i][1]);
582 G4double x2 = (G4double)(a[i][2]);
583 G4double x3 = (G4double)(a[i][3]);
584 G4double x4 = (G4double)(a[i][4]);
585 G4double slow = x1 * G4Exp(G4Log(T) * 0.45);
586 G4double shigh = G4Log( 1.0 + x3/T + x4*T ) * x2/T;
587 ionloss = slow*shigh*fac / (slow + shigh);
588
589 ionloss = std::max(ionloss, 0.0);
590
591 return ionloss;
592}
593
594//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
595
597{
598 G4double eloss = 0.0;
599
600 // check DB
603 baseMaterial = material->GetBaseMaterial()
604 ? material->GetBaseMaterial() : material;
605 iPSTAR = -1;
606 iMolecula = -1;
608
609 if(iICRU90 < 0) {
611 if(iPSTAR < 0) { HasMaterial(baseMaterial); }
612 }
613 //G4cout << "%%% " <<material->GetName() << " iMolecula= "
614 // << iMolecula << " iPSTAR= " << iPSTAR
615 // << " iICRU90= " << iICRU90<< G4endl;
616 }
617
618 // ICRU90 parameterisation
619 if(iICRU90 >= 0) {
620 return fICRU90->GetElectronicDEDXforProton(iICRU90, kineticEnergy)
621 *material->GetDensity();
622 }
623 // PSTAR parameterisation
624 if( iPSTAR >= 0 ) {
625 return fPSTAR->GetElectronicDEDX(iPSTAR, kineticEnergy)
626 *material->GetDensity();
627
628 }
629 const G4int numberOfElements = material->GetNumberOfElements();
630 const G4double* theAtomicNumDensityVector =
631 material->GetAtomicNumDensityVector();
632
633
634 if(iMolecula >= 0) {
635 eloss = StoppingPower(baseMaterial, kineticEnergy)*
636 material->GetDensity()/amu;
637
638 // Pure material ICRU49 paralmeterisation
639 } else if(1 == numberOfElements) {
640
641 G4double z = material->GetZ();
642 eloss = ElectronicStoppingPower(z, kineticEnergy)
643 * (material->GetTotNbOfAtomsPerVolume());
644
645
646 // Experimental data exist only for kinetic energy 125 keV
647 } else if( MolecIsInZiegler1988(material) ) {
648
649 // Loop over elements - calculation based on Bragg's rule
650 G4double eloss125 = 0.0 ;
651 const G4ElementVector* theElementVector =
652 material->GetElementVector();
653
654 // Loop for the elements in the material
655 for (G4int i=0; i<numberOfElements; ++i) {
656 const G4Element* element = (*theElementVector)[i] ;
657 G4double z = element->GetZ() ;
658 eloss += ElectronicStoppingPower(z,kineticEnergy)
659 * theAtomicNumDensityVector[i] ;
660 eloss125 += ElectronicStoppingPower(z,125.0*keV)
661 * theAtomicNumDensityVector[i] ;
662 }
663
664 // Chemical factor is taken into account
665 eloss *= ChemicalFactor(kineticEnergy, eloss125) ;
666
667 // Brugg's rule calculation
668 } else {
669 const G4ElementVector* theElementVector =
670 material->GetElementVector() ;
671
672 // loop for the elements in the material
673 for (G4int i=0; i<numberOfElements; ++i)
674 {
675 const G4Element* element = (*theElementVector)[i] ;
676 eloss += ElectronicStoppingPower(element->GetZ(), kineticEnergy)
677 * theAtomicNumDensityVector[i];
678 }
679 }
680 return eloss*theZieglerFactor;
681}
682
683//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
684
686{
687 // The list of molecules from
688 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
689 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
690
691 G4String myFormula = G4String(" ") ;
692 const G4String chFormula = material->GetChemicalFormula() ;
693 if (myFormula == chFormula ) { return false; }
694
695 // There are no evidence for difference of stopping power depended on
696 // phase of the compound except for water. The stopping power of the
697 // water in gas phase can be predicted using Bragg's rule.
698 //
699 // No chemical factor for water-gas
700
701 myFormula = G4String("H_2O") ;
702 const G4State theState = material->GetState() ;
703 if( theState == kStateGas && myFormula == chFormula) return false ;
704
705
706 // The coffecient from Table.4 of Ziegler & Manoyan
707 static const G4float HeEff = 2.8735f;
708
709 static const size_t numberOfMolecula = 53;
710 static const G4String nameOfMol[numberOfMolecula] = {
711 "H_2O", "C_2H_4O", "C_3H_6O", "C_2H_2", "C_H_3OH",
712 "C_2H_5OH", "C_3H_7OH", "C_3H_4", "NH_3", "C_14H_10",
713 "C_6H_6", "C_4H_10", "C_4H_6", "C_4H_8O", "CCl_4",
714 "CF_4", "C_6H_8", "C_6H_12", "C_6H_10O", "C_6H_10",
715 "C_8H_16", "C_5H_10", "C_5H_8", "C_3H_6-Cyclopropane","C_2H_4F_2",
716 "C_2H_2F_2", "C_4H_8O_2", "C_2H_6", "C_2F_6", "C_2H_6O",
717 "C_3H_6O", "C_4H_10O", "C_2H_4", "C_2H_4O", "C_2H_4S",
718 "SH_2", "CH_4", "CCLF_3", "CCl_2F_2", "CHCl_2F",
719 "(CH_3)_2S", "N_2O", "C_5H_10O", "C_8H_6", "(CH_2)_N",
720 "(C_3H_6)_N","(C_8H_8)_N", "C_3H_8", "C_3H_6-Propylene", "C_3H_6O",
721 "C_3H_6S", "C_4H_4S", "C_7H_8"
722 };
723
724 static const G4float expStopping[numberOfMolecula] = {
725 66.1f, 190.4f, 258.7f, 42.2f, 141.5f,
726 210.9f, 279.6f, 198.8f, 31.0f, 267.5f,
727 122.8f, 311.4f, 260.3f, 328.9f, 391.3f,
728 206.6f, 374.0f, 422.0f, 432.0f, 398.0f,
729 554.0f, 353.0f, 326.0f, 74.6f, 220.5f,
730 197.4f, 362.0f, 170.0f, 330.5f, 211.3f,
731 262.3f, 349.6f, 51.3f, 187.0f, 236.9f,
732 121.9f, 35.8f, 247.0f, 292.6f, 268.0f,
733 262.3f, 49.0f, 398.9f, 444.0f, 22.91f,
734 68.0f, 155.0f, 84.0f, 74.2f, 254.7f,
735 306.8f, 324.4f, 420.0f
736 } ;
737
738 static const G4float expCharge[53] = {
739 HeEff, HeEff, HeEff, 1.0f, HeEff,
740 HeEff, HeEff, HeEff, 1.0f, 1.0f,
741 1.0f, HeEff, HeEff, HeEff, HeEff,
742 HeEff, HeEff, HeEff, HeEff, HeEff,
743 HeEff, HeEff, HeEff, 1.0f, HeEff,
744 HeEff, HeEff, HeEff, HeEff, HeEff,
745 HeEff, HeEff, 1.0f, HeEff, HeEff,
746 HeEff, 1.0f, HeEff, HeEff, HeEff,
747 HeEff, 1.0f, HeEff, HeEff, 1.0f,
748 1.0f, 1.0f, 1.0f, 1.0f, HeEff,
749 HeEff, HeEff, HeEff
750 } ;
751
752 static const G4int numberOfAtomsPerMolecula[53] = {
753 3, 7, 10, 4, 6,
754 9, 12, 7, 4, 24,
755 12,14, 10, 13, 5,
756 5, 14, 18, 17, 17,
757 24,15, 13, 9, 8,
758 6, 14, 8, 8, 9,
759 10,15, 6, 7, 7,
760 3, 5, 5, 5, 5,
761 9, 3, 16, 14, 3,
762 9, 16, 11, 9, 10,
763 10, 9, 15};
764
765 // Search for the compaund in the table
766 for (size_t i=0; i<numberOfMolecula; ++i) {
767 if(chFormula == nameOfMol[i]) {
768 expStopPower125 = ((G4double)expStopping[i])
769 * (material->GetTotNbOfAtomsPerVolume()) /
770 ((G4double)(expCharge[i] * numberOfAtomsPerMolecula[i]));
771 return true;
772 }
773 }
774 return false;
775}
776
777//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
778
780 G4double eloss125) const
781{
782 // Approximation of Chemical Factor according to
783 // J.F.Ziegler and J.M.Manoyan, The stopping of ions in compaunds,
784 // Nucl. Inst. & Meth. in Phys. Res. B35 (1988) 215-228.
785
786 static const G4double gamma25 = 1.0 + 25.0*keV /proton_mass_c2;
787 static const G4double gamma125 = 1.0 + 125.0*keV/proton_mass_c2;
788 static const G4double beta25 = std::sqrt(1.0 - 1.0/(gamma25*gamma25));
789 static const G4double beta125 = std::sqrt(1.0 - 1.0/(gamma125*gamma125));
790 static const G4double f12525 = 1.0 + G4Exp( 1.48*(beta125/beta25 - 7.0) );
791
792 G4double gamma = 1.0 + kineticEnergy/proton_mass_c2;
793 G4double beta = std::sqrt(1.0 - 1.0/(gamma*gamma));
794
795 G4double factor = 1.0 + (expStopPower125/eloss125 - 1.0) * f12525/
796 (1.0 + G4Exp( 1.48 * ( beta/beta25 - 7.0 ) ) );
797
798 return factor ;
799}
800
801//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
802
std::vector< const G4Element * > G4ElementVector
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4State
Definition: G4Material.hh:111
@ kStateGas
Definition: G4Material.hh:111
static const G4double fac
static const G4int numberOfMolecula
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double keV
Definition: G4SIunits.hh:202
float G4float
Definition: G4Types.hh:84
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
G4ParticleChangeForLoss * fParticleChange
G4double protonMassAMU
G4double ratio
G4double ElectronicStoppingPower(G4double z, G4double kineticEnergy) const
~G4BraggModel() override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
const G4Material * baseMaterial
G4double massRate
G4double lowestKinEnergy
G4BraggModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")
Definition: G4BraggModel.cc:84
G4bool MolecIsInZiegler1988(const G4Material *material)
G4double theZieglerFactor
G4double StoppingPower(const G4Material *material, G4double kineticEnergy)
void SetParticle(const G4ParticleDefinition *p)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4EmCorrections * corr
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4ParticleDefinition * theElectron
G4ICRU90StoppingData * fICRU90
G4double mass
G4double DEDX(const G4Material *material, G4double kineticEnergy)
G4double chargeSquare
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
const G4ParticleDefinition * particle
G4double GetChargeSquareRatio() const
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
const G4Material * currentMaterial
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double expStopPower125
G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const
static G4PSTARStopping * fPSTAR
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void HasMaterial(const G4Material *material)
G4double spin
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double GetZ() const
Definition: G4Element.hh:131
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4EmParameters * Instance()
G4double GetElectronicDEDXforProton(const G4Material *, G4double kinEnergy) const
G4int GetIndex(const G4Material *) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
const G4String & GetChemicalFormula() const
Definition: G4Material.hh:174
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
G4ICRU90StoppingData * GetICRU90StoppingData()
static G4NistManager * Instance()
G4int GetIndex(const G4Material *) const
G4double GetElectronicDEDX(G4int idx, G4double energy) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleType() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:614
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:621
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:297
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:823
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:718
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:520
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:108
static constexpr double GeV
static constexpr double keV
static constexpr double MeV
static constexpr double cm2
static constexpr double twopi_mc2_rcl2
static constexpr double eV
G4double energy(const ThreeVector &p, const G4double m)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
string material
Definition: eplot.py:19
string pname
Definition: eplot.py:33
float electron_mass_c2
Definition: hepunit.py:273
float proton_mass_c2
Definition: hepunit.py:274
int G4lrint(double ad)
Definition: templates.hh:134