Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmCalculator.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: G4EmCalculator
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 28.06.2004
37//
38//
39// Class Description: V.Ivanchenko & M.Novak
40//
41// -------------------------------------------------------------------
42//
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45
46#include "G4EmCalculator.hh"
47#include "G4SystemOfUnits.hh"
48#include "G4LossTableManager.hh"
49#include "G4EmParameters.hh"
50#include "G4NistManager.hh"
51#include "G4DynamicParticle.hh"
52#include "G4VEmProcess.hh"
55#include "G4Material.hh"
58#include "G4ParticleTable.hh"
59#include "G4IonTable.hh"
60#include "G4PhysicsTable.hh"
62#include "G4ProcessManager.hh"
64#include "G4RegionStore.hh"
65#include "G4Element.hh"
66#include "G4EmCorrections.hh"
67#include "G4GenericIon.hh"
68#include "G4ProcessVector.hh"
69#include "G4Gamma.hh"
70#include "G4Electron.hh"
71#include "G4Positron.hh"
72#include "G4EmUtility.hh"
73
74namespace {
75 constexpr G4double length = CLHEP::nm;
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81{
84 theParameters = G4EmParameters::Instance();
85 corr = manager->EmCorrections();
86 cutenergy[0] = cutenergy[1] = cutenergy[2] = DBL_MAX;
87 theGenericIon = G4GenericIon::GenericIon();
88 ionEffCharge = new G4ionEffectiveCharge();
89 dynParticle = new G4DynamicParticle();
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96{
97 delete ionEffCharge;
98 delete dynParticle;
99 for (G4int i=0; i<nLocalMaterials; ++i) {
100 delete localCouples[i];
101 }
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107 const G4ParticleDefinition* p,
108 const G4Material* mat,
109 const G4Region* region)
110{
111 G4double res = 0.0;
112 const G4MaterialCutsCouple* couple = FindCouple(mat, region);
113 if (nullptr != couple && UpdateParticle(p, kinEnergy) ) {
114 if (FindEmModel(p, currentProcessName, kinEnergy)) {
115 if (isIon) {
116 G4double q2 = currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
117 res = manager->GetDEDX(p, kinEnergy, couple) * q2;
118 G4double eloss = res*length;
119 auto const v =
121 G4double cut = (*v)[couple->GetIndex()];
122 currentModel->CorrectionsAlongStep(mat, p, kinEnergy, cut, length, eloss);
123 res = std::max(eloss/length, 0.0);
124 } else {
125 res = manager->GetDEDX(p, kinEnergy, couple);
126 }
127 }
128 }
129 if (verbose > 1) {
130 G4cout << "G4EmCalculator::GetDEDX: E(MeV)= " << kinEnergy/MeV
131 << " DEDX(MeV/mm)= " << res*mm/MeV
132 << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
133 << " " << p->GetParticleName()
134 << " in " << mat->GetName()
135 << " isIon= " << isIon
136 << G4endl;
137 }
138 return res;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144 const G4ParticleDefinition* p,
145 const G4Material* mat,
146 const G4Region* region)
147{
148 G4double res = 0.0;
149 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
150 if(couple && UpdateParticle(p, kinEnergy)) {
151 res = manager->GetRangeFromRestricteDEDX(p, kinEnergy, couple);
152 if(verbose>1) {
153 G4cout << " G4EmCalculator::GetRangeFromRestrictedDEDX: E(MeV)= "
154 << kinEnergy/MeV
155 << " range(mm)= " << res/mm
156 << " " << p->GetParticleName()
157 << " in " << mat->GetName()
158 << G4endl;
159 }
160 }
161 return res;
162}
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165
167 const G4ParticleDefinition* p,
168 const G4Material* mat,
169 const G4Region* region)
170{
171 G4double res = 0.0;
172 if(!theParameters->BuildCSDARange()) {
174 ed << "G4EmCalculator::GetCSDARange: CSDA table is not built; "
175 << " use UI command: /process/eLoss/CSDARange true";
176 G4Exception("G4EmCalculator::GetCSDARange", "em0077",
177 JustWarning, ed);
178 return res;
179 }
180
181 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
182 if(nullptr != couple && UpdateParticle(p, kinEnergy)) {
183 res = manager->GetCSDARange(p, kinEnergy, couple);
184 if(verbose>1) {
185 G4cout << " G4EmCalculator::GetCSDARange: E(MeV)= " << kinEnergy/MeV
186 << " range(mm)= " << res/mm
187 << " " << p->GetParticleName()
188 << " in " << mat->GetName()
189 << G4endl;
190 }
191 }
192 return res;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
198 const G4ParticleDefinition* p,
199 const G4Material* mat,
200 const G4Region* region)
201{
202 G4double res = 0.0;
203 if(theParameters->BuildCSDARange()) {
204 res = GetCSDARange(kinEnergy, p, mat, region);
205 } else {
206 res = GetRangeFromRestricteDEDX(kinEnergy, p, mat, region);
207 }
208 return res;
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
212
214 const G4ParticleDefinition* p,
215 const G4Material* mat,
216 const G4Region* region)
217{
218 G4double res = 0.0;
219 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
220 if(nullptr != couple && UpdateParticle(p, 1.0*GeV)) {
221 res = manager->GetEnergy(p, range, couple);
222 if(verbose>0) {
223 G4cout << "G4EmCalculator::GetKinEnergy: Range(mm)= " << range/mm
224 << " KinE(MeV)= " << res/MeV
225 << " " << p->GetParticleName()
226 << " in " << mat->GetName()
227 << G4endl;
228 }
229 }
230 return res;
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234
236 const G4ParticleDefinition* p,
237 const G4String& processName,
238 const G4Material* mat,
239 const G4Region* region)
240{
241 G4double res = 0.0;
242 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
243
244 if(nullptr != couple && UpdateParticle(p, kinEnergy)) {
245 if(FindEmModel(p, processName, kinEnergy)) {
246 G4int idx = couple->GetIndex();
247 G4int procType = -1;
248 FindLambdaTable(p, processName, kinEnergy, procType);
249
250 G4VEmProcess* emproc = FindDiscreteProcess(p, processName);
251 if(nullptr != emproc) {
252 res = emproc->GetCrossSection(kinEnergy, couple);
253 } else if(currentLambda) {
254 // special tables are built for Msc models
255 // procType is set in FindLambdaTable
256 if(procType==2) {
257 auto mscM = static_cast<G4VMscModel*>(currentModel);
258 mscM->SetCurrentCouple(couple);
259 G4double tr1Mfp = mscM->GetTransportMeanFreePath(p, kinEnergy);
260 if (tr1Mfp<DBL_MAX) {
261 res = 1./tr1Mfp;
262 }
263 } else {
264 G4double e = kinEnergy*massRatio;
265 res = (((*currentLambda)[idx])->Value(e))*chargeSquare;
266 }
267 } else {
268 res = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, kinEnergy);
269 }
270 if(verbose>0) {
271 G4cout << "G4EmCalculator::GetXSPerVolume: E(MeV)= " << kinEnergy/MeV
272 << " cross(cm-1)= " << res*cm
273 << " " << p->GetParticleName()
274 << " in " << mat->GetName();
275 if(verbose>1)
276 G4cout << " idx= " << idx << " Escaled((MeV)= "
277 << kinEnergy*massRatio
278 << " q2= " << chargeSquare;
279 G4cout << G4endl;
280 }
281 }
282 }
283 return res;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
289 const G4String& particle,
290 G4int Z,
292 G4double kinEnergy)
293{
294 G4double res = 0.0;
295 const G4ParticleDefinition* p = FindParticle(particle);
296 G4VAtomDeexcitation* ad = manager->AtomDeexcitation();
297 if(nullptr != p && nullptr != ad) {
298 res = ad->GetShellIonisationCrossSectionPerAtom(p, Z, shell, kinEnergy);
299 }
300 return res;
301}
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
304
306 const G4ParticleDefinition* p,
307 const G4String& processName,
308 const G4Material* mat,
309 const G4Region* region)
310{
311 G4double res = DBL_MAX;
312 G4double x = GetCrossSectionPerVolume(kinEnergy,p, processName, mat,region);
313 if(x > 0.0) { res = 1.0/x; }
314 if(verbose>1) {
315 G4cout << "G4EmCalculator::GetMeanFreePath: E(MeV)= " << kinEnergy/MeV
316 << " MFP(mm)= " << res/mm
317 << " " << p->GetParticleName()
318 << " in " << mat->GetName()
319 << G4endl;
320 }
321 return res;
322}
323
324//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325
327{
328 const G4VEnergyLossProcess* elp = manager->GetEnergyLossProcess(p);
329 G4cout << "##### DEDX Table for " << p->GetParticleName() << G4endl;
330 if(nullptr != elp) G4cout << *(elp->DEDXTable()) << G4endl;
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334
336{
337 const G4VEnergyLossProcess* elp = manager->GetEnergyLossProcess(p);
338 G4cout << "##### Range Table for " << p->GetParticleName() << G4endl;
339 if(nullptr != elp) G4cout << *(elp->RangeTableForLoss()) << G4endl;
340}
341
342//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343
345{
346 const G4VEnergyLossProcess* elp = manager->GetEnergyLossProcess(p);
347 G4cout << "### G4EmCalculator: Inverse Range Table for "
348 << p->GetParticleName() << G4endl;
349 if(nullptr != elp) G4cout << *(elp->InverseRangeTable()) << G4endl;
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353
355 const G4ParticleDefinition* p,
356 const G4String& processName,
357 const G4Material* mat,
358 G4double cut)
359{
360 SetupMaterial(mat);
361 G4double res = 0.0;
362 if (verbose > 1) {
363 G4cout << "### G4EmCalculator::ComputeDEDX: " << p->GetParticleName()
364 << " in " << currentMaterialName
365 << " e(MeV)= " << kinEnergy/MeV << " cut(MeV)= " << cut/MeV
366 << G4endl;
367 }
368 if (UpdateParticle(p, kinEnergy)) {
369 if (FindEmModel(p, processName, kinEnergy)) {
370 G4double escaled = kinEnergy*massRatio;
371 chargeSquare = currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
372 if (nullptr != baseParticle && !isIon) {
373 res = currentModel->ComputeDEDXPerVolume(mat, baseParticle, escaled, cut)
374 * chargeSquare;
375 if (verbose > 1) {
376 G4cout << "Particle: " << p->GetParticleName()
377 << " E(MeV)=" << kinEnergy
378 << " Base particle: " << baseParticle->GetParticleName()
379 << " Escaled(MeV)= " << escaled
380 << " q2=" << chargeSquare
381 << " dEdx(MeV/mm)=" << res << " << cut=" << cut << G4endl;
382 }
383 } else {
384 res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
385 if(verbose > 1) {
386 G4cout << p->GetParticleName() << " E(MeV)=" << kinEnergy
387 << " dE/dx(MeV/mm)=" << res << G4endl;
388 }
389 }
390 if (verbose > 1) {
391 G4cout << currentModel->GetName() << ": E(MeV)=" << kinEnergy
392 << " cut(MeV)=" << cut << " DEDX(MeV/mm)=" << res
393 << " q2=" << chargeSquare << G4endl;
394 }
395 // emulate smoothing procedure
396 if (applySmoothing && nullptr != loweModel && loweModel != currentModel) {
397 G4double eth = loweModel->HighEnergyLimit();
398 G4double res0 = 0.0;
399 G4double res1 = 0.0;
400 if (nullptr != baseParticle && !isIon) {
401 res1 = chargeSquare*
402 currentModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut);
403 res0 = chargeSquare*
404 loweModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut);
405 } else {
406 res1 = currentModel->ComputeDEDXPerVolume(mat, p, eth, cut);
407 res0 = loweModel->ComputeDEDXPerVolume(mat, p, eth, cut);
408 }
409 if(res1 > 0.0 && escaled > 0.0) {
410 res *= (1.0 + (res0/res1 - 1.0)*eth/escaled);
411 }
412 if (verbose > 1) {
413 G4cout << "At boundary energy(MeV)= " << eth/MeV
414 << " DEDX(MeV/mm)= " << res0*mm/MeV << " " << res1*mm/MeV
415 << " after correction DEDX(MeV/mm)=" << res*mm/MeV << G4endl;
416 }
417 }
418 if (isIon) {
419 G4double eloss = res*length;
420 currentModel->CorrectionsAlongStep(mat, p, kinEnergy, cut, length, eloss);
421 res = std::max(eloss/length, 0.0);
422 }
423 if (verbose > 0) {
424 G4cout << "## E(MeV)= " << kinEnergy/MeV
425 << " DEDX(MeV/mm)= " << res*mm/MeV
426 << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
427 << " cut(MeV)= " << cut/MeV
428 << " " << p->GetParticleName()
429 << " in " << currentMaterialName
430 << " Zi^2= " << chargeSquare
431 << " isIon=" << isIon
432 << G4endl;
433 }
434 }
435 }
436 return res;
437}
438
439//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
440
442 const G4ParticleDefinition* part,
443 const G4Material* mat,
444 G4double cut)
445{
446 SetupMaterial(mat);
447 G4double dedx = 0.0;
448 if(UpdateParticle(part, kinEnergy)) {
449
451 const std::vector<G4VEnergyLossProcess*> vel =
452 lManager->GetEnergyLossProcessVector();
453 std::size_t n = vel.size();
454
455 //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
456 // << " n= " << n << G4endl;
457
458 for(std::size_t i=0; i<n; ++i) {
459 if(vel[i]) {
460 auto p = static_cast<G4VProcess*>(vel[i]);
461 if(ActiveForParticle(part, p)) {
462 //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
463 // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
464 dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),mat,cut);
465 }
466 }
467 }
468 }
469 return dedx;
470}
471
472//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473
476 const G4ParticleDefinition* part,
477 const G4Material* mat,
478 G4double rangecut)
479{
480 SetupMaterial(mat);
481 G4double dedx = 0.0;
482 if(UpdateParticle(part, kinEnergy)) {
483
485 const std::vector<G4VEnergyLossProcess*> vel =
486 lManager->GetEnergyLossProcessVector();
487 std::size_t n = vel.size();
488
489 if(mat != cutMaterial) {
490 cutMaterial = mat;
491 cutenergy[0] =
493 cutenergy[1] =
495 cutenergy[2] =
497 }
498
499 //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
500 // << " n= " << n << G4endl;
501
502 for(std::size_t i=0; i<n; ++i) {
503 if(vel[i]) {
504 auto p = static_cast<G4VProcess*>(vel[i]);
505 if(ActiveForParticle(part, p)) {
506 //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
507 // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
508 const G4ParticleDefinition* sec = (vel[i])->SecondaryParticle();
509 std::size_t idx = 0;
510 if(sec == G4Electron::Electron()) { idx = 1; }
511 else if(sec == G4Positron::Positron()) { idx = 2; }
512
513 dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),
514 mat,cutenergy[idx]);
515 }
516 }
517 }
518 }
519 return dedx;
520}
521
522//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523
525 const G4ParticleDefinition* part,
526 const G4Material* mat,
527 G4double cut)
528{
529 G4double dedx = ComputeElectronicDEDX(kinEnergy,part,mat,cut);
530 if(mass > 700.*MeV) { dedx += ComputeNuclearDEDX(kinEnergy,part,mat); }
531 return dedx;
532}
533
534//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
535
537 const G4ParticleDefinition* p,
538 const G4Material* mat)
539{
540 G4double res = 0.0;
541 G4VEmProcess* nucst = FindDiscreteProcess(p, "nuclearStopping");
542 if(nucst) {
543 G4VEmModel* mod = nucst->EmModel();
544 if(mod) {
545 mod->SetFluctuationFlag(false);
546 res = mod->ComputeDEDXPerVolume(mat, p, kinEnergy);
547 }
548 }
549
550 if(verbose > 1) {
551 G4cout << p->GetParticleName() << " E(MeV)= " << kinEnergy/MeV
552 << " NuclearDEDX(MeV/mm)= " << res*mm/MeV
553 << " NuclearDEDX(MeV*cm^2/g)= "
554 << res*gram/(MeV*cm2*mat->GetDensity())
555 << G4endl;
556 }
557 return res;
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
563 G4double kinEnergy,
564 const G4ParticleDefinition* p,
565 const G4String& processName,
566 const G4Material* mat,
567 G4double cut)
568{
569 SetupMaterial(mat);
570 G4double res = 0.0;
571 if(UpdateParticle(p, kinEnergy)) {
572 if(FindEmModel(p, processName, kinEnergy)) {
573 G4double e = kinEnergy;
574 G4double aCut = std::max(cut, theParameters->LowestElectronEnergy());
575 chargeSquare = currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
576 if (baseParticle && !isIon) {
577 e *= kinEnergy*massRatio;
578 res = currentModel->CrossSectionPerVolume(
579 mat, baseParticle, e, aCut, e) * chargeSquare;
580 } else {
581 res = currentModel->CrossSectionPerVolume(mat, p, e, aCut, e);
582 }
583 if(verbose>0) {
584 G4cout << "G4EmCalculator::ComputeXSPerVolume: E(MeV)= "
585 << kinEnergy/MeV
586 << " cross(cm-1)= " << res*cm
587 << " cut(keV)= " << aCut/keV
588 << " " << p->GetParticleName()
589 << " in " << mat->GetName()
590 << G4endl;
591 }
592 }
593 }
594 return res;
595}
596
597//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
598
601 const G4ParticleDefinition* p,
602 const G4String& processName,
604 G4double cut)
605{
606 G4double res = 0.0;
607 if(UpdateParticle(p, kinEnergy)) {
608 G4int iz = G4lrint(Z);
609 CheckMaterial(iz);
610 if(FindEmModel(p, processName, kinEnergy)) {
611 G4double e = kinEnergy;
612 G4double aCut = std::max(cut, theParameters->LowestElectronEnergy());
613 if(baseParticle) {
614 e *= kinEnergy*massRatio;
615 currentModel->InitialiseForElement(baseParticle, iz);
616 res = currentModel->ComputeCrossSectionPerAtom(
617 baseParticle, e, Z, A, aCut) * chargeSquare;
618 } else {
619 currentModel->InitialiseForElement(p, iz);
620 res = currentModel->ComputeCrossSectionPerAtom(p, e, Z, A, aCut);
621 }
622 if(verbose>0) {
623 G4cout << "E(MeV)= " << kinEnergy/MeV
624 << " cross(barn)= " << res/barn
625 << " " << p->GetParticleName()
626 << " Z= " << Z << " A= " << A/(g/mole) << " g/mole"
627 << " cut(keV)= " << aCut/keV
628 << G4endl;
629 }
630 }
631 }
632 return res;
633}
634
635//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
636
639 const G4ParticleDefinition* p,
640 const G4String& processName,
641 G4int Z, G4int shellIdx,
642 G4double cut)
643{
644 G4double res = 0.0;
645 if(UpdateParticle(p, kinEnergy)) {
646 CheckMaterial(Z);
647 if(FindEmModel(p, processName, kinEnergy)) {
648 G4double e = kinEnergy;
649 G4double aCut = std::max(cut, theParameters->LowestElectronEnergy());
650 if(nullptr != baseParticle) {
651 e *= kinEnergy*massRatio;
652 currentModel->InitialiseForElement(baseParticle, Z);
653 res =
654 currentModel->ComputeCrossSectionPerShell(baseParticle, Z, shellIdx,
655 e, aCut) * chargeSquare;
656 } else {
657 currentModel->InitialiseForElement(p, Z);
658 res = currentModel->ComputeCrossSectionPerAtom(p, Z, shellIdx, e, aCut);
659 }
660 if(verbose>0) {
661 G4cout << "E(MeV)= " << kinEnergy/MeV
662 << " cross(barn)= " << res/barn
663 << " " << p->GetParticleName()
664 << " Z= " << Z << " shellIdx= " << shellIdx
665 << " cut(keV)= " << aCut/keV
666 << G4endl;
667 }
668 }
669 }
670 return res;
671}
672
673//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
674
677 const G4Material* mat)
678{
679 G4double res = 0.0;
680 const G4ParticleDefinition* gamma = G4Gamma::Gamma();
681 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "conv", mat, 0.0);
682 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "compt", mat, 0.0);
683 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "phot", mat, 0.0);
684 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "Rayl", mat, 0.0);
685 if(res > 0.0) { res = 1.0/res; }
686 return res;
687}
688
689//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
690
692 const G4String& particle,
693 G4int Z,
695 G4double kinEnergy,
696 const G4Material* mat)
697{
698 G4double res = 0.0;
699 const G4ParticleDefinition* p = FindParticle(particle);
700 G4VAtomDeexcitation* ad = manager->AtomDeexcitation();
701 if(p && ad) {
702 res = ad->ComputeShellIonisationCrossSectionPerAtom(p, Z, shell,
703 kinEnergy, mat);
704 }
705 return res;
706}
707
708//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
709
711 const G4ParticleDefinition* p,
712 const G4String& processName,
713 const G4Material* mat,
714 G4double cut)
715{
716 G4double mfp = DBL_MAX;
717 G4double x =
718 ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, cut);
719 if(x > 0.0) { mfp = 1.0/x; }
720 if(verbose>1) {
721 G4cout << "E(MeV)= " << kinEnergy/MeV
722 << " MFP(mm)= " << mfp/mm
723 << " " << p->GetParticleName()
724 << " in " << mat->GetName()
725 << G4endl;
726 }
727 return mfp;
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
733 G4double range,
734 const G4ParticleDefinition* part,
735 const G4Material* mat)
736{
738 ConvertRangeToEnergy(part, mat, range);
739}
740
741//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742
743G4bool G4EmCalculator::UpdateParticle(const G4ParticleDefinition* p,
744 G4double kinEnergy)
745{
746 if(p != currentParticle) {
747
748 // new particle
749 currentParticle = p;
750 dynParticle->SetDefinition(const_cast<G4ParticleDefinition*>(p));
751 dynParticle->SetKineticEnergy(kinEnergy);
752 baseParticle = nullptr;
753 currentParticleName = p->GetParticleName();
754 massRatio = 1.0;
755 mass = p->GetPDGMass();
756 chargeSquare = 1.0;
757 currentProcess = manager->GetEnergyLossProcess(p);
758 currentProcessName = "";
759 isIon = false;
760
761 // ionisation process exist
762 if(nullptr != currentProcess) {
763 currentProcessName = currentProcess->GetProcessName();
764 baseParticle = currentProcess->BaseParticle();
765 if (currentProcessName == "ionIoni" && !(p->GetParticleName() == "alpha")) {
766 baseParticle = theGenericIon;
767 isIon = true;
768 }
769
770 // base particle is used
771 if(nullptr != baseParticle) {
772 massRatio = baseParticle->GetPDGMass()/p->GetPDGMass();
773 G4double q = p->GetPDGCharge()/baseParticle->GetPDGCharge();
774 chargeSquare = q*q;
775 }
776 }
777 }
778 // Effective charge for ions
779 if(isIon && nullptr != currentProcess) {
780 chargeSquare =
781 corr->EffectiveChargeSquareRatio(p, currentMaterial, kinEnergy);
782 currentProcess->SetDynamicMassCharge(massRatio,chargeSquare);
783 if(verbose>1) {
784 G4cout <<"\n NewIon: massR= "<< massRatio << " q2= "
785 << chargeSquare << " " << currentProcess << G4endl;
786 }
787 }
788 return true;
789}
790
791//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
792
794{
795 const G4ParticleDefinition* p = nullptr;
796 if(name != currentParticleName) {
798 if(nullptr == p) {
799 G4cout << "### WARNING: G4EmCalculator::FindParticle fails to find "
800 << name << G4endl;
801 }
802 } else {
803 p = currentParticle;
804 }
805 return p;
806}
807
808//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
809
811{
812 const G4ParticleDefinition* p = ionTable->GetIon(Z,A,0);
813 return p;
814}
815
816//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
817
819{
820 if(name != currentMaterialName) {
822 if(nullptr == currentMaterial) {
823 G4cout << "### WARNING: G4EmCalculator::FindMaterial fails to find "
824 << name << G4endl;
825 }
826 }
827 return currentMaterial;
828}
829
830//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
831
833{
834 return G4EmUtility::FindRegion(reg);
835}
836
837//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
838
840 const G4Material* material,
841 const G4Region* region)
842{
843 const G4MaterialCutsCouple* couple = nullptr;
844 SetupMaterial(material);
845 if(nullptr != currentMaterial) {
846 // Access to materials
847 const G4ProductionCutsTable* theCoupleTable=
849 const G4Region* r = region;
850 if(nullptr != r) {
851 couple = theCoupleTable->GetMaterialCutsCouple(material,
852 r->GetProductionCuts());
853 } else {
855 std::size_t nr = store->size();
856 if(0 < nr) {
857 for(std::size_t i=0; i<nr; ++i) {
858 couple = theCoupleTable->GetMaterialCutsCouple(
859 material, ((*store)[i])->GetProductionCuts());
860 if(nullptr != couple) { break; }
861 }
862 }
863 }
864 }
865 if(nullptr == couple) {
867 ed << "G4EmCalculator::FindCouple: fail for material <"
868 << currentMaterialName << ">";
869 if(region) { ed << " and region " << region->GetName(); }
870 G4Exception("G4EmCalculator::FindCouple", "em0078",
871 FatalException, ed);
872 }
873 return couple;
874}
875
876//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
877
878G4bool G4EmCalculator::UpdateCouple(const G4Material* material, G4double cut)
879{
880 SetupMaterial(material);
881 if (nullptr == currentMaterial) { return false; }
882 for (G4int i=0; i<nLocalMaterials; ++i) {
883 G4cout << material->GetName()
884 << " " << cut << " " << localCuts[i] << G4endl;
885 if(material == localMaterials[i] && cut == localCuts[i]) {
886 currentCouple = localCouples[i];
887 currentCoupleIndex = currentCouple->GetIndex();
888 currentCut = cut;
889 return true;
890 }
891 }
892 auto prodcut = new G4ProductionCuts();
893 prodcut->SetProductionCut(cut, 0);
894 prodcut->SetProductionCut(cut, 1);
895 auto const cc = new G4MaterialCutsCouple(material, prodcut);
896 localMaterials.push_back(material);
897 localCouples.push_back(cc);
898 localCuts.push_back(cut);
899 ++nLocalMaterials;
900 currentCouple = cc;
901 currentCoupleIndex = currentCouple->GetIndex();
902 currentCut = cut;
903 return true;
904}
905
906//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
907
908void G4EmCalculator::FindLambdaTable(const G4ParticleDefinition* p,
909 const G4String& processName,
910 G4double kinEnergy, G4int& proctype)
911{
912 // Search for the process
913 if (!currentLambda || p != lambdaParticle || processName != lambdaName) {
914 lambdaName = processName;
915 currentLambda = nullptr;
916 lambdaParticle = p;
917 isApplicable = false;
918
919 const G4ParticleDefinition* part = (isIon) ? theGenericIon : p;
920
921 // Search for energy loss process
922 currentName = processName;
923 currentModel = nullptr;
924 loweModel = nullptr;
925
926 G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
927 if(nullptr != elproc) {
928 currentLambda = elproc->LambdaTable();
929 proctype = 0;
930 if(nullptr != currentLambda) {
931 isApplicable = true;
932 if(verbose>1) {
933 G4cout << "G4VEnergyLossProcess is found out: " << currentName
934 << G4endl;
935 }
936 }
937 curProcess = elproc;
938 return;
939 }
940
941 // Search for discrete process
942 G4VEmProcess* proc = FindDiscreteProcess(part, processName);
943 if(nullptr != proc) {
944 currentLambda = proc->LambdaTable();
945 proctype = 1;
946 if(nullptr != currentLambda) {
947 isApplicable = true;
948 if(verbose>1) {
949 G4cout << "G4VEmProcess is found out: " << currentName << G4endl;
950 }
951 }
952 curProcess = proc;
953 return;
954 }
955
956 // Search for msc process
957 G4VMultipleScattering* msc = FindMscProcess(part, processName);
958 if(nullptr != msc) {
959 currentModel = msc->SelectModel(kinEnergy,0);
960 proctype = 2;
961 if(nullptr != currentModel) {
962 currentLambda = currentModel->GetCrossSectionTable();
963 if(nullptr != currentLambda) {
964 isApplicable = true;
965 if(verbose>1) {
966 G4cout << "G4VMultipleScattering is found out: " << currentName
967 << G4endl;
968 }
969 }
970 }
971 curProcess = msc;
972 }
973 }
974}
975
976//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
977
978G4bool G4EmCalculator::FindEmModel(const G4ParticleDefinition* p,
979 const G4String& processName,
980 G4double kinEnergy)
981{
982 isApplicable = false;
983 if(nullptr == p || nullptr == currentMaterial) {
984 G4cout << "G4EmCalculator::FindEmModel WARNING: no particle"
985 << " or materail defined; particle: " << p << G4endl;
986 return isApplicable;
987 }
988 G4String partname = p->GetParticleName();
989 G4double scaledEnergy = kinEnergy*massRatio;
990 const G4ParticleDefinition* part = (isIon) ? theGenericIon : p;
991
992 if(verbose > 1) {
993 G4cout << "## G4EmCalculator::FindEmModel for " << partname
994 << " (type= " << p->GetParticleType()
995 << ") and " << processName << " at E(MeV)= " << scaledEnergy
996 << G4endl;
997 if(p != part) { G4cout << " GenericIon is the base particle" << G4endl; }
998 }
999
1000 // Search for energy loss process
1001 currentName = processName;
1002 currentModel = nullptr;
1003 loweModel = nullptr;
1004 std::size_t idx = 0;
1005
1006 G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
1007 if(nullptr != elproc) {
1008 currentModel = elproc->SelectModelForMaterial(scaledEnergy, idx);
1009 currentModel->InitialiseForMaterial(part, currentMaterial);
1010 currentModel->SetupForMaterial(part, currentMaterial, kinEnergy);
1011 G4double eth = currentModel->LowEnergyLimit();
1012 if(eth > 0.0) {
1013 loweModel = elproc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1014 if(loweModel == currentModel) { loweModel = nullptr; }
1015 else {
1016 loweModel->InitialiseForMaterial(part, currentMaterial);
1017 loweModel->SetupForMaterial(part, currentMaterial, eth - CLHEP::eV);
1018 }
1019 }
1020 }
1021
1022 // Search for discrete process
1023 if(nullptr == currentModel) {
1024 G4VEmProcess* proc = FindDiscreteProcess(part, processName);
1025 if(nullptr != proc) {
1026 currentModel = proc->SelectModelForMaterial(kinEnergy, idx);
1027 currentModel->InitialiseForMaterial(part, currentMaterial);
1028 currentModel->SetupForMaterial(part, currentMaterial, kinEnergy);
1029 G4double eth = currentModel->LowEnergyLimit();
1030 if(eth > 0.0) {
1031 loweModel = proc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1032 if(loweModel == currentModel) { loweModel = nullptr; }
1033 else {
1034 loweModel->InitialiseForMaterial(part, currentMaterial);
1035 loweModel->SetupForMaterial(part, currentMaterial, eth - CLHEP::eV);
1036 }
1037 }
1038 }
1039 }
1040
1041 // Search for msc process
1042 if(nullptr == currentModel) {
1043 G4VMultipleScattering* proc = FindMscProcess(part, processName);
1044 if(nullptr != proc) {
1045 currentModel = proc->SelectModel(kinEnergy, idx);
1046 loweModel = nullptr;
1047 }
1048 }
1049 if(nullptr != currentModel) {
1050 if(loweModel == currentModel) { loweModel = nullptr; }
1051 isApplicable = true;
1052 currentModel->InitialiseForMaterial(part, currentMaterial);
1053 if(loweModel) {
1054 loweModel->InitialiseForMaterial(part, currentMaterial);
1055 }
1056 if(verbose > 1) {
1057 G4cout << " Model <" << currentModel->GetName()
1058 << "> Emin(MeV)= " << currentModel->LowEnergyLimit()/MeV
1059 << " for " << part->GetParticleName();
1060 if(nullptr != elproc) {
1061 G4cout << " and " << elproc->GetProcessName() << " " << elproc
1062 << G4endl;
1063 }
1064 if(nullptr != loweModel) {
1065 G4cout << " LowEnergy model <" << loweModel->GetName() << ">";
1066 }
1067 G4cout << G4endl;
1068 }
1069 }
1070 return isApplicable;
1071}
1072
1073//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1074
1076G4EmCalculator::FindEnLossProcess(const G4ParticleDefinition* part,
1077 const G4String& processName)
1078{
1079 G4VEnergyLossProcess* proc = nullptr;
1080 const std::vector<G4VEnergyLossProcess*> v =
1081 manager->GetEnergyLossProcessVector();
1082 std::size_t n = v.size();
1083 for(std::size_t i=0; i<n; ++i) {
1084 if((v[i])->GetProcessName() == processName) {
1085 auto p = static_cast<G4VProcess*>(v[i]);
1086 if(ActiveForParticle(part, p)) {
1087 proc = v[i];
1088 break;
1089 }
1090 }
1091 }
1092 return proc;
1093}
1094
1095//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1096
1098G4EmCalculator::FindDiscreteProcess(const G4ParticleDefinition* part,
1099 const G4String& processName)
1100{
1101 G4VEmProcess* proc = nullptr;
1102 auto v = manager->GetEmProcessVector();
1103 std::size_t n = v.size();
1104 for(std::size_t i=0; i<n; ++i) {
1105 const G4String& pName = v[i]->GetProcessName();
1106 if(pName == "GammaGeneralProc") {
1107 proc = v[i]->GetEmProcess(processName);
1108 break;
1109 } else if(pName == processName) {
1110 const auto p = static_cast<G4VProcess*>(v[i]);
1111 if(ActiveForParticle(part, p)) {
1112 proc = v[i];
1113 break;
1114 }
1115 }
1116 }
1117 return proc;
1118}
1119
1120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1121
1123G4EmCalculator::FindMscProcess(const G4ParticleDefinition* part,
1124 const G4String& processName)
1125{
1126 G4VMultipleScattering* proc = nullptr;
1127 const std::vector<G4VMultipleScattering*> v =
1128 manager->GetMultipleScatteringVector();
1129 std::size_t n = v.size();
1130 for(std::size_t i=0; i<n; ++i) {
1131 if((v[i])->GetProcessName() == processName) {
1132 auto p = static_cast<G4VProcess*>(v[i]);
1133 if(ActiveForParticle(part, p)) {
1134 proc = v[i];
1135 break;
1136 }
1137 }
1138 }
1139 return proc;
1140}
1141
1142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1143
1145 const G4String& processName)
1146{
1147 G4VProcess* proc = nullptr;
1148 const G4ProcessManager* procman = part->GetProcessManager();
1149 G4ProcessVector* pv = procman->GetProcessList();
1150 G4int nproc = (G4int)pv->size();
1151 for(G4int i=0; i<nproc; ++i) {
1152 if(processName == (*pv)[i]->GetProcessName()) {
1153 proc = (*pv)[i];
1154 break;
1155 }
1156 }
1157 return proc;
1158}
1159
1160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1161
1162G4bool G4EmCalculator::ActiveForParticle(const G4ParticleDefinition* part,
1163 G4VProcess* proc)
1164{
1165 G4ProcessManager* pm = part->GetProcessManager();
1166 G4ProcessVector* pv = pm->GetProcessList();
1167 G4int n = (G4int)pv->size();
1168 G4bool res = false;
1169 for(G4int i=0; i<n; ++i) {
1170 if((*pv)[i] == proc) {
1171 if(pm->GetProcessActivation(i)) { res = true; }
1172 break;
1173 }
1174 }
1175 return res;
1176}
1177
1178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1179
1181{
1182 if (nullptr != mat) {
1183 currentMaterial = mat;
1184 currentMaterialName = mat->GetName();
1185 } else {
1186 currentMaterial = nullptr;
1187 currentMaterialName = "";
1188 }
1189}
1190
1191//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1192
1194{
1195 SetupMaterial(nist->FindOrBuildMaterial(mname));
1196}
1197
1198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1199
1200void G4EmCalculator::CheckMaterial(G4int Z)
1201{
1202 G4bool isFound = false;
1203 if(nullptr != currentMaterial) {
1204 G4int nn = (G4int)currentMaterial->GetNumberOfElements();
1205 for(G4int i=0; i<nn; ++i) {
1206 if(Z == currentMaterial->GetElement(i)->GetZasInt()) {
1207 isFound = true;
1208 break;
1209 }
1210 }
1211 }
1212 if(!isFound) {
1213 SetupMaterial(nist->FindOrBuildSimpleMaterial(Z));
1214 }
1215}
1216
1217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1218
1220{
1221 verbose = verb;
1222}
1223
1224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1225
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Electron()
Definition G4Electron.cc:91
G4int GetZasInt() const
Definition G4Element.hh:120
const G4ParticleDefinition * FindParticle(const G4String &)
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
void SetVerbose(G4int val)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double ComputeCrossSectionPerShell(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4int Z, G4int shellIdx, G4double cut=0.0)
G4double ComputeShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy, const G4Material *mat=nullptr)
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=nullptr)
void SetupMaterial(const G4Material *)
void PrintRangeTable(const G4ParticleDefinition *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double rangecut=DBL_MAX)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
const G4Material * FindMaterial(const G4String &)
static G4EmParameters * Instance()
static const G4Region * FindRegion(const G4String &regionName, const G4int verbose=0)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4double GetDensity() const
const G4Element * GetElement(G4int iel) const
std::size_t GetNumberOfElements() const
const G4String & GetName() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
static G4NistManager * Instance()
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4Positron * Positron()
Definition G4Positron.cc:90
G4bool GetProcessActivation(G4VProcess *aProcess) const
G4ProcessVector * GetProcessList() const
std::size_t size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4RegionStore is a singleton class, acting as container for all geometrical regions,...
static G4RegionStore * GetInstance()
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
virtual G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
void SetCurrentCouple(const G4MaterialCutsCouple *)
void SetFluctuationFlag(G4bool val)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
G4VEmModel * EmModel(std::size_t index=0) const
G4PhysicsTable * LambdaTable() const
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, std::size_t idxCouple) const
G4double GetCrossSection(const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
const G4ParticleDefinition * BaseParticle() const
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, std::size_t &idxCouple) const
G4PhysicsTable * LambdaTable() const
G4PhysicsTable * DEDXTable() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
const G4String & GetProcessName() const
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62