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