Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEmProcess.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// GEANT4 Class file
29//
30//
31// File name: G4VEmProcess
32//
33// Author: Vladimir Ivanchenko on base of Laszlo Urban code
34//
35// Creation date: 01.10.2003
36//
37// Modifications: by V.Ivanchenko
38//
39// Class Description: based class for discrete and rest/discrete EM processes
40//
41
42// -------------------------------------------------------------------
43//
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
47#include "G4VEmProcess.hh"
49#include "G4SystemOfUnits.hh"
50#include "G4ProcessManager.hh"
51#include "G4LossTableManager.hh"
52#include "G4LossTableBuilder.hh"
53#include "G4Step.hh"
55#include "G4VEmModel.hh"
56#include "G4DataVector.hh"
57#include "G4PhysicsTable.hh"
58#include "G4EmDataHandler.hh"
59#include "G4PhysicsLogVector.hh"
60#include "G4VParticleChange.hh"
62#include "G4Region.hh"
63#include "G4Gamma.hh"
64#include "G4Electron.hh"
65#include "G4Positron.hh"
67#include "G4EmBiasingManager.hh"
68#include "G4GenericIon.hh"
69#include "G4Log.hh"
70#include <iostream>
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75 G4VDiscreteProcess(name, type),
76 secondaryParticle(nullptr),
77 buildLambdaTable(true),
78 numberOfModels(0),
79 theLambdaTable(nullptr),
80 theLambdaTablePrim(nullptr),
81 integral(false),
82 applyCuts(false),
83 startFromNull(false),
84 splineFlag(true),
85 isIon(false),
86 currentCouple(nullptr),
87 isTheMaster(true),
88 masterProc(nullptr),
89 theData(nullptr),
90 currentModel(nullptr),
91 particle(nullptr),
92 currentParticle(nullptr)
93{
96
97 // Size of tables assuming spline
98 minKinEnergy = 0.1*keV;
99 maxKinEnergy = 100.0*TeV;
100 nLambdaBins = 84;
101 minKinEnergyPrim = DBL_MAX;
102 actBinning = actSpline = actMinKinEnergy = actMaxKinEnergy = false;
103
104 // default lambda factor
105 lambdaFactor = 0.8;
106 logLambdaFactor = G4Log(lambdaFactor);
107
108 // default limit on polar angle
109 biasFactor = fFactor = 1.0;
110
111 // particle types
114 thePositron = G4Positron::Positron();
115
116 theCuts = theCutsGamma = theCutsElectron = theCutsPositron = nullptr;
117
120 secParticles.reserve(5);
121
122 baseMaterial = currentMaterial = nullptr;
123
127 massRatio = 1.0;
128
130
131 modelManager = new G4EmModelManager();
132 biasManager = nullptr;
133 biasFlag = false;
134 weightFlag = false;
136 lManager->Register(this);
140
141 secID = fluoID = augerID = biasID = -1;
142 mainSecondaries = 100;
143 if("phot" == GetProcessName() || "compt" == GetProcessName()
144 || "e-_G4DNAIonisation" == GetProcessName()
145 || "hydrogen_G4DNAIonisation" == GetProcessName()
146 || "helium_G4DNAIonisation" == GetProcessName()
147 || "alpha_G4DNAIonisation" == GetProcessName()
148 || "alpha+_G4DNAIonisation" == GetProcessName()
149 || "proton_G4DNAIonisation" == GetProcessName()
150 || "GenericIon_G4DNAIonisation" == GetProcessName() )
151 {
152 mainSecondaries = 1;
153 }
154}
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157
159{
160 /*
161 if(1 < verboseLevel) {
162 G4cout << "G4VEmProcess destruct " << GetProcessName()
163 << " " << this << " " << theLambdaTable <<G4endl;
164 }
165 */
166 if(isTheMaster) {
167 delete theData;
168 theData = nullptr;
169 }
170 delete modelManager;
171 delete biasManager;
172 lManager->DeRegister(this);
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
177void G4VEmProcess::Clear()
178{
179 currentCouple = nullptr;
180 preStepLambda = 0.0;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184
186 const G4Material*)
187{
188 return 0.0;
189}
190
191//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192
194 const G4Region* region)
195{
196 G4VEmFluctuationModel* fm = nullptr;
197 modelManager->AddEmModel(order, p, fm, region);
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
202
204{
205 for(auto & em : emModels) { if(em == ptr) { return; } }
206 emModels.push_back(ptr);
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210
212{
213 return (index < emModels.size()) ? emModels[index] : nullptr;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217
219 G4double emin, G4double emax)
220{
221 modelManager->UpdateEmModel(nam, emin, emax);
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
225
227{
228 return modelManager->NumberOfModels();
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232
234{
235 return modelManager->NumberOfRegionModels(couple_index);
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
240G4VEmModel* G4VEmProcess::GetRegionModel(G4int idx, size_t couple_index) const
241{
242 return modelManager->GetRegionModel(idx, couple_index);
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246
248{
249 return modelManager->GetModel(idx, ver);
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253
255{
257
258 if(!particle) { SetParticle(&part); }
259
260 if(part.GetParticleType() == "nucleus" &&
261 part.GetParticleSubType() == "generic") {
262
263 G4String pname = part.GetParticleName();
264 if(pname != "deuteron" && pname != "triton" &&
265 pname != "alpha" && pname != "He3" &&
266 pname != "alpha+" && pname != "helium" &&
267 pname != "hydrogen") {
268
269 particle = G4GenericIon::GenericIon();
270 isIon = true;
271 }
272 }
273
274 if(1 < verboseLevel) {
275 G4cout << "G4VEmProcess::PreparePhysicsTable() for "
276 << GetProcessName()
277 << " and particle " << part.GetParticleName()
278 << " local particle " << particle->GetParticleName()
279 << G4endl;
280 }
281
282 if(particle != &part) { return; }
283
285
287
288 Clear();
289 InitialiseProcess(particle);
290
291 const G4ProductionCutsTable* theCoupleTable=
293 size_t n = theCoupleTable->GetTableSize();
294
295 theEnergyOfCrossSectionMax.resize(n, 0.0);
296 theCrossSectionMax.resize(n, DBL_MAX);
297
298 // initialisation of the process
299 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
300 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
301 if(!actSpline) { splineFlag = theParameters->Spline(); }
302
303 if(isTheMaster) {
305 if(!theData) { theData = new G4EmDataHandler(2); }
306 } else {
308 }
309 applyCuts = theParameters->ApplyCuts();
310 lambdaFactor = theParameters->LambdaFactor();
311 logLambdaFactor = G4Log(lambdaFactor);
313
314 // initialisation of models
315 numberOfModels = modelManager->NumberOfModels();
316 for(G4int i=0; i<numberOfModels; ++i) {
317 G4VEmModel* mod = modelManager->GetModel(i);
318 if(0 == i) { currentModel = mod; }
321 if(mod->HighEnergyLimit() > maxKinEnergy) {
322 mod->SetHighEnergyLimit(maxKinEnergy);
323 }
324 }
325
326 if(lManager->AtomDeexcitation()) { modelManager->SetFluoFlag(true); }
327 theCuts = modelManager->Initialise(particle,secondaryParticle,
328 2.,verboseLevel);
329 theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
330 theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
331 theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
332
333 // prepare tables
334 if(buildLambdaTable && isTheMaster){
335 theLambdaTable = theData->MakeTable(0);
336 bld->InitialiseBaseMaterials(theLambdaTable);
337 }
338 // high energy table
339 if(isTheMaster && minKinEnergyPrim < maxKinEnergy){
340 theLambdaTablePrim = theData->MakeTable(1);
341 bld->InitialiseBaseMaterials(theLambdaTablePrim);
342 }
344 // forced biasing
345 if(biasManager) {
347 biasFlag = false;
348 }
349 // defined ID of secondary particles
350 G4String nam1 = GetProcessName();
352 if(100 > mainSecondaries) {
353 G4String nam2 = nam1 + "_fluo" ;
354 G4String nam3 = nam1 + "_auger";
355 G4String nam4 = nam1 + "_split";
359 }
360}
361
362//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
363
365{
366 if(!masterProc) {
367 if(isTheMaster) { masterProc = this; }
368 else { masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());}
369 }
370
371 G4String num = part.GetParticleName();
372 if(1 < verboseLevel) {
373 G4cout << "### G4VEmProcess::BuildPhysicsTable() for "
374 << GetProcessName()
375 << " and particle " << num
376 << " buildLambdaTable= " << buildLambdaTable
377 << " isTheMaster= " << isTheMaster
378 << " " << masterProc
379 << G4endl;
380 }
381
382 if(particle == &part) {
383
384 // worker initialisation
385 if(!isTheMaster) {
386 theLambdaTable = masterProc->LambdaTable();
387 theLambdaTablePrim = masterProc->LambdaTablePrim();
388
389 if(theLambdaTable) { FindLambdaMax(); }
390
391 // local initialisation of models
392 G4bool printing = true;
393 numberOfModels = modelManager->NumberOfModels();
394 for(G4int i=0; i<numberOfModels; ++i) {
395 G4VEmModel* mod = GetModelByIndex(i, printing);
396 G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing);
397 //G4cout << i << ". " << mod << " " << mod0 << " "
398 // << particle->GetParticleName() << G4endl;
399 mod->InitialiseLocal(particle, mod0);
400 }
401 // master thread
402 } else {
403 if(buildLambdaTable || minKinEnergyPrim < maxKinEnergy) {
404 BuildLambdaTable();
405 }
406 }
407 }
408
409 // explicitly defined printout by particle name
410 if(1 < verboseLevel ||
411 (0 < verboseLevel && (num == "gamma" || num == "e-" ||
412 num == "e+" || num == "mu+" ||
413 num == "mu-" || num == "proton"||
414 num == "pi+" || num == "pi-" ||
415 num == "kaon+" || num == "kaon-" ||
416 num == "alpha" || num == "anti_proton" ||
417 num == "GenericIon"|| num == "alpha++" ||
418 num == "alpha+" || num == "helium" ||
419 num == "hydrogen")))
420 {
421 StreamInfo(G4cout, part);
422 }
423
424 if(1 < verboseLevel) {
425 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
426 << GetProcessName()
427 << " and particle " << num
428 << G4endl;
429 }
430}
431
432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
433
434void G4VEmProcess::BuildLambdaTable()
435{
436 if(1 < verboseLevel) {
437 G4cout << "G4EmProcess::BuildLambdaTable() for process "
438 << GetProcessName() << " and particle "
439 << particle->GetParticleName() << " " << this
440 << G4endl;
441 }
442
443 // Access to materials
444 const G4ProductionCutsTable* theCoupleTable=
446 size_t numOfCouples = theCoupleTable->GetTableSize();
447
449
450 G4PhysicsLogVector* aVector = nullptr;
451 G4PhysicsLogVector* aVectorPrim = nullptr;
452 G4PhysicsLogVector* bVectorPrim = nullptr;
453
455 G4int nbin =
456 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
457 scale = G4Log(scale);
458 if(actBinning) { nbin = std::max(nbin, nLambdaBins); }
459 G4double emax1 = std::min(maxKinEnergy, minKinEnergyPrim);
460
461 for(size_t i=0; i<numOfCouples; ++i) {
462
463 if (bld->GetFlag(i)) {
464
465 // create physics vector and fill it
466 const G4MaterialCutsCouple* couple =
467 theCoupleTable->GetMaterialCutsCouple(i);
468
469 // build main table
470 if(buildLambdaTable) {
471 delete (*theLambdaTable)[i];
472
473 // if start from zero then change the scale
474 G4double emin = minKinEnergy;
475 G4bool startNull = false;
476 if(startFromNull) {
477 G4double e = MinPrimaryEnergy(particle,couple->GetMaterial());
478 if(e >= emin) {
479 emin = e;
480 startNull = true;
481 }
482 }
483 G4double emax = emax1;
484 if(emax <= emin) { emax = 2*emin; }
485 G4int bin = G4lrint(nbin*G4Log(emax/emin)/scale);
486 if(bin < 3) { bin = 3; }
487 aVector = new G4PhysicsLogVector(emin, emax, bin);
488 aVector->SetSpline(splineFlag);
489 modelManager->FillLambdaVector(aVector, couple, startNull);
490 if(splineFlag) { aVector->FillSecondDerivatives(); }
491 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector);
492 }
493 // build high energy table
494 if(minKinEnergyPrim < maxKinEnergy) {
495 delete (*theLambdaTablePrim)[i];
496
497 // start not from zero
498 if(!bVectorPrim) {
499 G4int bin = G4lrint(nbin*G4Log(maxKinEnergy/minKinEnergyPrim)/scale);
500 if(bin < 3) { bin = 3; }
501 aVectorPrim =
502 new G4PhysicsLogVector(minKinEnergyPrim, maxKinEnergy, bin);
503 bVectorPrim = aVectorPrim;
504 } else {
505 aVectorPrim = new G4PhysicsLogVector(*bVectorPrim);
506 }
507 // always use spline
508 aVectorPrim->SetSpline(splineFlag);
509 modelManager->FillLambdaVector(aVectorPrim, couple, false,
511 aVectorPrim->FillSecondDerivatives();
512 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTablePrim, i,
513 aVectorPrim);
514 }
515 }
516 }
517
518 if(buildLambdaTable) { FindLambdaMax(); }
519
520 if(1 < verboseLevel) {
521 G4cout << "Lambda table is built for "
522 << particle->GetParticleName()
523 << G4endl;
524 }
525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528
529void G4VEmProcess::StreamInfo(std::ostream& out,
530 const G4ParticleDefinition& part, G4bool rst) const
531{
532 G4String indent = (rst ? " " : "");
533 out << std::setprecision(6);
534 out << G4endl << indent << GetProcessName() << ": ";
535 if (!rst) {
536 out << " for " << part.GetParticleName();
537 if (integral) { out << ","; }
538 }
539 if(integral) { out << " integral:1 "; }
540 if(applyCuts) { out << " applyCuts:1 "; }
541 out << " SubType=" << GetProcessSubType();
542 if(biasFactor != 1.0) { out << " BiasingFactor= " << biasFactor; }
543 out << " BuildTable=" << buildLambdaTable << G4endl;
544 if(buildLambdaTable) {
545 if(particle == &part) {
546 size_t length = theLambdaTable->length();
547 for(size_t i=0; i<length; ++i) {
548 G4PhysicsVector* v = (*theLambdaTable)[i];
549 if(v) {
550 out << " Lambda table from ";
551 G4double emin = v->Energy(0);
552 G4double emax = v->GetMaxEnergy();
553 G4int nbin = v->GetVectorLength() - 1;
554 if(emin > minKinEnergy) { out << "threshold "; }
555 else { out << G4BestUnit(emin,"Energy"); }
556 out << " to "
557 << G4BestUnit(emax,"Energy")
558 << ", " << G4lrint(nbin/std::log10(emax/emin))
559 << " bins/decade, spline: "
560 << splineFlag << G4endl;
561 break;
562 }
563 }
564 } else {
565 out << " Used Lambda table of "
566 << particle->GetParticleName() << G4endl;
567 }
568 }
569 if(minKinEnergyPrim < maxKinEnergy) {
570 if(particle == &part) {
571 size_t length = theLambdaTablePrim->length();
572 for(size_t i=0; i<length; ++i) {
573 G4PhysicsVector* v = (*theLambdaTablePrim)[i];
574 if(v) {
575 out << " LambdaPrime table from "
576 << G4BestUnit(v->Energy(0),"Energy")
577 << " to "
578 << G4BestUnit(v->GetMaxEnergy(),"Energy")
579 << " in " << v->GetVectorLength()-1
580 << " bins " << G4endl;
581 break;
582 }
583 }
584 } else {
585 out << " Used LambdaPrime table of "
586 << particle->GetParticleName() << G4endl;
587 }
588 }
590 modelManager->DumpModelList(out, verboseLevel);
591
592 if(verboseLevel > 2 && buildLambdaTable) {
593 out << " LambdaTable address= " << theLambdaTable << G4endl;
594 if(theLambdaTable && particle == &part) {
595 out << (*theLambdaTable) << G4endl;
596 }
597 }
598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
603{
604 // reset parameters for the new track
605 currentParticle = track->GetParticleDefinition();
608
609 if(isIon) { massRatio = proton_mass_c2/currentParticle->GetPDGMass(); }
610
611 // forced biasing only for primary particles
612 if(biasManager) {
613 if(0 == track->GetParentID()) {
614 // primary particle
615 biasFlag = true;
617 }
618 }
619}
620
621//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
622
624 const G4Track& track,
625 G4double previousStepSize,
627{
629 G4double x = DBL_MAX;
630
634 G4double scaledEnergy = preStepKinEnergy*massRatio;
635 SelectModel(scaledEnergy, currentCoupleIndex);
636
637 if(!currentModel->IsActive(scaledEnergy)) {
640 return x;
641 }
642
643 // forced biasing only for primary particles
644 if(biasManager) {
645 if(0 == track.GetParentID()) {
646 if(biasFlag &&
648 return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
649 }
650 }
651 }
652
653 // compute mean free path
655 if (integral) {
656 ComputeIntegralLambda(preStepKinEnergy, preStepLogKinEnergy);
657 } else {
659 }
660
661 // zero cross section
662 if(preStepLambda <= 0.0) {
665 }
666 }
667
668 // non-zero cross section
669 if(preStepLambda > 0.0) {
670
672
673 // beggining of tracking (or just after DoIt of this process)
676
677 } else if(currentInteractionLength < DBL_MAX) {
678
680 previousStepSize/currentInteractionLength;
683 }
684
685 // new mean free path and step limit for the next step
688 }
689 return x;
690}
691
692//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
693
694void G4VEmProcess::ComputeIntegralLambda(G4double e, G4double loge)
695{
696 // condition to skip recomputation of cross section
697 const G4double epeak = theEnergyOfCrossSectionMax[currentCoupleIndex];
698 if(e <= epeak && e/lambdaFactor >= mfpKinEnergy) { return; }
699
700 // recomputation is needed
701 if (e <= epeak) {
702 preStepLambda = GetCurrentLambda(e, loge);
703 mfpKinEnergy = e;
704 } else {
705 const G4double e1 = e*lambdaFactor;
706 if (e1 > epeak) {
707 preStepLambda = GetCurrentLambda(e, loge);
708 mfpKinEnergy = e;
709 const G4double preStepLambda1 = GetCurrentLambda(e1,loge+logLambdaFactor);
710 if (preStepLambda1 > preStepLambda) {
711 mfpKinEnergy = e1;
712 preStepLambda = preStepLambda1;
713 }
714 } else {
715 preStepLambda = fFactor*theCrossSectionMax[currentCoupleIndex];
716 mfpKinEnergy = epeak;
717 }
718 }
719}
720
721//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
722
724 const G4Step& step)
725{
726 // In all cases clear number of interaction lengths
729
731
732 // Do not make anything if particle is stopped, the annihilation then
733 // should be performed by the AtRestDoIt!
734 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
735
736 const G4double finalT = track.GetKineticEnergy();
737 const G4double logFinalT = track.GetDynamicParticle()->GetLogKineticEnergy();
738
739 // forced process - should happen only once per track
740 if(biasFlag) {
742 biasFlag = false;
743 }
744 }
745
746 // Integral approach
747 if (integral) {
748 G4double lx = GetLambda(finalT, currentCouple, logFinalT);
749 if(preStepLambda<lx && 1 < verboseLevel) {
750 G4cout << "WARNING: for " << currentParticle->GetParticleName()
751 << " and " << GetProcessName()
752 << " E(MeV)= " << finalT/MeV
753 << " preLambda= " << preStepLambda << " < "
754 << lx << " (postLambda) "
755 << G4endl;
756 }
757
758 if(preStepLambda*G4UniformRand() > lx) {
760 return &fParticleChange;
761 }
762 }
763
764 G4double scaledEnergy = finalT*massRatio;
765 SelectModel(scaledEnergy, currentCoupleIndex);
766 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
767
768 // define new weight for primary and secondaries
770 if(weightFlag) {
771 weight /= biasFactor;
773 }
774
775
776 if(1 < verboseLevel) {
777 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
778 << finalT/MeV
779 << " MeV; model= (" << currentModel->LowEnergyLimit()
780 << ", " << currentModel->HighEnergyLimit() << ")"
781 << G4endl;
782 }
783
784
785 // sample secondaries
786 secParticles.clear();
787 currentModel->SampleSecondaries(&secParticles,
789 track.GetDynamicParticle(),
790 (*theCuts)[currentCoupleIndex]);
791
792 G4int num0 = secParticles.size();
793
794 // splitting or Russian roulette
795 if(biasManager) {
797 G4double eloss = 0.0;
799 secParticles, track, currentModel, &fParticleChange, eloss,
801 step.GetPostStepPoint()->GetSafety());
802 if(eloss > 0.0) {
805 }
806 }
807 }
808
809 // save secondaries
810 G4int num = secParticles.size();
811 if(num > 0) {
812
815 G4double time = track.GetGlobalTime();
816
817 for (G4int i=0; i<num; ++i) {
818 if (secParticles[i]) {
821 G4double e = dp->GetKineticEnergy();
822 G4bool good = true;
823 if(applyCuts) {
824 if (p == theGamma) {
825 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
826
827 } else if (p == theElectron) {
828 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
829
830 } else if (p == thePositron) {
831 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
832 e < (*theCutsPositron)[currentCoupleIndex]) {
833 good = false;
834 e += 2.0*electron_mass_c2;
835 }
836 }
837 // added secondary if it is good
838 }
839 if (good) {
840 G4Track* t = new G4Track(dp, time, track.GetPosition());
842 if (biasManager) {
843 t->SetWeight(weight * biasManager->GetWeight(i));
844 } else {
845 t->SetWeight(weight);
846 }
848
849 // define type of secondary
851 else if(i < num0) {
852 if(p == theGamma) {
854 } else {
856 }
857 } else {
859 }
860 /*
861 G4cout << "Secondary(post step) has weight " << t->GetWeight()
862 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
863 << GetProcessName() << " fluoID= " << fluoID
864 << " augerID= " << augerID <<G4endl;
865 */
866 } else {
867 delete dp;
868 edep += e;
869 }
870 }
871 }
873 }
874
877 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
880 }
881
882 return &fParticleChange;
883}
884
885//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
886
888 const G4String& directory,
889 G4bool ascii)
890{
891 G4bool yes = true;
892 if(!isTheMaster) { return yes; }
893
894 if ( theLambdaTable && part == particle) {
895 const G4String& nam =
896 GetPhysicsTableFileName(part,directory,"Lambda",ascii);
897 yes = theLambdaTable->StorePhysicsTable(nam,ascii);
898
899 if ( yes ) {
900 G4cout << "Physics table is stored for " << particle->GetParticleName()
901 << " and process " << GetProcessName()
902 << " in the directory <" << directory
903 << "> " << G4endl;
904 } else {
905 G4cout << "Fail to store Physics Table for "
906 << particle->GetParticleName()
907 << " and process " << GetProcessName()
908 << " in the directory <" << directory
909 << "> " << G4endl;
910 }
911 }
912 if ( theLambdaTablePrim && part == particle) {
913 const G4String& name =
914 GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
915 yes = theLambdaTablePrim->StorePhysicsTable(name,ascii);
916
917 if ( yes ) {
918 G4cout << "Physics table prim is stored for "
919 << particle->GetParticleName()
920 << " and process " << GetProcessName()
921 << " in the directory <" << directory
922 << "> " << G4endl;
923 } else {
924 G4cout << "Fail to store Physics Table Prim for "
925 << particle->GetParticleName()
926 << " and process " << GetProcessName()
927 << " in the directory <" << directory
928 << "> " << G4endl;
929 }
930 }
931 return yes;
932}
933
934//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
935
937 const G4String& directory,
938 G4bool ascii)
939{
940 if(1 < verboseLevel) {
941 G4cout << "G4VEmProcess::RetrievePhysicsTable() for "
942 << part->GetParticleName() << " and process "
943 << GetProcessName() << G4endl;
944 }
945 G4bool yes = true;
946
947 if((!buildLambdaTable && minKinEnergyPrim > maxKinEnergy)
948 || particle != part) { return yes; }
949
950 const G4String particleName = part->GetParticleName();
951
952 if(buildLambdaTable) {
953 const G4String& filename =
954 GetPhysicsTableFileName(part,directory,"Lambda",ascii);
956 filename,ascii);
957 if ( yes ) {
958 if (0 < verboseLevel) {
959 G4cout << "Lambda table for " << particleName
960 << " is Retrieved from <"
961 << filename << ">"
962 << G4endl;
963 }
964 if(theParameters->Spline()) {
965 size_t n = theLambdaTable->length();
966 for(size_t i=0; i<n; ++i) {
967 if((* theLambdaTable)[i]) {
968 (* theLambdaTable)[i]->SetSpline(true);
969 }
970 }
971 }
972 } else {
973 if (1 < verboseLevel) {
974 G4cout << "Lambda table for " << particleName << " in file <"
975 << filename << "> is not exist"
976 << G4endl;
977 }
978 }
979 }
980 if(minKinEnergyPrim < maxKinEnergy) {
981 const G4String& filename =
982 GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
983 yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTablePrim,
984 filename,ascii);
985 if ( yes ) {
986 if (0 < verboseLevel) {
987 G4cout << "Lambda table prim for " << particleName
988 << " is Retrieved from <"
989 << filename << ">"
990 << G4endl;
991 }
992 if(theParameters->Spline()) {
993 size_t n = theLambdaTablePrim->length();
994 for(size_t i=0; i<n; ++i) {
995 if((* theLambdaTablePrim)[i]) {
996 (* theLambdaTablePrim)[i]->SetSpline(true);
997 }
998 }
999 }
1000 } else {
1001 if (1 < verboseLevel) {
1002 G4cout << "Lambda table prim for " << particleName << " in file <"
1003 << filename << "> is not exist"
1004 << G4endl;
1005 }
1006 }
1007 }
1008 return yes;
1009}
1010
1011//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1012
1013G4double
1015 const G4MaterialCutsCouple* couple,
1016 G4double logKinEnergy)
1017{
1018 // Cross section per atom is calculated
1019 DefineMaterial(couple);
1020 G4double cross = 0.0;
1021 if(buildLambdaTable) {
1022 cross = GetCurrentLambda(kineticEnergy,
1023 (logKinEnergy < DBL_MAX) ? logKinEnergy : G4Log(kineticEnergy));
1024 } else {
1025 SelectModel(kineticEnergy, currentCoupleIndex);
1026 if(currentModel) {
1027 cross = fFactor*currentModel->CrossSectionPerVolume(currentMaterial,
1028 currentParticle,
1029 kineticEnergy);
1030 }
1031 }
1032 return std::max(cross, 0.0);
1033}
1034
1035//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1036
1038 G4double,
1040{
1042 return G4VEmProcess::MeanFreePath(track);
1043}
1044
1045//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1046
1048{
1049 const G4double kinEnergy = track.GetKineticEnergy();
1050 CurrentSetup(track.GetMaterialCutsCouple(), kinEnergy);
1051 const G4double xs = GetCurrentLambda(kinEnergy,
1053 return (0.0 < xs) ? 1.0/xs : DBL_MAX;
1054}
1055
1056//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1057
1058G4double
1060 G4double Z, G4double A, G4double cut)
1061{
1062 SelectModel(kinEnergy, currentCoupleIndex);
1063 return (currentModel) ?
1064 currentModel->ComputeCrossSectionPerAtom(currentParticle, kinEnergy,
1065 Z, A, cut) : 0.0;
1066}
1067
1068//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1069
1070void G4VEmProcess::FindLambdaMax()
1071{
1072 if(1 < verboseLevel) {
1073 G4cout << "### G4VEmProcess::FindLambdaMax: "
1074 << particle->GetParticleName()
1075 << " and process " << GetProcessName() << " " << G4endl;
1076 }
1077 size_t n = theLambdaTable->length();
1078 G4PhysicsVector* pv;
1079 G4double e, ss, emax, smax;
1080
1081 size_t i;
1082
1083 // first loop on existing vectors
1084 for (i=0; i<n; ++i) {
1085 pv = (*theLambdaTable)[i];
1086 if(pv) {
1087 size_t nb = pv->GetVectorLength();
1088 emax = DBL_MAX;
1089 smax = 0.0;
1090 if(nb > 0) {
1091 for (size_t j=0; j<nb; ++j) {
1092 e = pv->Energy(j);
1093 ss = (*pv)(j);
1094 if(ss > smax) {
1095 smax = ss;
1096 emax = e;
1097 }
1098 }
1099 }
1100 theEnergyOfCrossSectionMax[i] = emax;
1101 theCrossSectionMax[i] = smax;
1102 if(1 < verboseLevel) {
1103 G4cout << "For " << particle->GetParticleName()
1104 << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
1105 << " lambda= " << smax << G4endl;
1106 }
1107 }
1108 }
1109 // second loop using base materials
1110 for (i=0; i<n; ++i) {
1111 pv = (*theLambdaTable)[i];
1112 if(!pv){
1113 G4int j = (*theDensityIdx)[i];
1114 theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
1115 theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
1116 }
1117 }
1118}
1119
1120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1121
1124{
1125 G4PhysicsVector* v =
1126 new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nLambdaBins);
1128 return v;
1129}
1130
1131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1132
1134{
1135 const G4Element* elm =
1136 (currentModel) ? currentModel->GetCurrentElement() : nullptr;
1137 return elm;
1138}
1139
1140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1141
1143{
1144 if(f > 0.0) {
1145 biasFactor = f;
1146 weightFlag = flag;
1147 if(1 < verboseLevel) {
1148 G4cout << "### SetCrossSectionBiasingFactor: for "
1149 << particle->GetParticleName()
1150 << " and process " << GetProcessName()
1151 << " biasFactor= " << f << " weightFlag= " << flag
1152 << G4endl;
1153 }
1154 }
1155}
1156
1157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1158
1159void
1161 G4bool flag)
1162{
1164 if(1 < verboseLevel) {
1165 G4cout << "### ActivateForcedInteraction: for "
1166 << particle->GetParticleName()
1167 << " and process " << GetProcessName()
1168 << " length(mm)= " << length/mm
1169 << " in G4Region <" << r
1170 << "> weightFlag= " << flag
1171 << G4endl;
1172 }
1173 weightFlag = flag;
1175}
1176
1177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1178
1179void
1181 G4double factor,
1182 G4double energyLimit)
1183{
1184 if (0.0 <= factor) {
1185
1186 // Range cut can be applied only for e-
1187 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1188 { return; }
1189
1191 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1192 if(1 < verboseLevel) {
1193 G4cout << "### ActivateSecondaryBiasing: for "
1194 << " process " << GetProcessName()
1195 << " factor= " << factor
1196 << " in G4Region <" << region
1197 << "> energyLimit(MeV)= " << energyLimit/MeV
1198 << G4endl;
1199 }
1200 }
1201}
1202
1203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1204
1206{
1207 if(5 < n && n < 10000000) {
1208 nLambdaBins = n;
1209 actBinning = true;
1210 } else {
1211 G4double e = (G4double)n;
1212 PrintWarning("SetLambdaBinning", e);
1213 }
1214}
1215
1216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1217
1219{
1220 if(1.e-3*eV < e && e < maxKinEnergy) {
1221 nLambdaBins = G4lrint(nLambdaBins*G4Log(maxKinEnergy/e)
1222 /G4Log(maxKinEnergy/minKinEnergy));
1223 minKinEnergy = e;
1224 actMinKinEnergy = true;
1225 } else { PrintWarning("SetMinKinEnergy", e); }
1226}
1227
1228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1229
1231{
1232 if(minKinEnergy < e && e < 1.e+6*TeV) {
1233 nLambdaBins = G4lrint(nLambdaBins*G4Log(e/minKinEnergy)
1234 /G4Log(maxKinEnergy/minKinEnergy));
1235 maxKinEnergy = e;
1236 actMaxKinEnergy = true;
1237 } else { PrintWarning("SetMaxKinEnergy", e); }
1238}
1239
1240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1241
1243{
1244 if(theParameters->MinKinEnergy() <= e &&
1245 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
1246 else { PrintWarning("SetMinKinEnergyPrim", e); }
1247}
1248
1249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1250
1252{
1253 return (nam == GetProcessName()) ? this : nullptr;
1254}
1255
1256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1257
1258G4double
1260{
1261 CurrentSetup(couple, kinEnergy);
1262 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
1263}
1264
1265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1266
1267void G4VEmProcess::PrintWarning(G4String tit, G4double val)
1268{
1269 G4String ss = "G4VEmProcess::" + tit;
1271 ed << "Parameter is out of range: " << val
1272 << " it will have no effect!\n" << " Process "
1273 << GetProcessName() << " nbins= " << theParameters->NumberOfBins()
1274 << " Emin(keV)= " << theParameters->MinKinEnergy()/keV
1275 << " Emax(GeV)= " << theParameters->MaxKinEnergy()/GeV;
1276 G4Exception(ss, "em0044", JustWarning, ed);
1277}
1278
1279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1280
1281void G4VEmProcess::ProcessDescription(std::ostream& out) const
1282{
1283 if(particle) {
1284 StreamInfo(out, *particle, true);
1285 }
1286}
1287
1288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double A(double temperature)
@ fIsCrossSectionPrim
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ProcessType
@ idxG4ElectronCut
@ idxG4GammaCut
@ idxG4PositronCut
#define G4BestUnit(a, b)
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetLogKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4bool ForcedInteractionRegion(G4int coupleIdx)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
G4double GetWeight(G4int i)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
G4PhysicsTable * MakeTable(size_t idx)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
void UpdateEmModel(const G4String &model_name, G4double emin, G4double emax)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
void SetFluoFlag(G4bool val)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
G4int NumberOfRegionModels(size_t index_couple) const
G4VEmModel * GetRegionModel(G4int idx, size_t index_couple)
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
static G4EmParameters * Instance()
G4int NumberOfBins() const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MscThetaLimit() const
void DefineRegParamForEM(G4VEmProcess *) const
G4int Verbose() const
G4int WorkerVerbose() const
G4bool ApplyCuts() const
G4double MaxKinEnergy() const
G4bool Spline() const
G4double LambdaFactor() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
G4bool GetFlag(size_t idx)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
void InitializeForPostStep(const G4Track &)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static G4int Register(const G4String &)
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
std::size_t length() const
G4double Energy(std::size_t index) const
G4double GetMaxEnergy() const
void FillSecondDerivatives()
void SetSpline(G4bool)
std::size_t GetVectorLength() const
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetSafety() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetCreatorModelIndex(G4int idx)
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:792
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:359
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:729
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:495
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:220
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double MeanFreePath(const G4Track &track)
G4int mainSecondaries
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
virtual void StreamProcessInfo(std::ostream &) const
Definition: G4VEmProcess.hh:95
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:74
G4LossTableManager * lManager
G4EmBiasingManager * biasManager
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
G4VEmModel * EmModel(size_t index=0) const
G4double preStepLogKinEnergy
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4PhysicsTable * LambdaTable() const
void SetMinKinEnergy(G4double e)
size_t basedCoupleIndex
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
virtual void StartTracking(G4Track *) override
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4double mfpKinEnergy
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
void SetLambdaBinning(G4int nbins)
const std::vector< G4double > * theDensityFactor
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
virtual void ProcessDescription(std::ostream &outFile) const override
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void UpdateEmModel(const G4String &, G4double, G4double)
virtual G4VEmProcess * GetEmProcess(const G4String &name)
G4double MaxKinEnergy() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t index)
G4bool isTheMaster
std::vector< G4DynamicParticle * > secParticles
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
const G4ParticleDefinition * theElectron
G4EmParameters * theParameters
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * theGamma
G4PhysicsTable * LambdaTablePrim() const
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4double preStepKinEnergy
size_t currentCoupleIndex
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4ParticleChangeForGamma fParticleChange
const std::vector< G4int > * theDensityIdx
void SetParticle(const G4ParticleDefinition *p)
G4VEmModel * GetRegionModel(G4int idx, size_t couple_index) const
void SetMinKinEnergyPrim(G4double e)
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMaxKinEnergy(G4double e)
G4int GetNumberOfRegionModels(size_t couple_index) const
const G4Material * currentMaterial
const G4Element * GetCurrentElement() const
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4int GetNumberOfModels() const
virtual ~G4VEmProcess()
G4double GetParentWeight() const
void ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void ProposeWeight(G4double finalWeight)
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4TrackStatus GetTrackStatus() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define LOG_EKIN_MIN
Definition: templates.hh:98
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62