Geant4 11.2.2
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 "G4EmParameters.hh"
69#include "G4EmProcessSubType.hh"
70#include "G4EmTableUtil.hh"
71#include "G4EmUtility.hh"
72#include "G4DNAModelSubType.hh"
73#include "G4GenericIon.hh"
74#include "G4Log.hh"
75#include <iostream>
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
80 G4VDiscreteProcess(name, type)
81{
82 theParameters = G4EmParameters::Instance();
84
85 // Size of tables
86 minKinEnergy = 0.1*CLHEP::keV;
87 maxKinEnergy = 100.0*CLHEP::TeV;
88
89 // default lambda factor
90 invLambdaFactor = 1.0/lambdaFactor;
91
92 // particle types
93 theGamma = G4Gamma::Gamma();
94 theElectron = G4Electron::Electron();
95 thePositron = G4Positron::Positron();
96
99 secParticles.reserve(5);
100
101 modelManager = new G4EmModelManager();
102 lManager = G4LossTableManager::Instance();
103 lManager->Register(this);
104 isTheMaster = lManager->IsMaster();
105 G4LossTableBuilder* bld = lManager->GetTableBuilder();
106 theDensityFactor = bld->GetDensityFactors();
107 theDensityIdx = bld->GetCoupleIndexes();
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
113{
114 if(isTheMaster) {
115 delete theData;
117 }
118 delete modelManager;
119 delete biasManager;
120 lManager->DeRegister(this);
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124
126 const G4Region* region)
127{
128 if(nullptr == ptr) { return; }
129 G4VEmFluctuationModel* fm = nullptr;
130 modelManager->AddEmModel(order, ptr, fm, region);
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135
137{
138 if(nullptr == ptr) { return; }
139 if(!emModels.empty()) {
140 for(auto & em : emModels) { if(em == ptr) { return; } }
141 }
142 emModels.push_back(ptr);
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146
148{
149 if(nullptr == particle) { SetParticle(&part); }
150
151 if(part.GetParticleType() == "nucleus" &&
152 part.GetParticleSubType() == "generic") {
153
154 G4String pname = part.GetParticleName();
155 if(pname != "deuteron" && pname != "triton" &&
156 pname != "alpha" && pname != "alpha+" &&
157 pname != "helium" && pname != "hydrogen") {
158
159 particle = G4GenericIon::GenericIon();
160 isIon = true;
161 }
162 }
163 if(particle != &part) { return; }
164
165 lManager->PreparePhysicsTable(&part, this);
166
167 // for new run
168 currentCouple = nullptr;
169 preStepLambda = 0.0;
170 fLambdaEnergy = 0.0;
171
172 InitialiseProcess(particle);
173
174 G4LossTableBuilder* bld = lManager->GetTableBuilder();
175 const G4ProductionCutsTable* theCoupleTable=
177 theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
178 theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
179 theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
180
181 // initialisation of the process
182 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
183 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
184
185 applyCuts = theParameters->ApplyCuts();
186 lambdaFactor = theParameters->LambdaFactor();
187 invLambdaFactor = 1.0/lambdaFactor;
188 theParameters->DefineRegParamForEM(this);
189
190 // integral option may be disabled
191 if(!theParameters->Integral()) { fXSType = fEmNoIntegral; }
192
193 // prepare tables
194 if(isTheMaster) {
195 if(nullptr == theData) { theData = new G4EmDataHandler(2); }
196
197 if(buildLambdaTable) {
198 theLambdaTable = theData->MakeTable(0);
199 bld->InitialiseBaseMaterials(theLambdaTable);
200 }
201 // high energy table
202 if(minKinEnergyPrim < maxKinEnergy) {
203 theLambdaTablePrim = theData->MakeTable(1);
204 bld->InitialiseBaseMaterials(theLambdaTablePrim);
205 }
206 }
207 // models
209 numberOfModels = modelManager->NumberOfModels();
210 currentModel = modelManager->GetModel(0);
211 if(nullptr != lManager->AtomDeexcitation()) {
212 modelManager->SetFluoFlag(true);
213 }
214 // forced biasing
215 if(nullptr != biasManager) {
217 biasFlag = false;
218 }
219
220 theCuts =
221 G4EmTableUtil::PrepareEmProcess(this, particle, secondaryParticle,
222 modelManager, maxKinEnergy,
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
230{
231 if(nullptr == masterProc) {
232 if(isTheMaster) { masterProc = this; }
233 else { masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());}
234 }
235 G4int nModels = modelManager->NumberOfModels();
236 G4bool isLocked = theParameters->IsPrintLocked();
237 G4bool toBuild = (buildLambdaTable || minKinEnergyPrim < maxKinEnergy);
238
239 G4EmTableUtil::BuildEmProcess(this, masterProc, particle, &part,
240 nModels, verboseLevel, isTheMaster,
241 isLocked, toBuild, baseMat);
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245
247{
248 G4double scale = theParameters->MaxKinEnergy()/theParameters->MinKinEnergy();
249 G4int nbin =
250 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
251 if(actBinning) { nbin = std::max(nbin, nLambdaBins); }
252 scale = nbin/G4Log(scale);
253
254 G4LossTableBuilder* bld = lManager->GetTableBuilder();
255 G4EmTableUtil::BuildLambdaTable(this, particle, modelManager,
256 bld, theLambdaTable, theLambdaTablePrim,
257 minKinEnergy, minKinEnergyPrim,
258 maxKinEnergy, scale, verboseLevel,
259 startFromNull, splineFlag);
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
263
264void G4VEmProcess::StreamInfo(std::ostream& out,
265 const G4ParticleDefinition& part, G4bool rst) const
266{
267 G4String indent = (rst ? " " : "");
268 out << std::setprecision(6);
269 out << G4endl << indent << GetProcessName() << ": ";
270 if (!rst) {
271 out << " for " << part.GetParticleName();
272 }
273 if(fXSType != fEmNoIntegral) { out << " XStype:" << fXSType; }
274 if(applyCuts) { out << " applyCuts:1 "; }
275 out << " SubType=" << GetProcessSubType();
276 if(biasFactor != 1.0) { out << " BiasingFactor= " << biasFactor; }
277 out << " BuildTable=" << buildLambdaTable << G4endl;
278 if(buildLambdaTable) {
279 if(particle == &part) {
280 for(auto & v : *theLambdaTable) {
281 if(nullptr != v) {
282 out << " Lambda table from ";
283 G4double emin = v->Energy(0);
284 G4double emax = v->GetMaxEnergy();
285 G4int nbin = G4int(v->GetVectorLength() - 1);
286 if(emin > minKinEnergy) { out << "threshold "; }
287 else { out << G4BestUnit(emin,"Energy"); }
288 out << " to "
289 << G4BestUnit(emax,"Energy")
290 << ", " << G4lrint(nbin/std::log10(emax/emin))
291 << " bins/decade, spline: "
292 << splineFlag << G4endl;
293 break;
294 }
295 }
296 } else {
297 out << " Used Lambda table of "
298 << particle->GetParticleName() << G4endl;
299 }
300 }
301 if(minKinEnergyPrim < maxKinEnergy) {
302 if(particle == &part) {
303 for(auto & v : *theLambdaTablePrim) {
304 if(nullptr != v) {
305 out << " LambdaPrime table from "
306 << G4BestUnit(v->Energy(0),"Energy")
307 << " to "
308 << G4BestUnit(v->GetMaxEnergy(),"Energy")
309 << " in " << v->GetVectorLength()-1
310 << " bins " << G4endl;
311 break;
312 }
313 }
314 } else {
315 out << " Used LambdaPrime table of "
316 << particle->GetParticleName() << G4endl;
317 }
318 }
320 modelManager->DumpModelList(out, verboseLevel);
321
322 if(verboseLevel > 2 && buildLambdaTable) {
323 out << " LambdaTable address= " << theLambdaTable << G4endl;
324 if(theLambdaTable && particle == &part) {
325 out << (*theLambdaTable) << G4endl;
326 }
327 }
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
331
333{
334 // reset parameters for the new track
335 currentParticle = track->GetParticleDefinition();
338 preStepLambda = 0.0;
339
340 if(isIon) { massRatio = proton_mass_c2/currentParticle->GetPDGMass(); }
341
342 // forced biasing only for primary particles
343 if(biasManager) {
344 if(0 == track->GetParentID()) {
345 // primary particle
346 biasFlag = true;
348 }
349 }
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353
355 const G4Track& track,
356 G4double previousStepSize,
358{
360 G4double x = DBL_MAX;
361
364 const G4double scaledEnergy = preStepKinEnergy*massRatio;
365 SelectModel(scaledEnergy, currentCoupleIndex);
366 /*
367 G4cout << "PostStepGetPhysicalInteractionLength: idx= " << currentCoupleIndex
368 << " couple: " << currentCouple << G4endl;
369 */
370 if(!currentModel->IsActive(scaledEnergy)) {
374 preStepLambda = 0.0;
375 return x;
376 }
377
378 // forced biasing only for primary particles
379 if(biasManager) {
380 if(0 == track.GetParentID()) {
381 if(biasFlag &&
383 return biasManager->GetStepLimit((G4int)currentCoupleIndex, previousStepSize);
384 }
385 }
386 }
387
388 // compute mean free path
389
390 ComputeIntegralLambda(preStepKinEnergy, track);
391
392 // zero cross section
393 if(preStepLambda <= 0.0) {
396
397 } else {
398
399 // non-zero cross section
401
402 // beggining of tracking (or just after DoIt of this process)
405
406 } else {
407
409 previousStepSize/currentInteractionLength;
412 }
413
414 // new mean free path and step limit for the next step
417 }
418 return x;
419}
420
421//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
422
423void G4VEmProcess::ComputeIntegralLambda(G4double e, const G4Track& track)
424{
425 if (fXSType == fEmNoIntegral) {
426 preStepLambda = GetCurrentLambda(e, LogEkin(track));
427
428 } else if (fXSType == fEmIncreasing) {
429 if(e*invLambdaFactor < mfpKinEnergy) {
430 preStepLambda = GetCurrentLambda(e, LogEkin(track));
431 mfpKinEnergy = (preStepLambda > 0.0) ? e : 0.0;
432 }
433
434 } else if(fXSType == fEmDecreasing) {
435 if(e < mfpKinEnergy) {
436 const G4double e1 = e*lambdaFactor;
437 preStepLambda = GetCurrentLambda(e1);
438 mfpKinEnergy = e1;
439 }
440
441 } else if(fXSType == fEmOnePeak) {
442 const G4double epeak = (*theEnergyOfCrossSectionMax)[currentCoupleIndex];
443 if(e <= epeak) {
444 if(e*invLambdaFactor < mfpKinEnergy) {
445 preStepLambda = GetCurrentLambda(e, LogEkin(track));
446 mfpKinEnergy = (preStepLambda > 0.0) ? e : 0.0;
447 }
448 } else if(e < mfpKinEnergy) {
449 const G4double e1 = std::max(epeak, e*lambdaFactor);
450 preStepLambda = GetCurrentLambda(e1);
451 mfpKinEnergy = e1;
452 }
453 } else {
454 preStepLambda = GetCurrentLambda(e, LogEkin(track));
455 }
456}
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459
461 const G4Step& step)
462{
463 // clear number of interaction lengths in any case
466
468
469 // Do not make anything if particle is stopped, the annihilation then
470 // should be performed by the AtRestDoIt!
471 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
472
473 const G4double finalT = track.GetKineticEnergy();
474
475 // forced process - should happen only once per track
476 if(biasFlag) {
478 biasFlag = false;
479 }
480 }
481
482 // check active and select model
483 const G4double scaledEnergy = finalT*massRatio;
484 SelectModel(scaledEnergy, currentCoupleIndex);
485 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
486
487 // Integral approach
488 if (fXSType != fEmNoIntegral) {
489 const G4double logFinalT =
491 const G4double lx = std::max(GetCurrentLambda(finalT, logFinalT), 0.0);
492#ifdef G4VERBOSE
493 if(preStepLambda < lx && 1 < verboseLevel) {
494 G4cout << "WARNING: for " << currentParticle->GetParticleName()
495 << " and " << GetProcessName() << " E(MeV)= " << finalT/MeV
496 << " preLambda= " << preStepLambda
497 << " < " << lx << " (postLambda) " << G4endl;
498 }
499#endif
500 // if false interaction then use new cross section value
501 // if both values are zero - no interaction
502 if(preStepLambda*G4UniformRand() >= lx) {
503 return &fParticleChange;
504 }
505 }
506
507 // define new weight for primary and secondaries
509 if(weightFlag) {
510 weight /= biasFactor;
512 }
513
514#ifdef G4VERBOSE
515 if(1 < verboseLevel) {
516 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
517 << finalT/MeV
518 << " MeV; model= (" << currentModel->LowEnergyLimit()
519 << ", " << currentModel->HighEnergyLimit() << ")"
520 << G4endl;
521 }
522#endif
523
524 // sample secondaries
525 secParticles.clear();
526 currentModel->SampleSecondaries(&secParticles,
528 track.GetDynamicParticle(),
529 (*theCuts)[currentCoupleIndex]);
530
531 G4int num0 = (G4int)secParticles.size();
532
533 // splitting or Russian roulette
534 if(biasManager) {
536 G4double eloss = 0.0;
538 secParticles, track, currentModel, &fParticleChange, eloss,
540 step.GetPostStepPoint()->GetSafety());
541 if(eloss > 0.0) {
544 }
545 }
546 }
547
548 // save secondaries
549 G4int num = (G4int)secParticles.size();
550 if(num > 0) {
551
554 G4double time = track.GetGlobalTime();
555
556 G4int n1(0), n2(0);
557 if(num > mainSecondaries) {
558 currentModel->FillNumberOfSecondaries(n1, n2);
559 }
560
561 for (G4int i=0; i<num; ++i) {
563 if (nullptr != dp) {
565 G4double e = dp->GetKineticEnergy();
566 G4bool good = true;
567 if(applyCuts) {
568 if (p == theGamma) {
569 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
570
571 } else if (p == theElectron) {
572 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
573
574 } else if (p == thePositron) {
575 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
576 e < (*theCutsPositron)[currentCoupleIndex]) {
577 good = false;
578 e += 2.0*electron_mass_c2;
579 }
580 }
581 // added secondary if it is good
582 }
583 if (good) {
584 G4Track* t = new G4Track(dp, time, track.GetPosition());
586 if (biasManager) {
587 t->SetWeight(weight * biasManager->GetWeight(i));
588 } else {
589 t->SetWeight(weight);
590 }
592
593 // define type of secondary
594 if(i < mainSecondaries) {
596 if(GetProcessSubType() == fComptonScattering && p == theGamma) {
598 }
599 } else if(i < mainSecondaries + n1) {
601 } else if(i < mainSecondaries + n1 + n2) {
603 } else {
604 if(i < num0) {
605 if(p == theGamma) {
607 } else {
609 }
610 } else {
612 }
613 }
614 /*
615 G4cout << "Secondary(post step) has weight " << t->GetWeight()
616 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
617 << GetProcessName() << " fluoID= " << fluoID
618 << " augerID= " << augerID <<G4endl;
619 */
620 } else {
621 delete dp;
622 edep += e;
623 }
624 }
625 }
627 }
628
631 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
634 }
635
636 return &fParticleChange;
637}
638
639//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
640
642 const G4String& directory,
643 G4bool ascii)
644{
645 if(!isTheMaster || part != particle) { return true; }
646 if(G4EmTableUtil::StoreTable(this, part, theLambdaTable,
647 directory, "Lambda",
648 verboseLevel, ascii) &&
649 G4EmTableUtil::StoreTable(this, part, theLambdaTablePrim,
650 directory, "LambdaPrim",
651 verboseLevel, ascii)) {
652 return true;
653 }
654 return false;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
658
660 const G4String& dir,
661 G4bool ascii)
662{
663 if(!isTheMaster || part != particle) { return true; }
664 G4bool yes = true;
665 if(buildLambdaTable) {
666 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTable, dir,
667 "Lambda", verboseLevel,
668 ascii, splineFlag);
669 }
670 if(yes && minKinEnergyPrim < maxKinEnergy) {
671 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTablePrim, dir,
672 "LambdaPrim", verboseLevel,
673 ascii, splineFlag);
674 }
675 return yes;
676}
677
678//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
679
681 const G4MaterialCutsCouple* couple)
682{
683 CurrentSetup(couple, kinEnergy);
684 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
685}
686
687//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
688
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
701 G4double Z, G4double A, G4double cut)
702{
704 return (currentModel) ?
705 currentModel->ComputeCrossSectionPerAtom(currentParticle, kinEnergy,
706 Z, A, cut) : 0.0;
707}
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710
713{
714 DefineMaterial(couple);
715 G4PhysicsVector* newv = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy,
716 nLambdaBins, splineFlag);
717 return newv;
718}
719
720//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
721
723{
724 return (nullptr != currentModel) ?
725 currentModel->GetCurrentElement(currentMaterial) : nullptr;
726}
727
728//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
729
731{
732 return (nullptr != currentModel) ?
733 currentModel->GetCurrentElement(currentMaterial) : nullptr;
734}
735
736//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
737
739{
740 return (nullptr != currentModel) ?
741 currentModel->GetCurrentIsotope(GetCurrentElement()) : nullptr;
742}
743
744//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745
747{
748 if(f > 0.0) {
749 biasFactor = f;
750 weightFlag = flag;
751 if(1 < verboseLevel) {
752 G4cout << "### SetCrossSectionBiasingFactor: for "
753 << particle->GetParticleName()
754 << " and process " << GetProcessName()
755 << " biasFactor= " << f << " weightFlag= " << flag
756 << G4endl;
757 }
758 }
759}
760
761//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
762
763void
765 G4bool flag)
766{
767 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
768 if(1 < verboseLevel) {
769 G4cout << "### ActivateForcedInteraction: for "
770 << particle->GetParticleName()
771 << " and process " << GetProcessName()
772 << " length(mm)= " << length/mm
773 << " in G4Region <" << r
774 << "> weightFlag= " << flag
775 << G4endl;
776 }
777 weightFlag = flag;
779}
780
781//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
782
783void
785 G4double factor,
786 G4double energyLimit)
787{
788 if (0.0 <= factor) {
789
790 // Range cut can be applied only for e-
791 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
792 { return; }
793
795 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
796 if(1 < verboseLevel) {
797 G4cout << "### ActivateSecondaryBiasing: for "
798 << " process " << GetProcessName()
799 << " factor= " << factor
800 << " in G4Region <" << region
801 << "> energyLimit(MeV)= " << energyLimit/MeV
802 << G4endl;
803 }
804 }
805}
806
807//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808
810{
811 if(5 < n && n < 10000000) {
812 nLambdaBins = n;
813 actBinning = true;
814 } else {
815 G4double e = (G4double)n;
816 PrintWarning("SetLambdaBinning", e);
817 }
818}
819
820//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
821
823{
824 if(1.e-3*eV < e && e < maxKinEnergy) {
825 nLambdaBins = G4lrint(nLambdaBins*G4Log(maxKinEnergy/e)
826 /G4Log(maxKinEnergy/minKinEnergy));
827 minKinEnergy = e;
828 actMinKinEnergy = true;
829 } else { PrintWarning("SetMinKinEnergy", e); }
830}
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
833
835{
836 if(minKinEnergy < e && e < 1.e+6*TeV) {
837 nLambdaBins = G4lrint(nLambdaBins*G4Log(e/minKinEnergy)
838 /G4Log(maxKinEnergy/minKinEnergy));
839 maxKinEnergy = e;
840 actMaxKinEnergy = true;
841 } else { PrintWarning("SetMaxKinEnergy", e); }
842}
843
844//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
845
847{
848 if(theParameters->MinKinEnergy() <= e &&
849 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
850 else { PrintWarning("SetMinKinEnergyPrim", e); }
851}
852
853//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
854
856{
857 return (nam == GetProcessName()) ? this : nullptr;
858}
859
860//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
861
863{
864 return theParameters->MscThetaLimit();
865}
866
867//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
868
869void G4VEmProcess::PrintWarning(G4String tit, G4double val)
870{
871 G4String ss = "G4VEmProcess::" + tit;
873 ed << "Parameter is out of range: " << val
874 << " it will have no effect!\n" << " Process "
875 << GetProcessName() << " nbins= " << theParameters->NumberOfBins()
876 << " Emin(keV)= " << theParameters->MinKinEnergy()/keV
877 << " Emax(GeV)= " << theParameters->MaxKinEnergy()/GeV;
878 G4Exception(ss, "em0044", JustWarning, ed);
879}
880
881//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
882
883void G4VEmProcess::ProcessDescription(std::ostream& out) const
884{
885 if(nullptr != particle) {
886 StreamInfo(out, *particle, true);
887 }
888}
889
890//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fComptonScattering
@ fEmOnePeak
@ fEmDecreasing
@ fEmNoIntegral
@ fEmIncreasing
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition G4Log.hh:227
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
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
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:91
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 *fm, const G4Region *r)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
void SetFluoFlag(G4bool val)
G4bool IsPrintLocked() const
static G4EmParameters * Instance()
G4int NumberOfBins() const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MscThetaLimit() const
void DefineRegParamForEM(G4VEmProcess *) const
G4bool ApplyCuts() const
G4double MaxKinEnergy() const
G4bool Integral() const
G4double LambdaFactor() const
static void BuildEmProcess(G4VEmProcess *proc, const G4VEmProcess *masterProc, const G4ParticleDefinition *firstPart, const G4ParticleDefinition *part, const G4int nModels, const G4int verb, const G4bool master, const G4bool isLocked, const G4bool toBuild, G4bool &baseMat)
static G4bool RetrieveTable(G4VProcess *ptr, const G4ParticleDefinition *part, G4PhysicsTable *aTable, const G4String &dir, const G4String &tname, const G4int verb, const G4bool ascii, const G4bool spline)
static const G4DataVector * PrepareEmProcess(G4VEmProcess *proc, const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4EmModelManager *modelManager, const G4double &maxKinEnergy, G4int &secID, G4int &tripletID, G4int &mainSec, const G4int &verb, const G4bool &master)
static G4bool StoreTable(G4VProcess *, const G4ParticleDefinition *, G4PhysicsTable *, const G4String &dir, const G4String &tname, G4int verb, G4bool ascii)
static void BuildLambdaTable(G4VEmProcess *proc, const G4ParticleDefinition *part, G4EmModelManager *modelManager, G4LossTableBuilder *bld, G4PhysicsTable *theLambdaTable, G4PhysicsTable *theLambdaTablePrim, const G4double minKinEnergy, const G4double minKinEnergyPrim, const G4double maxKinEnergy, const G4double scale, const G4int verbose, const G4bool startFromNull, const G4bool splineFlag)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4VAtomDeexcitation * AtomDeexcitation()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void InitializeForPostStep(const G4Track &)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static G4Positron * Positron()
Definition G4Positron.cc:90
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetSafety() const
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
void SetCreatorModelID(const G4int id)
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4double LowEnergyLimit() const
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
G4double HighEnergyLimit() const
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
G4bool IsActive(G4double kinEnergy) const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4double MeanFreePath(const G4Track &track)
G4int mainSecondaries
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
virtual void StreamProcessInfo(std::ostream &) const
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
G4EmBiasingManager * biasManager
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
std::vector< G4double > * theEnergyOfCrossSectionMax
void SetMinKinEnergy(G4double e)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void StartTracking(G4Track *) override
G4double GetCrossSection(const G4double kinEnergy, const G4MaterialCutsCouple *couple) 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 StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
~G4VEmProcess() override
void SetLambdaBinning(G4int nbins)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void ProcessDescription(std::ostream &outFile) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4Element * GetTargetElement() const
virtual G4VEmProcess * GetEmProcess(const G4String &name)
G4double MaxKinEnergy() const
const G4Isotope * GetTargetIsotope() const
std::vector< G4DynamicParticle * > secParticles
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
const G4MaterialCutsCouple * currentCouple
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4double preStepKinEnergy
G4double PolarAngleLimit() const
size_t currentCoupleIndex
G4ParticleChangeForGamma fParticleChange
void SetParticle(const G4ParticleDefinition *p)
void SetMinKinEnergyPrim(G4double e)
void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMaxKinEnergy(G4double e)
const G4Material * currentMaterial
const G4Element * GetCurrentElement() const
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void BuildLambdaTable()
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
G4double theInitialNumberOfInteractionLength
void SetVerboseLevel(G4int value)
const G4VProcess * GetMasterProcess() const
G4int verboseLevel
G4double theNumberOfInteractionLengthLeft
G4VParticleChange * pParticleChange
G4int GetProcessSubType() const
const G4String & GetProcessName() const
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62