Geant4 11.4.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 "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
98 fParticleChange.SetSecondaryWeightByProcess(true);
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) {
116 }
117 delete modelManager;
118 delete biasManager;
119 lManager->DeRegister(this);
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
125 const G4Region* region)
126{
127 if(nullptr == ptr) { return; }
128 G4VEmFluctuationModel* fm = nullptr;
129 modelManager->AddEmModel(order, ptr, fm, region);
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134
136{
137 if(nullptr == ptr) { return; }
138 if(!emModels.empty()) {
139 for(auto & em : emModels) { if(em == ptr) { return; } }
140 }
141 emModels.push_back(ptr);
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145
147{
148 if(nullptr == particle) { SetParticle(&part); }
149
150 if(part.GetParticleType() == "nucleus" &&
151 part.GetParticleSubType() == "generic") {
152
153 G4String pname = part.GetParticleName();
154 if(pname != "deuteron" && pname != "triton" &&
155 pname != "He3" && pname != "alpha" && pname != "alpha+" &&
156 pname != "helium" && pname != "hydrogen") {
157
158 particle = G4GenericIon::GenericIon();
159 isIon = true;
160 }
161 }
162 if(particle != &part) { return; }
163
164 lManager->PreparePhysicsTable(&part, this);
165
166 // for new run
167 currentCouple = nullptr;
168 preStepLambda = 0.0;
169 fLambdaEnergy = 0.0;
170
171 InitialiseProcess(particle);
172
173 G4LossTableBuilder* bld = lManager->GetTableBuilder();
174 const G4ProductionCutsTable* theCoupleTable=
176 theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
177 theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
178 theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
179
180 // initialisation of the process
181 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
182 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
183
184 applyCuts = theParameters->ApplyCuts();
185 lambdaFactor = theParameters->LambdaFactor();
186 invLambdaFactor = 1.0/lambdaFactor;
187 theParameters->DefineRegParamForEM(this);
188
189 // integral option may be disabled
190 if(!theParameters->Integral()) { fXSType = fEmNoIntegral; }
191
192 // prepare tables
193 if(isTheMaster) {
194 if(buildLambdaTable) {
195 if (nullptr == theData) { theData = new G4EmDataHandler(2, particle->GetParticleName()); }
196 theLambdaTable = theData->MakeTable(theLambdaTable, 0);
197 bld->InitialiseBaseMaterials(theLambdaTable);
198 }
199 // high energy table
200 if(minKinEnergyPrim < maxKinEnergy) {
201 if (nullptr == theData) { theData = new G4EmDataHandler(2, particle->GetParticleName()); }
202 theLambdaTablePrim = theData->MakeTable(theLambdaTablePrim, 1);
203 bld->InitialiseBaseMaterials(theLambdaTablePrim);
204 }
205 }
206 // models
208 numberOfModels = modelManager->NumberOfModels();
209 currentModel = modelManager->GetModel(0);
210 if(nullptr != lManager->AtomDeexcitation()) {
211 modelManager->SetFluoFlag(true);
212 }
213 // forced biasing
214 if(nullptr != biasManager) {
215 biasManager->Initialise(part, GetProcessName(), verboseLevel);
216 biasFlag = false;
217 }
218
219 theCuts =
220 G4EmTableUtil::PrepareEmProcess(this, particle, secondaryParticle,
221 modelManager, maxKinEnergy,
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
229{
230 if(nullptr == masterProc) {
231 if(isTheMaster) { masterProc = this; }
232 else { masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());}
233 }
234 G4int nModels = modelManager->NumberOfModels();
235 G4bool isLocked = theParameters->IsPrintLocked();
236 G4bool toBuild = (buildLambdaTable || minKinEnergyPrim < maxKinEnergy);
237
238 G4EmTableUtil::BuildEmProcess(this, masterProc, particle, &part,
239 nModels, verboseLevel, isTheMaster,
240 isLocked, toBuild, baseMat);
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244
246{
247 G4double scale = theParameters->MaxKinEnergy()/theParameters->MinKinEnergy();
248 G4int nbin =
249 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
250 if(actBinning) { nbin = std::max(nbin, nLambdaBins); }
251 scale = nbin/G4Log(scale);
252
253 G4LossTableBuilder* bld = lManager->GetTableBuilder();
254 G4EmTableUtil::BuildLambdaTable(this, particle, modelManager,
255 bld, theLambdaTable, theLambdaTablePrim,
256 minKinEnergy, minKinEnergyPrim,
257 maxKinEnergy, scale, verboseLevel,
258 startFromNull, splineFlag);
259}
260
261//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
262
263void G4VEmProcess::StreamInfo(std::ostream& out,
264 const G4ParticleDefinition& part, G4bool rst) const
265{
266 G4String indent = (rst ? " " : "");
267 out << std::setprecision(6);
268 out << G4endl << indent << GetProcessName() << ": ";
269 if (!rst) {
270 out << " for " << part.GetParticleName();
271 }
272 if(fXSType != fEmNoIntegral) { out << " XStype:" << fXSType; }
273 if(applyCuts) { out << " applyCuts:1 "; }
274 G4int subtype = GetProcessSubType();
275 out << " SubType=" << subtype;
276 if (subtype == fAnnihilation) {
277 G4int mod = theParameters->PositronAtRestModelType();
278 const G4String namp[2] = {"Simple", "Allison"};
279 out << " AtRestModel:" << namp[mod];
280 }
281 if(biasFactor != 1.0) { out << " BiasingFactor=" << biasFactor; }
282 out << " BuildTable=" << buildLambdaTable << G4endl;
283 if(buildLambdaTable) {
284 if(particle == &part) {
285 for(auto & v : *theLambdaTable) {
286 if(nullptr != v) {
287 out << " Lambda table from ";
288 G4double emin = v->Energy(0);
289 G4double emax = v->GetMaxEnergy();
290 G4int nbin = G4int(v->GetVectorLength() - 1);
291 if(emin > minKinEnergy) { out << "threshold "; }
292 else { out << G4BestUnit(emin,"Energy"); }
293 out << " to "
294 << G4BestUnit(emax,"Energy")
295 << ", " << G4lrint(nbin/std::log10(emax/emin))
296 << " bins/decade, spline: "
297 << splineFlag << G4endl;
298 break;
299 }
300 }
301 } else {
302 out << " Used Lambda table of "
303 << particle->GetParticleName() << G4endl;
304 }
305 }
306 if(minKinEnergyPrim < maxKinEnergy) {
307 if(particle == &part) {
308 for(auto & v : *theLambdaTablePrim) {
309 if(nullptr != v) {
310 out << " LambdaPrime table from "
311 << G4BestUnit(v->Energy(0),"Energy")
312 << " to "
313 << G4BestUnit(v->GetMaxEnergy(),"Energy")
314 << " in " << v->GetVectorLength()-1
315 << " bins " << G4endl;
316 break;
317 }
318 }
319 } else {
320 out << " Used LambdaPrime table of "
321 << particle->GetParticleName() << G4endl;
322 }
323 }
325 modelManager->DumpModelList(out, verboseLevel);
326
327 if(verboseLevel > 2 && buildLambdaTable) {
328 out << " LambdaTable address= " << theLambdaTable << G4endl;
329 if(theLambdaTable && particle == &part) {
330 out << (*theLambdaTable) << G4endl;
331 }
332 }
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
336
338{
339 // reset parameters for the new track
340 currentParticle = track->GetParticleDefinition();
343 preStepLambda = 0.0;
344
345 if(isIon) { massRatio = proton_mass_c2/currentParticle->GetPDGMass(); }
346
347 // forced biasing only for primary particles
348 if(biasManager) {
349 if(0 == track->GetParentID()) {
350 // primary particle
351 biasFlag = true;
352 biasManager->ResetForcedInteraction();
353 }
354 }
355 for (G4int i=0; i<numberOfModels; ++i) {
356 auto ptr = GetModelByIndex(i);
357 ptr->StartTracking(track);
358 }
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
362
364 const G4Track& track,
365 G4double previousStepSize,
367{
369 G4double x = DBL_MAX;
370
373 const G4double scaledEnergy = preStepKinEnergy*massRatio;
374 SelectModel(scaledEnergy, currentCoupleIndex);
375
376 // In models applied to ions the dynamic charge is needed
377 if (isIon) { currentModel->ChargeSquareRatio(track); }
378 /*
379 G4cout << "PostStepGetPhysicalInteractionLength: idx= " << currentCoupleIndex
380 << " couple: " << currentCouple << G4endl;
381 */
382 if(!currentModel->IsActive(scaledEnergy)) {
386 preStepLambda = 0.0;
387 return x;
388 }
389
390 // forced biasing only for primary particles
391 if(biasManager) {
392 if(0 == track.GetParentID()) {
393 if(biasFlag &&
394 biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
395 return biasManager->GetStepLimit((G4int)currentCoupleIndex, previousStepSize);
396 }
397 }
398 }
399
400 // compute mean free path
401
402 ComputeIntegralLambda(preStepKinEnergy, track);
403
404 // zero cross section
405 if(preStepLambda <= 0.0) {
408
409 } else {
410
411 // non-zero cross section
413
414 // beggining of tracking (or just after DoIt of this process)
417
418 } else {
419
421 previousStepSize/currentInteractionLength;
424 }
425
426 // new mean free path and step limit for the next step
429 }
430 return x;
431}
432
433//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
434
435void G4VEmProcess::ComputeIntegralLambda(G4double e, const G4Track& track)
436{
437 if (fXSType == fEmNoIntegral) {
438 preStepLambda = GetCurrentLambda(e, LogEkin(track));
439
440 } else if (fXSType == fEmIncreasing) {
441 if(e*invLambdaFactor < mfpKinEnergy) {
442 preStepLambda = GetCurrentLambda(e, LogEkin(track));
443 mfpKinEnergy = (preStepLambda > 0.0) ? e : 0.0;
444 }
445
446 } else if(fXSType == fEmDecreasing) {
447 if(e < mfpKinEnergy) {
448 const G4double e1 = e*lambdaFactor;
449 preStepLambda = GetCurrentLambda(e1);
450 mfpKinEnergy = e1;
451 }
452
453 } else if(fXSType == fEmOnePeak) {
454 const G4double epeak = (*theEnergyOfCrossSectionMax)[currentCoupleIndex];
455 if(e <= epeak) {
456 if(e*invLambdaFactor < mfpKinEnergy) {
457 preStepLambda = GetCurrentLambda(e, LogEkin(track));
458 mfpKinEnergy = (preStepLambda > 0.0) ? e : 0.0;
459 }
460 } else if(e < mfpKinEnergy) {
461 const G4double e1 = std::max(epeak, e*lambdaFactor);
462 preStepLambda = GetCurrentLambda(e1);
463 mfpKinEnergy = e1;
464 }
465 } else {
466 preStepLambda = GetCurrentLambda(e, LogEkin(track));
467 }
468}
469
470//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
471
473 const G4Step& step)
474{
475 // clear number of interaction lengths in any case
478
479 fParticleChange.InitializeForPostStep(track);
480
481 // Do not make anything if particle is stopped, the annihilation then
482 // should be performed by the AtRestDoIt!
483 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
484
485 const G4double finalT = track.GetKineticEnergy();
486
487 // forced process - should happen only once per track
488 if(biasFlag) {
489 if(biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
490 biasFlag = false;
491 }
492 }
493
494 // check active and select model
495 const G4double scaledEnergy = finalT*massRatio;
496 SelectModel(scaledEnergy, currentCoupleIndex);
497 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
498
499 // Integral approach
500 if (fXSType != fEmNoIntegral) {
501 const G4double logFinalT =
503 const G4double lx = std::max(GetCurrentLambda(finalT, logFinalT), 0.0);
504#ifdef G4VERBOSE
505 if(preStepLambda < lx && 1 < verboseLevel) {
506 G4cout << "WARNING: for " << currentParticle->GetParticleName()
507 << " and " << GetProcessName() << " E(MeV)= " << finalT/MeV
508 << " preLambda= " << preStepLambda
509 << " < " << lx << " (postLambda) " << G4endl;
510 }
511#endif
512 // if false interaction then use new cross section value
513 // if both values are zero - no interaction
514 if(preStepLambda*G4UniformRand() >= lx) {
515 return &fParticleChange;
516 }
517 }
518
519 // define new weight for primary and secondaries
520 G4double weight = fParticleChange.GetParentWeight();
521 if(weightFlag) {
522 weight /= biasFactor;
523 fParticleChange.ProposeWeight(weight);
524 }
525
526#ifdef G4VERBOSE
527 if(1 < verboseLevel) {
528 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
529 << finalT/MeV
530 << " MeV; model= (" << currentModel->LowEnergyLimit()
531 << ", " << currentModel->HighEnergyLimit() << ")"
532 << G4endl;
533 }
534#endif
535
536 // sample secondaries
537 secParticles.clear();
538 currentModel->SampleSecondaries(&secParticles,
540 track.GetDynamicParticle(),
541 (*theCuts)[currentCoupleIndex]);
542
543 G4int num0 = (G4int)secParticles.size();
544
545 // splitting or Russian roulette
546 if(biasManager) {
547 if(biasManager->SecondaryBiasingRegion((G4int)currentCoupleIndex)) {
548 G4double eloss = 0.0;
549 weight *= biasManager->ApplySecondaryBiasing(
550 secParticles, track, currentModel, &fParticleChange, eloss,
552 step.GetPostStepPoint()->GetSafety());
553 if(eloss > 0.0) {
554 eloss += fParticleChange.GetLocalEnergyDeposit();
555 fParticleChange.ProposeLocalEnergyDeposit(eloss);
556 }
557 }
558 }
559
560 // save secondaries
561 G4int num = (G4int)secParticles.size();
562 if(num > 0) {
563
564 fParticleChange.SetNumberOfSecondaries(num);
565 G4double edep = fParticleChange.GetLocalEnergyDeposit();
566 G4double time = track.GetGlobalTime();
567
568 G4int n1(0), n2(0);
569 if(num0 > mainSecondaries) {
570 currentModel->FillNumberOfSecondaries(n1, n2);
571 }
572
573 for (G4int i=0; i<num; ++i) {
575 if (nullptr != dp) {
577 G4double e = dp->GetKineticEnergy();
578 G4bool good = true;
579 if(applyCuts) {
580 if (p == theGamma) {
581 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
582
583 } else if (p == theElectron) {
584 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
585
586 } else if (p == thePositron) {
587 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
588 e < (*theCutsPositron)[currentCoupleIndex]) {
589 good = false;
590 e += 2.0*electron_mass_c2;
591 }
592 }
593 // added secondary if it is good
594 }
595 if (good) {
596 G4Track* t = new G4Track(dp, time, track.GetPosition());
598 if (biasManager) {
599 t->SetWeight(weight * biasManager->GetWeight(i));
600 } else {
601 t->SetWeight(weight);
602 }
603 pParticleChange->AddSecondary(t);
604
605 // define type of secondary
606 if(i < mainSecondaries) {
608 if(GetProcessSubType() == fComptonScattering && p == theGamma) {
610 }
611 } else if(i < mainSecondaries + n1) {
613 } else if(i < mainSecondaries + n1 + n2) {
615 } else {
616 if(i < num0) {
617 if(p == theGamma) {
619 } else {
621 }
622 } else {
624 }
625 }
626 /*
627 G4cout << "Secondary(post step) has weight " << t->GetWeight()
628 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
629 << GetProcessName() << " fluoID= " << fluoID
630 << " augerID= " << augerID <<G4endl;
631 */
632 } else {
633 delete dp;
634 edep += e;
635 }
636 }
637 }
638 fParticleChange.ProposeLocalEnergyDeposit(edep);
639 }
640
641 if(0.0 == fParticleChange.GetProposedKineticEnergy() &&
642 fAlive == fParticleChange.GetTrackStatus()) {
643 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
644 { fParticleChange.ProposeTrackStatus(fStopButAlive); }
645 else { fParticleChange.ProposeTrackStatus(fStopAndKill); }
646 }
647
648 return &fParticleChange;
649}
650
651//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
652
654 const G4String& directory,
655 G4bool ascii)
656{
657 if(!isTheMaster || part != particle) { return true; }
658 if(G4EmTableUtil::StoreTable(this, part, theLambdaTable,
659 directory, "Lambda",
660 verboseLevel, ascii) &&
661 G4EmTableUtil::StoreTable(this, part, theLambdaTablePrim,
662 directory, "LambdaPrim",
663 verboseLevel, ascii)) {
664 return true;
665 }
666 return false;
667}
668
669//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
670
672 const G4String& dir,
673 G4bool ascii)
674{
675 if(!isTheMaster || part != particle) { return true; }
676 G4bool yes = true;
677 if(buildLambdaTable) {
678 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTable, dir,
679 "Lambda", verboseLevel,
680 ascii, splineFlag);
681 }
682 if(yes && minKinEnergyPrim < maxKinEnergy) {
683 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTablePrim, dir,
684 "LambdaPrim", verboseLevel,
685 ascii, splineFlag);
686 }
687 return yes;
688}
689
690//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691
693 const G4MaterialCutsCouple* couple)
694{
695 CurrentSetup(couple, kinEnergy);
696 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
697}
698
699//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
700
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710
713 G4double Z, G4double A, G4double cut)
714{
716 return (currentModel) ?
717 currentModel->ComputeCrossSectionPerAtom(currentParticle, kinEnergy,
718 Z, A, cut) : 0.0;
719}
720
721//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
722
725{
726 DefineMaterial(couple);
727 G4PhysicsVector* newv = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy,
728 nLambdaBins, splineFlag);
729 return newv;
730}
731
732//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
733
735{
736 return (nullptr != currentModel) ?
737 currentModel->GetCurrentElement(currentMaterial) : nullptr;
738}
739
740//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
741
743{
744 return (nullptr != currentModel) ?
745 currentModel->GetCurrentElement(currentMaterial) : nullptr;
746}
747
748//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
749
751{
752 return (nullptr != currentModel) ?
753 currentModel->GetCurrentIsotope(GetCurrentElement()) : nullptr;
754}
755
756//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
757
759{
760 if(f > 0.0) {
761 biasFactor = f;
762 weightFlag = flag;
763 if(1 < verboseLevel) {
764 G4cout << "### SetCrossSectionBiasingFactor: for "
765 << particle->GetParticleName()
766 << " and process " << GetProcessName()
767 << " biasFactor= " << f << " weightFlag= " << flag
768 << G4endl;
769 }
770 }
771}
772
773//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
774
775void
777 G4bool flag)
778{
779 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
780 if(1 < verboseLevel) {
781 G4cout << "### ActivateForcedInteraction: for "
782 << particle->GetParticleName()
783 << " and process " << GetProcessName()
784 << " length(mm)= " << length/mm
785 << " in G4Region <" << r
786 << "> weightFlag= " << flag
787 << G4endl;
788 }
789 weightFlag = flag;
790 biasManager->ActivateForcedInteraction(length, r);
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
795void
797 G4double factor,
798 G4double energyLimit)
799{
800 if (0.0 <= factor) {
801
802 // Range cut can be applied only for e-
803 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
804 { return; }
805
807 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
808 if(1 < verboseLevel) {
809 G4cout << "### ActivateSecondaryBiasing: for "
810 << " process " << GetProcessName()
811 << " factor= " << factor
812 << " in G4Region <" << region
813 << "> energyLimit(MeV)= " << energyLimit/MeV
814 << G4endl;
815 }
816 }
817}
818
819//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
820
822{
823 if(5 < n && n < 10000000) {
824 nLambdaBins = n;
825 actBinning = true;
826 } else {
827 G4double e = (G4double)n;
828 PrintWarning("SetLambdaBinning", e);
829 }
830}
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
833
835{
836 if(1.e-3*eV < e && e < maxKinEnergy) {
837 nLambdaBins = G4lrint(nLambdaBins*G4Log(maxKinEnergy/e)
838 /G4Log(maxKinEnergy/minKinEnergy));
839 minKinEnergy = e;
840 actMinKinEnergy = true;
841 } else { PrintWarning("SetMinKinEnergy", e); }
842}
843
844//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
845
847{
848 if(minKinEnergy < e && e < 1.e+6*TeV) {
849 nLambdaBins = G4lrint(nLambdaBins*G4Log(e/minKinEnergy)
850 /G4Log(maxKinEnergy/minKinEnergy));
851 maxKinEnergy = e;
852 actMaxKinEnergy = true;
853 } else { PrintWarning("SetMaxKinEnergy", e); }
854}
855
856//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
857
859{
860 if(theParameters->MinKinEnergy() <= e &&
861 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
862 else { PrintWarning("SetMinKinEnergyPrim", e); }
863}
864
865//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
866
868{
869 return (nam == GetProcessName()) ? this : nullptr;
870}
871
872//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
873
875{
876 return theParameters->MscThetaLimit();
877}
878
879//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
880
881void G4VEmProcess::PrintWarning(G4String tit, G4double val)
882{
883 G4String ss = "G4VEmProcess::" + tit;
885 ed << "Parameter is out of range: " << val
886 << " it will have no effect!\n" << " Process "
887 << GetProcessName() << " nbins= " << theParameters->NumberOfBins()
888 << " Emin(keV)= " << theParameters->MinKinEnergy()/keV
889 << " Emax(GeV)= " << theParameters->MaxKinEnergy()/GeV;
890 G4Exception(ss, "em0044", JustWarning, ed);
891}
892
893//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
894
895void G4VEmProcess::ProcessDescription(std::ostream& out) const
896{
897 if(nullptr != particle) {
898 StreamInfo(out, *particle, true);
899 }
900}
901
902//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fComptonScattering
@ fAnnihilation
@ 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:169
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
static G4EmParameters * Instance()
G4int NumberOfBins() const
G4double MinKinEnergy() const
G4double MaxKinEnergy() 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()
static G4bool GetBaseMaterialFlag()
static const std::vector< G4double > * GetDensityFactors()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
static const std::vector< G4int > * GetCoupleIndexes()
static G4LossTableManager * Instance()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static G4Positron * Positron()
Definition G4Positron.cc:90
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
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
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
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)
virtual void StreamProcessInfo(std::ostream &) const
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
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
std::size_t currentCoupleIndex
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
G4ParticleChangeForGamma fParticleChange
G4VEmModel * SelectModel(G4double kinEnergy, std::size_t)
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 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