Geant4-11
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
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
76{
81 cutenergy[0] = cutenergy[1] = cutenergy[2] = DBL_MAX;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
91{
92 delete ionEffCharge;
93 delete dynParticle;
94 for (G4int i=0; i<nLocalMaterials; ++i) {
95 delete localCouples[i];
96 }
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
102 const G4ParticleDefinition* p,
103 const G4Material* mat,
104 const G4Region* region)
105{
106 G4double res = 0.0;
107 const G4MaterialCutsCouple* couple = FindCouple(mat, region);
108 if(couple && UpdateParticle(p, kinEnergy) ) {
109 res = manager->GetDEDX(p, kinEnergy, couple);
110
111 if(isIon) {
112 if(FindEmModel(p, currentProcessName, kinEnergy)) {
113 G4double length = CLHEP::nm;
114 G4double eloss = res*length;
115 //G4cout << "### GetDEDX: E= " << kinEnergy << " dedx0= " << res
116 // << " de= " << eloss << G4endl;;
117 dynParticle->SetKineticEnergy(kinEnergy);
118 currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
119 currentModel->CorrectionsAlongStep(couple,dynParticle,length,eloss);
120 res = eloss/length;
121 //G4cout << " de1= " << eloss << " res1= " << res
122 // << " " << p->GetParticleName() <<G4endl;;
123 }
124 }
125
126 if(verbose>0) {
127 G4cout << "G4EmCalculator::GetDEDX: E(MeV)= " << kinEnergy/MeV
128 << " DEDX(MeV/mm)= " << res*mm/MeV
129 << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
130 << " " << p->GetParticleName()
131 << " in " << mat->GetName()
132 << " isIon= " << isIon
133 << G4endl;
134 }
135 }
136 return res;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
142 const G4ParticleDefinition* p,
143 const G4Material* mat,
144 const G4Region* region)
145{
146 G4double res = 0.0;
147 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
148 if(couple && UpdateParticle(p, kinEnergy)) {
149 res = manager->GetRangeFromRestricteDEDX(p, kinEnergy, couple);
150 if(verbose>1) {
151 G4cout << " G4EmCalculator::GetRangeFromRestrictedDEDX: E(MeV)= "
152 << kinEnergy/MeV
153 << " range(mm)= " << res/mm
154 << " " << p->GetParticleName()
155 << " in " << mat->GetName()
156 << G4endl;
157 }
158 }
159 return res;
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163
165 const G4ParticleDefinition* p,
166 const G4Material* mat,
167 const G4Region* region)
168{
169 G4double res = 0.0;
172 ed << "G4EmCalculator::GetCSDARange: CSDA table is not built; "
173 << " use UI command: /process/eLoss/CSDARange true";
174 G4Exception("G4EmCalculator::GetCSDARange", "em0077",
175 JustWarning, ed);
176 return res;
177 }
178
179 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
180 if(couple && UpdateParticle(p, kinEnergy)) {
181 res = manager->GetCSDARange(p, kinEnergy, couple);
182 if(verbose>1) {
183 G4cout << " G4EmCalculator::GetCSDARange: E(MeV)= " << kinEnergy/MeV
184 << " range(mm)= " << res/mm
185 << " " << p->GetParticleName()
186 << " in " << mat->GetName()
187 << G4endl;
188 }
189 }
190 return res;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
196 const G4ParticleDefinition* p,
197 const G4Material* mat,
198 const G4Region* region)
199{
200 G4double res = 0.0;
202 res = GetCSDARange(kinEnergy, p, mat, region);
203 } else {
204 res = GetRangeFromRestricteDEDX(kinEnergy, p, mat, region);
205 }
206 return res;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210
212 const G4ParticleDefinition* p,
213 const G4Material* mat,
214 const G4Region* region)
215{
216 G4double res = 0.0;
217 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
218 if(couple && UpdateParticle(p, 1.0*GeV)) {
219 res = manager->GetEnergy(p, range, couple);
220 if(verbose>0) {
221 G4cout << "G4EmCalculator::GetKinEnergy: Range(mm)= " << range/mm
222 << " KinE(MeV)= " << res/MeV
223 << " " << p->GetParticleName()
224 << " in " << mat->GetName()
225 << G4endl;
226 }
227 }
228 return res;
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232
234 const G4ParticleDefinition* p,
235 const G4String& processName,
236 const G4Material* mat,
237 const G4Region* region)
238{
239 G4double res = 0.0;
240 const G4MaterialCutsCouple* couple = FindCouple(mat,region);
241
242 if(couple && UpdateParticle(p, kinEnergy)) {
243 if(FindEmModel(p, processName, kinEnergy)) {
244 G4int idx = couple->GetIndex();
245 G4int procType = -1;
246 FindLambdaTable(p, processName, kinEnergy, procType);
247
248 G4VEmProcess* emproc = FindDiscreteProcess(p, processName);
249 if(emproc) {
250 res = emproc->CrossSectionPerVolume(kinEnergy, couple);
251 } else if(currentLambda) {
252 // special tables are built for Msc models (procType is set in FindLambdaTable
253 if(procType==2) {
254 G4VMscModel* mscM = static_cast<G4VMscModel*>(currentModel);
255 mscM->SetCurrentCouple(couple);
256 G4double tr1Mfp = mscM->GetTransportMeanFreePath(p, kinEnergy);
257 if (tr1Mfp<DBL_MAX) {
258 res = 1./tr1Mfp;
259 }
260 } else {
261 G4double e = kinEnergy*massRatio;
262 res = (((*currentLambda)[idx])->Value(e))*chargeSquare;
263 }
264 } else {
265 res = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, kinEnergy);
266 }
267 if(verbose>0) {
268 G4cout << "G4EmCalculator::GetXSPerVolume: E(MeV)= " << kinEnergy/MeV
269 << " cross(cm-1)= " << res*cm
270 << " " << p->GetParticleName()
271 << " in " << mat->GetName();
272 if(verbose>1)
273 G4cout << " idx= " << idx << " Escaled((MeV)= "
274 << kinEnergy*massRatio
275 << " q2= " << chargeSquare;
276 G4cout << G4endl;
277 }
278 }
279 }
280 return res;
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
284
286 const G4String& particle,
287 G4int Z,
289 G4double kinEnergy)
290{
291 G4double res = 0.0;
292 const G4ParticleDefinition* p = FindParticle(particle);
294 if(p && ad) {
295 res = ad->GetShellIonisationCrossSectionPerAtom(p, Z, shell, kinEnergy);
296 }
297 return res;
298}
299
300//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
301
303 const G4ParticleDefinition* p,
304 const G4String& processName,
305 const G4Material* mat,
306 const G4Region* region)
307{
308 G4double res = DBL_MAX;
309 G4double x = GetCrossSectionPerVolume(kinEnergy,p, processName, mat,region);
310 if(x > 0.0) { res = 1.0/x; }
311 if(verbose>1) {
312 G4cout << "G4EmCalculator::GetMeanFreePath: E(MeV)= " << kinEnergy/MeV
313 << " MFP(mm)= " << res/mm
314 << " " << p->GetParticleName()
315 << " in " << mat->GetName()
316 << G4endl;
317 }
318 return res;
319}
320
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
322
324{
326 G4cout << "##### DEDX Table for " << p->GetParticleName() << G4endl;
327 if(elp) G4cout << *(elp->DEDXTable()) << G4endl;
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
331
333{
335 G4cout << "##### Range Table for " << p->GetParticleName() << G4endl;
336 if(elp) G4cout << *(elp->RangeTableForLoss()) << G4endl;
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340
342{
344 G4cout << "### G4EmCalculator: Inverse Range Table for "
345 << p->GetParticleName() << G4endl;
346 if(elp) G4cout << *(elp->InverseRangeTable()) << G4endl;
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350
352 const G4ParticleDefinition* p,
353 const G4String& processName,
354 const G4Material* mat,
355 G4double cut)
356{
357 SetupMaterial(mat);
358 G4double res = 0.0;
359 if(verbose > 1) {
360 G4cout << "### G4EmCalculator::ComputeDEDX: " << p->GetParticleName()
361 << " in " << currentMaterialName
362 << " e(MeV)= " << kinEnergy/MeV << " cut(MeV)= " << cut/MeV
363 << G4endl;
364 }
365 if(UpdateParticle(p, kinEnergy)) {
366 if(FindEmModel(p, processName, kinEnergy)) {
367
368 // Special case of ICRU'73 model
369 const G4String& mname = currentModel->GetName();
370 if(mname == "ParamICRU73" || mname == "LinhardSorensen" || mname == "Atima") {
371 res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
372 if(verbose > 1) {
373 G4cout << mname << " ion E(MeV)= " << kinEnergy << " ";
374 G4cout << currentModel->GetName() << ": DEDX(MeV/mm)= " << res*mm/MeV
375 << " DEDX(MeV*cm^2/g)= "
376 << res*gram/(MeV*cm2*mat->GetDensity())
377 << G4endl;
378 }
379 } else {
380
381 G4double escaled = kinEnergy*massRatio;
382 if(baseParticle) {
384 mat, baseParticle, escaled, cut) * chargeSquare;
385 if(verbose > 1) {
387 << " Escaled(MeV)= " << escaled;
388 }
389 } else {
390 res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
391 if(verbose > 1) { G4cout << " no basePart E(MeV)= " << kinEnergy << " "; }
392 }
393 if(verbose > 1) {
394 G4cout << currentModel->GetName() << ": DEDX(MeV/mm)= " << res*mm/MeV
395 << " DEDX(MeV*cm^2/g)= "
396 << res*gram/(MeV*cm2*mat->GetDensity())
397 << G4endl;
398 }
399
400 // emulate smoothing procedure
402 // G4cout << "massRatio= " << massRatio << " eth= " << eth << G4endl;
403 if(loweModel) {
404 G4double res0 = 0.0;
405 G4double res1 = 0.0;
406 if(baseParticle) {
407 res1 = currentModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
408 * chargeSquare;
409 res0 = loweModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
410 * chargeSquare;
411 } else {
412 res1 = currentModel->ComputeDEDXPerVolume(mat, p, eth, cut);
413 res0 = loweModel->ComputeDEDXPerVolume(mat, p, eth, cut);
414 }
415 if(verbose > 1) {
416 G4cout << "At boundary energy(MeV)= " << eth/MeV
417 << " DEDX(MeV/mm)= " << res1*mm/MeV
418 << G4endl;
419 }
420
421 //G4cout << "eth= " << eth << " escaled= " << escaled
422 // << " res0= " << res0 << " res1= "
423 // << res1 << " q2= " << chargeSquare << G4endl;
424
425 if(res1 > 0.0 && escaled > 0.0) {
426 res *= (1.0 + (res0/res1 - 1.0)*eth/escaled);
427 }
428 }
429
430 // low energy correction for ions
431 if(isIon) {
432 G4double length = CLHEP::nm;
433 const G4Region* r = 0;
434 const G4MaterialCutsCouple* couple = FindCouple(mat, r);
435 G4double eloss = res*length;
436 dynParticle->SetKineticEnergy(kinEnergy);
437 currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
438 currentModel->CorrectionsAlongStep(couple,dynParticle,length,eloss);
439 res = eloss/length;
440
441 if(verbose > 1) {
442 G4cout << "After Corrections: DEDX(MeV/mm)= " << res*mm/MeV
443 << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
444 << G4endl;
445 }
446 }
447 }
448 }
449 if(verbose > 0) {
450 G4cout << "Sum: E(MeV)= " << kinEnergy/MeV
451 << " DEDX(MeV/mm)= " << res*mm/MeV
452 << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
453 << " cut(MeV)= " << cut/MeV
454 << " " << p->GetParticleName()
455 << " in " << currentMaterialName
456 << " Zi^2= " << chargeSquare
457 << " isIon=" << isIon
458 << G4endl;
459 }
460 }
461 return res;
462}
463
464//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
465
467 const G4ParticleDefinition* part,
468 const G4Material* mat,
469 G4double cut)
470{
471 SetupMaterial(mat);
472 G4double dedx = 0.0;
473 if(UpdateParticle(part, kinEnergy)) {
474
476 const std::vector<G4VEnergyLossProcess*> vel =
477 lManager->GetEnergyLossProcessVector();
478 G4int n = vel.size();
479
480 //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
481 // << " n= " << n << G4endl;
482
483 for(G4int i=0; i<n; ++i) {
484 if(vel[i]) {
485 G4VProcess* p = reinterpret_cast<G4VProcess*>(vel[i]);
486 if(ActiveForParticle(part, p)) {
487 //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
488 // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
489 dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),mat,cut);
490 }
491 }
492 }
493 }
494 return dedx;
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
498
500 const G4ParticleDefinition* part,
501 const G4Material* mat,
502 G4double rangecut)
503{
504 SetupMaterial(mat);
505 G4double dedx = 0.0;
506 if(UpdateParticle(part, kinEnergy)) {
507
509 const std::vector<G4VEnergyLossProcess*> vel =
510 lManager->GetEnergyLossProcessVector();
511 G4int n = vel.size();
512
513 if(mat != cutMaterial) {
514 cutMaterial = mat;
518 }
519
520 //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
521 // << " n= " << n << G4endl;
522
523 for(G4int i=0; i<n; ++i) {
524 if(vel[i]) {
525 G4VProcess* p = reinterpret_cast<G4VProcess*>(vel[i]);
526 if(ActiveForParticle(part, p)) {
527 //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
528 // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
529 const G4ParticleDefinition* sec = (vel[i])->SecondaryParticle();
530 G4int idx = 0;
531 if(sec == G4Electron::Electron()) { idx = 1; }
532 else if(sec == G4Positron::Positron()) { idx = 2; }
533
534 dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),
535 mat,cutenergy[idx]);
536 }
537 }
538 }
539 }
540 return dedx;
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
544
546 const G4ParticleDefinition* part,
547 const G4Material* mat,
548 G4double cut)
549{
550 G4double dedx = ComputeElectronicDEDX(kinEnergy,part,mat,cut);
551 if(mass > 700.*MeV) { dedx += ComputeNuclearDEDX(kinEnergy,part,mat); }
552 return dedx;
553}
554
555//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
556
558 const G4ParticleDefinition* p,
559 const G4Material* mat)
560{
561 G4double res = 0.0;
562 G4VEmProcess* nucst = FindDiscreteProcess(p, "nuclearStopping");
563 if(nucst) {
564 G4VEmModel* mod = nucst->EmModel();
565 if(mod) {
566 mod->SetFluctuationFlag(false);
567 res = mod->ComputeDEDXPerVolume(mat, p, kinEnergy);
568 }
569 }
570
571 if(verbose > 1) {
572 G4cout << p->GetParticleName() << " E(MeV)= " << kinEnergy/MeV
573 << " NuclearDEDX(MeV/mm)= " << res*mm/MeV
574 << " NuclearDEDX(MeV*cm^2/g)= "
575 << res*gram/(MeV*cm2*mat->GetDensity())
576 << G4endl;
577 }
578 return res;
579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
582
584 G4double kinEnergy,
585 const G4ParticleDefinition* p,
586 const G4String& processName,
587 const G4Material* mat,
588 G4double cut)
589{
590 SetupMaterial(mat);
591 G4double res = 0.0;
592 if(UpdateParticle(p, kinEnergy)) {
593 if(FindEmModel(p, processName, kinEnergy)) {
594 G4double e = kinEnergy;
596 if(baseParticle) {
597 e *= kinEnergy*massRatio;
599 mat, baseParticle, e, aCut, e) * chargeSquare;
600 } else {
601 res = currentModel->CrossSectionPerVolume(mat, p, e, aCut, e);
602 }
603 if(verbose>0) {
604 G4cout << "G4EmCalculator::ComputeXSPerVolume: E(MeV)= " << kinEnergy/MeV
605 << " cross(cm-1)= " << res*cm
606 << " cut(keV)= " << aCut/keV
607 << " " << p->GetParticleName()
608 << " in " << mat->GetName()
609 << G4endl;
610 }
611 }
612 }
613 return res;
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
617
619 G4double kinEnergy,
620 const G4ParticleDefinition* p,
621 const G4String& processName,
623 G4double cut)
624{
625 G4double res = 0.0;
626 if(UpdateParticle(p, kinEnergy)) {
627 G4int iz = G4lrint(Z);
628 CheckMaterial(iz);
629 if(FindEmModel(p, processName, kinEnergy)) {
630 G4double e = kinEnergy;
632 if(baseParticle) {
633 e *= kinEnergy*massRatio;
636 baseParticle, e, Z, A, aCut) * chargeSquare;
637 } else {
639 res = currentModel->ComputeCrossSectionPerAtom(p, e, Z, A, aCut);
640 }
641 if(verbose>0) {
642 G4cout << "E(MeV)= " << kinEnergy/MeV
643 << " cross(barn)= " << res/barn
644 << " " << p->GetParticleName()
645 << " Z= " << Z << " A= " << A/(g/mole) << " g/mole"
646 << " cut(keV)= " << aCut/keV
647 << G4endl;
648 }
649 }
650 }
651 return res;
652}
653
654//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
655
657 const G4ParticleDefinition* p,
658 const G4String& processName,
659 G4int Z, G4int shellIdx,
660 G4double cut)
661{
662 G4double res = 0.0;
663 if(UpdateParticle(p, kinEnergy)) {
665 if(FindEmModel(p, processName, kinEnergy)) {
666 G4double e = kinEnergy;
668 if(baseParticle) {
669 e *= kinEnergy*massRatio;
672 e, aCut) * chargeSquare;
673 } else {
675 res = currentModel->ComputeCrossSectionPerAtom(p, Z, shellIdx, e, aCut);
676 }
677 if(verbose>0) {
678 G4cout << "E(MeV)= " << kinEnergy/MeV
679 << " cross(barn)= " << res/barn
680 << " " << p->GetParticleName()
681 << " Z= " << Z << " shellIdx= " << shellIdx
682 << " cut(keV)= " << aCut/keV
683 << G4endl;
684 }
685 }
686 }
687 return res;
688}
689
690//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691
694 const G4Material* mat)
695{
696 G4double res = 0.0;
697 const G4ParticleDefinition* gamma = G4Gamma::Gamma();
698 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "conv", mat, 0.0);
699 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "compt", mat, 0.0);
700 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "phot", mat, 0.0);
701 res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "Rayl", mat, 0.0);
702 if(res > 0.0) { res = 1.0/res; }
703 return res;
704}
705
706//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
707
709 const G4String& particle,
710 G4int Z,
712 G4double kinEnergy,
713 const G4Material* mat)
714{
715 G4double res = 0.0;
716 const G4ParticleDefinition* p = FindParticle(particle);
718 if(p && ad) {
719 res = ad->ComputeShellIonisationCrossSectionPerAtom(p, Z, shell,
720 kinEnergy, mat);
721 }
722 return res;
723}
724
725//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
726
728 const G4ParticleDefinition* p,
729 const G4String& processName,
730 const G4Material* mat,
731 G4double cut)
732{
733 G4double mfp = DBL_MAX;
734 G4double x = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, cut);
735 if(x > 0.0) { mfp = 1.0/x; }
736 if(verbose>1) {
737 G4cout << "E(MeV)= " << kinEnergy/MeV
738 << " MFP(mm)= " << mfp/mm
739 << " " << p->GetParticleName()
740 << " in " << mat->GetName()
741 << G4endl;
742 }
743 return mfp;
744}
745
746//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
747
749 G4double range,
750 const G4ParticleDefinition* part,
751 const G4Material* mat)
752{
754 ConvertRangeToEnergy(part, mat, range);
755}
756
757//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
758
760 G4double kinEnergy)
761{
762 if(p != currentParticle) {
763
764 // new particle
765 currentParticle = p;
767 dynParticle->SetKineticEnergy(kinEnergy);
768 baseParticle = 0;
770 massRatio = 1.0;
771 mass = p->GetPDGMass();
772 chargeSquare = 1.0;
775 isIon = false;
776
777 // ionisation process exist
778 if(currentProcess) {
781
782 // base particle is used
783 if(baseParticle) {
786 chargeSquare = q*q;
787 }
788
789 if(p->GetParticleType() == "nucleus"
790 && currentParticleName != "deuteron"
791 && currentParticleName != "triton"
792 && currentParticleName != "alpha+"
793 && currentParticleName != "alpha"
794 ) {
795 isIon = true;
798 if(verbose>1) {
799 G4cout << "\n G4EmCalculator::UpdateParticle: isIon 1 "
800 << p->GetParticleName()
801 << " in " << currentMaterial->GetName()
802 << " e= " << kinEnergy << G4endl;
803 }
804 }
805 }
806 }
807
808 // Effective charge for ions
809 if(isIon) {
813 if(currentProcess) {
815 if(verbose>1) {
816 G4cout <<"\n NewIon: massR= "<< massRatio << " q2= "
817 << chargeSquare << " " << currentProcess << G4endl;
818 }
819 }
820 }
821 return true;
822}
823
824//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
825
827{
828 const G4ParticleDefinition* p = 0;
831 if(!p) {
832 G4cout << "### WARNING: G4EmCalculator::FindParticle fails to find "
833 << name << G4endl;
834 }
835 } else {
836 p = currentParticle;
837 }
838 return p;
839}
840
841//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
842
844{
845 const G4ParticleDefinition* p = ionTable->GetIon(Z,A,0);
846 return p;
847}
848
849//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
850
852{
855 if(!currentMaterial) {
856 G4cout << "### WARNING: G4EmCalculator::FindMaterial fails to find "
857 << name << G4endl;
858 }
859 }
860 return currentMaterial;
861}
862
863//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
864
866{
867 const G4Region* r = 0;
868 if(reg != "" && reg != "world") {
870 } else {
871 r = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld");
872 }
873 return r;
874}
875
876//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
877
879 const G4Material* material,
880 const G4Region* region)
881{
882 const G4MaterialCutsCouple* couple = nullptr;
884 if(currentMaterial) {
885 // Access to materials
886 const G4ProductionCutsTable* theCoupleTable=
888 const G4Region* r = region;
889 if(r) {
890 couple = theCoupleTable->GetMaterialCutsCouple(material,
891 r->GetProductionCuts());
892 } else {
894 size_t nr = store->size();
895 if(0 < nr) {
896 for(size_t i=0; i<nr; ++i) {
897 couple = theCoupleTable->GetMaterialCutsCouple(
898 material, ((*store)[i])->GetProductionCuts());
899 if(couple) { break; }
900 }
901 }
902 }
903 }
904 if(!couple) {
906 ed << "G4EmCalculator::FindCouple: fail for material <"
907 << currentMaterialName << ">";
908 if(region) { ed << " and region " << region->GetName(); }
909 G4Exception("G4EmCalculator::FindCouple", "em0078",
910 FatalException, ed);
911 }
912 return couple;
913}
914
915//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
916
918{
920 if(!currentMaterial) { return false; }
921 for (G4int i=0; i<nLocalMaterials; ++i) {
922 if(material == localMaterials[i] && cut == localCuts[i]) {
925 currentCut = cut;
926 return true;
927 }
928 }
930 localMaterials.push_back(material);
931 localCouples.push_back(cc);
932 localCuts.push_back(cut);
934 currentCouple = cc;
936 currentCut = cut;
937 return true;
938}
939
940//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
941
943 const G4String& processName,
944 G4double kinEnergy, G4int& proctype)
945{
946 // Search for the process
947 if (!currentLambda || p != lambdaParticle || processName != lambdaName) {
948 lambdaName = processName;
949 currentLambda = nullptr;
950 lambdaParticle = p;
951
952 const G4ParticleDefinition* part = p;
953 if(isIon) { part = theGenericIon; }
954
955 // Search for energy loss process
956 currentName = processName;
957 currentModel = nullptr;
958 loweModel = nullptr;
959
960 G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
961 if(elproc) {
962 currentLambda = elproc->LambdaTable();
963 proctype = 0;
964 if(currentLambda) {
965 isApplicable = true;
966 if(verbose>1) {
967 G4cout << "G4VEnergyLossProcess is found out: " << currentName
968 << G4endl;
969 }
970 }
971 curProcess = elproc;
972 return;
973 }
974
975 // Search for discrete process
976 G4VEmProcess* proc = FindDiscreteProcess(part, processName);
977 if(proc) {
978 currentLambda = proc->LambdaTable();
979 proctype = 1;
980 if(currentLambda) {
981 isApplicable = true;
982 if(verbose>1) {
983 G4cout << "G4VEmProcess is found out: " << currentName << G4endl;
984 }
985 }
986 curProcess = proc;
987 return;
988 }
989
990 // Search for msc process
991 G4VMultipleScattering* msc = FindMscProcess(part, processName);
992 if(msc) {
993 currentModel = msc->SelectModel(kinEnergy,0);
994 proctype = 2;
995 if(currentModel) {
997 if(currentLambda) {
998 isApplicable = true;
999 if(verbose>1) {
1000 G4cout << "G4VMultipleScattering is found out: " << currentName
1001 << G4endl;
1002 }
1003 }
1004 }
1005 curProcess = msc;
1006 }
1007 }
1008}
1009
1010//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1011
1013 const G4String& processName,
1014 G4double kinEnergy)
1015{
1016 isApplicable = false;
1017 if(!p || !currentMaterial) {
1018 G4cout << "G4EmCalculator::FindEmModel WARNING: no particle"
1019 << " or materail defined; particle: " << p << G4endl;
1020 return isApplicable;
1021 }
1022 G4String partname = p->GetParticleName();
1023 const G4ParticleDefinition* part = p;
1024 G4double scaledEnergy = kinEnergy*massRatio;
1025 if(isIon) { part = theGenericIon; }
1026
1027 if(verbose > 1) {
1028 G4cout << "## G4EmCalculator::FindEmModel for " << partname
1029 << " (type= " << p->GetParticleType()
1030 << ") and " << processName << " at E(MeV)= " << scaledEnergy
1031 << G4endl;
1032 if(p != part) { G4cout << " GenericIon is the base particle" << G4endl; }
1033 }
1034
1035 // Search for energy loss process
1036 currentName = processName;
1037 currentModel = nullptr;
1038 loweModel = nullptr;
1039 size_t idx = 0;
1040
1041 G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
1042 if(elproc) {
1043 currentModel = elproc->SelectModelForMaterial(scaledEnergy, idx);
1045 currentModel->SetupForMaterial(part, currentMaterial, scaledEnergy);
1047 if(eth > 0.0) {
1048 loweModel = elproc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1049 if(loweModel == currentModel) { loweModel = nullptr; }
1050 else {
1053 }
1054 }
1055 }
1056
1057 // Search for discrete process
1058 if(!currentModel) {
1059 G4VEmProcess* proc = FindDiscreteProcess(part, processName);
1060 if(proc) {
1061 currentModel = proc->SelectModelForMaterial(kinEnergy, idx);
1065 if(eth > 0.0) {
1066 loweModel = proc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1067 if(loweModel == currentModel) { loweModel = nullptr; }
1068 else {
1071 }
1072 }
1073 }
1074 }
1075
1076 // Search for msc process
1077 if(!currentModel) {
1078 G4VMultipleScattering* proc = FindMscProcess(part, processName);
1079 if(proc) {
1080 currentModel = proc->SelectModel(kinEnergy, idx);
1081 loweModel = nullptr;
1082 }
1083 }
1084 if(currentModel) {
1085 if(loweModel == currentModel) { loweModel = nullptr; }
1086 isApplicable = true;
1088 if(loweModel) {
1090 }
1091 if(verbose > 1) {
1092 G4cout << " Model <" << currentModel->GetName()
1093 << "> Emin(MeV)= " << currentModel->LowEnergyLimit()/MeV
1094 << " for " << part->GetParticleName();
1095 if(elproc) {
1096 G4cout << " and " << elproc->GetProcessName() << " " << elproc
1097 << G4endl;
1098 }
1099 if(loweModel) {
1100 G4cout << " LowEnergy model <" << loweModel->GetName() << ">";
1101 }
1102 G4cout << G4endl;
1103 }
1104 }
1105 return isApplicable;
1106}
1107
1108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1109
1111 const G4ParticleDefinition* p)
1112{
1113 G4VEnergyLossProcess* elp = nullptr;
1114 G4String partname = p->GetParticleName();
1115 const G4ParticleDefinition* part = p;
1116
1117 if(p->GetParticleType() == "nucleus"
1118 && currentParticleName != "deuteron"
1119 && currentParticleName != "triton"
1120 && currentParticleName != "alpha"
1121 && currentParticleName != "alpha+"
1122 ) { part = theGenericIon; }
1123
1124 elp = manager->GetEnergyLossProcess(part);
1125 /*
1126 G4cout << "\n G4EmCalculator::FindEnergyLossProcess: for " << p->GetParticleName()
1127 << " found " << elp->GetProcessName() << " of "
1128 << elp->Particle()->GetParticleName() << " " << elp << G4endl;
1129 */
1130 return elp;
1131}
1132
1133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1134
1137 const G4String& processName)
1138{
1139 G4VEnergyLossProcess* proc = 0;
1140 const std::vector<G4VEnergyLossProcess*> v =
1142 G4int n = v.size();
1143 for(G4int i=0; i<n; ++i) {
1144 if((v[i])->GetProcessName() == processName) {
1145 G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1146 if(ActiveForParticle(part, p)) {
1147 proc = v[i];
1148 break;
1149 }
1150 }
1151 }
1152 return proc;
1153}
1154
1155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1156
1159 const G4String& processName)
1160{
1161 G4VEmProcess* proc = nullptr;
1162 auto v = manager->GetEmProcessVector();
1163 G4int n = v.size();
1164 for(G4int i=0; i<n; ++i) {
1165 auto pName = v[i]->GetProcessName();
1166 if(pName == "GammaGeneralProc") {
1167 proc = v[i]->GetEmProcess(processName);
1168 break;
1169 } else if(pName == processName) {
1170 G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1171 if(ActiveForParticle(part, p)) {
1172 proc = v[i];
1173 break;
1174 }
1175 }
1176 }
1177 return proc;
1178}
1179
1180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1181
1184 const G4String& processName)
1185{
1186 G4VMultipleScattering* proc = 0;
1187 const std::vector<G4VMultipleScattering*> v =
1189 G4int n = v.size();
1190 for(G4int i=0; i<n; ++i) {
1191 if((v[i])->GetProcessName() == processName) {
1192 G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1193 if(ActiveForParticle(part, p)) {
1194 proc = v[i];
1195 break;
1196 }
1197 }
1198 }
1199 return proc;
1200}
1201
1202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1203
1205 const G4String& processName)
1206{
1207 G4VProcess* proc = 0;
1208 const G4ProcessManager* procman = part->GetProcessManager();
1209 G4ProcessVector* pv = procman->GetProcessList();
1210 G4int nproc = pv->size();
1211 for(G4int i=0; i<nproc; ++i) {
1212 if(processName == (*pv)[i]->GetProcessName()) {
1213 proc = (*pv)[i];
1214 break;
1215 }
1216 }
1217 return proc;
1218}
1219
1220//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1221
1223 G4VProcess* proc)
1224{
1225 G4ProcessManager* pm = part->GetProcessManager();
1226 G4ProcessVector* pv = pm->GetProcessList();
1227 G4int n = pv->size();
1228 G4bool res = false;
1229 for(G4int i=0; i<n; ++i) {
1230 if((*pv)[i] == proc) {
1231 if(pm->GetProcessActivation(i)) { res = true; }
1232 break;
1233 }
1234 }
1235 return res;
1236}
1237
1238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1239
1241{
1242 if(mat) {
1243 currentMaterial = mat;
1245 } else {
1246 currentMaterial = 0;
1248 }
1249}
1250
1251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1252
1254{
1256}
1257
1258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1259
1261{
1262 G4bool isFound = false;
1263 if(currentMaterial) {
1265 for(size_t i=0; i<nn; ++i) {
1266 if(Z == currentMaterial->GetElement(i)->GetZasInt()) {
1267 isFound = true;
1268 break;
1269 }
1270 }
1271 }
1272 if(!isFound) {
1274 }
1275}
1276
1277//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1278
1280{
1281 verbose = verb;
1282}
1283
1284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1285
G4AtomicShellEnumerator
static const G4double reg
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double gram
Definition: G4SIunits.hh:163
static constexpr double mole
Definition: G4SIunits.hh:279
static constexpr double barn
Definition: G4SIunits.hh:85
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double cm
Definition: G4SIunits.hh:99
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4int GetZasInt() const
Definition: G4Element.hh:132
const G4ParticleDefinition * FindParticle(const G4String &)
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4String currentName
G4VEmModel * loweModel
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
std::vector< const G4Material * > localMaterials
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
const G4ParticleDefinition * theGenericIon
G4double cutenergy[3]
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
const G4ParticleDefinition * currentParticle
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4LossTableManager * manager
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)
G4VProcess * curProcess
G4EmParameters * theParameters
G4bool UpdateParticle(const G4ParticleDefinition *, G4double kinEnergy)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4EmCorrections * corr
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)
G4VEnergyLossProcess * FindEnLossProcess(const G4ParticleDefinition *, const G4String &processName)
G4String currentProcessName
G4bool FindEmModel(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
G4VEmModel * currentModel
G4VMultipleScattering * FindMscProcess(const G4ParticleDefinition *, const G4String &processName)
void FindLambdaTable(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy, G4int &proctype)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
G4String currentMaterialName
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)
G4bool ActiveForParticle(const G4ParticleDefinition *part, G4VProcess *proc)
G4VEnergyLossProcess * currentProcess
G4DynamicParticle * dynParticle
const G4ParticleDefinition * lambdaParticle
const G4Material * currentMaterial
G4bool UpdateCouple(const G4Material *, G4double cut)
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 *)
G4String currentParticleName
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4IonTable * ionTable
G4VEmProcess * FindDiscreteProcess(const G4ParticleDefinition *, const G4String &processName)
void CheckMaterial(G4int Z)
const G4MaterialCutsCouple * currentCouple
std::vector< const G4MaterialCutsCouple * > localCouples
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)
const G4PhysicsTable * currentLambda
G4double chargeSquare
G4ionEffectiveCharge * ionEffCharge
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
const G4Material * FindMaterial(const G4String &)
G4NistManager * nist
const G4ParticleDefinition * baseParticle
std::vector< G4double > localCuts
const G4Material * cutMaterial
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4EmParameters * Instance()
G4bool BuildCSDARange() const
G4double LowestElectronEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4LossTableManager * Instance()
const std::vector< G4VEmProcess * > & GetEmProcessVector()
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4VAtomDeexcitation * AtomDeexcitation()
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4EmCorrections * EmCorrections()
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetDensity() const
Definition: G4Material.hh:176
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4String & GetName() const
Definition: G4Material.hh:173
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:686
G4Material * FindOrBuildSimpleMaterial(G4int Z, G4bool warning=false)
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4bool GetProcessActivation(G4VProcess *aProcess) const
G4ProcessVector * GetProcessList() const
std::size_t size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
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 G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:209
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
void SetFluctuationFlag(G4bool val)
Definition: G4VEmModel.hh:732
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:223
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:440
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:351
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:399
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:382
const G4String & GetName() const
Definition: G4VEmModel.hh:837
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:870
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
G4VEmModel * EmModel(size_t index=0) const
G4PhysicsTable * LambdaTable() const
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t idxCouple) const
virtual G4VEmProcess * GetEmProcess(const G4String &name)
const G4ParticleDefinition * BaseParticle() const
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idxCouple) const
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4PhysicsTable * LambdaTable() const
G4PhysicsTable * DEDXTable() const
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:319
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
static constexpr double eV
static constexpr double nm
Definition: SystemOfUnits.h:93
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
string material
Definition: eplot.py:19
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62