Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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, isTheMaster);
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
339 if(isIon) { massRatio = proton_mass_c2/currentParticle->GetPDGMass(); }
340
341 // forced biasing only for primary particles
342 if(biasManager) {
343 if(0 == track->GetParentID()) {
344 // primary particle
345 biasFlag = true;
347 }
348 }
349}
350
351//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
352
354 const G4Track& track,
355 G4double previousStepSize,
357{
359 G4double x = DBL_MAX;
360
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)) {
373 return x;
374 }
375
376 // forced biasing only for primary particles
377 if(biasManager) {
378 if(0 == track.GetParentID()) {
379 if(biasFlag &&
381 return biasManager->GetStepLimit((G4int)currentCoupleIndex, previousStepSize);
382 }
383 }
384 }
385
386 // compute mean free path
387 ComputeIntegralLambda(preStepKinEnergy, preStepLogKinEnergy);
388
389 // zero cross section
390 if(preStepLambda <= 0.0) {
393
394 } else {
395
396 // non-zero cross section
398
399 // beggining of tracking (or just after DoIt of this process)
402
403 } else {
404
406 previousStepSize/currentInteractionLength;
409 }
410
411 // new mean free path and step limit for the next step
414 }
415 return x;
416}
417
418//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
419
420void G4VEmProcess::ComputeIntegralLambda(G4double e, G4double loge)
421{
422 if(fXSType == fEmNoIntegral) {
423 preStepLambda = GetCurrentLambda(e, loge);
424
425 } else if(fXSType == fEmIncreasing) {
426 if(e*invLambdaFactor < mfpKinEnergy) {
427 mfpKinEnergy = e;
428 preStepLambda = GetCurrentLambda(e, loge);
429 }
430
431 } else if(fXSType == fEmDecreasing) {
432 if(e < mfpKinEnergy) {
433 const G4double e1 = e*lambdaFactor;
434 preStepLambda = GetCurrentLambda(e1);
435 mfpKinEnergy = e1;
436 }
437
438 } else if(fXSType == fEmOnePeak) {
439 const G4double epeak = (*theEnergyOfCrossSectionMax)[currentCoupleIndex];
440 if(e <= epeak) {
441 if(e*invLambdaFactor < mfpKinEnergy) {
442 mfpKinEnergy = e;
443 preStepLambda = GetCurrentLambda(e, loge);
444 }
445 } else if(e < mfpKinEnergy) {
446 const G4double e1 = std::max(epeak, e*lambdaFactor);
447 preStepLambda = GetCurrentLambda(e1);
448 mfpKinEnergy = e1;
449 }
450
451 } else {
452 preStepLambda = GetCurrentLambda(e, loge);
453 }
454}
455
456//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
457
459 const G4Step& step)
460{
461 // clear number of interaction lengths in any case
464
466
467 // Do not make anything if particle is stopped, the annihilation then
468 // should be performed by the AtRestDoIt!
469 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
470
471 const G4double finalT = track.GetKineticEnergy();
472
473 // forced process - should happen only once per track
474 if(biasFlag) {
476 biasFlag = false;
477 }
478 }
479
480 // check active and select model
481 const G4double scaledEnergy = finalT*massRatio;
482 SelectModel(scaledEnergy, currentCoupleIndex);
483 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
484
485 // Integral approach
486 if (fXSType != fEmNoIntegral) {
487 const G4double logFinalT =
489 const G4double lx = std::max(GetCurrentLambda(finalT, logFinalT), 0.0);
490#ifdef G4VERBOSE
491 if(preStepLambda < lx && 1 < verboseLevel) {
492 G4cout << "WARNING: for " << currentParticle->GetParticleName()
493 << " and " << GetProcessName() << " E(MeV)= " << finalT/MeV
494 << " preLambda= " << preStepLambda
495 << " < " << lx << " (postLambda) " << G4endl;
496 }
497#endif
498 // if false interaction then use new cross section value
499 // if both values are zero - no interaction
500 if(preStepLambda*G4UniformRand() >= lx) {
501 return &fParticleChange;
502 }
503 }
504
505 // define new weight for primary and secondaries
507 if(weightFlag) {
508 weight /= biasFactor;
510 }
511
512#ifdef G4VERBOSE
513 if(1 < verboseLevel) {
514 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
515 << finalT/MeV
516 << " MeV; model= (" << currentModel->LowEnergyLimit()
517 << ", " << currentModel->HighEnergyLimit() << ")"
518 << G4endl;
519 }
520#endif
521
522 // sample secondaries
523 secParticles.clear();
524 currentModel->SampleSecondaries(&secParticles,
526 track.GetDynamicParticle(),
527 (*theCuts)[currentCoupleIndex]);
528
529 G4int num0 = (G4int)secParticles.size();
530
531 // splitting or Russian roulette
532 if(biasManager) {
534 G4double eloss = 0.0;
536 secParticles, track, currentModel, &fParticleChange, eloss,
538 step.GetPostStepPoint()->GetSafety());
539 if(eloss > 0.0) {
542 }
543 }
544 }
545
546 // save secondaries
547 G4int num = (G4int)secParticles.size();
548 if(num > 0) {
549
552 G4double time = track.GetGlobalTime();
553
554 G4int n1(0), n2(0);
555 if(num > mainSecondaries) {
556 currentModel->FillNumberOfSecondaries(n1, n2);
557 }
558
559 for (G4int i=0; i<num; ++i) {
561 if (nullptr != dp) {
563 G4double e = dp->GetKineticEnergy();
564 G4bool good = true;
565 if(applyCuts) {
566 if (p == theGamma) {
567 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
568
569 } else if (p == theElectron) {
570 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
571
572 } else if (p == thePositron) {
573 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
574 e < (*theCutsPositron)[currentCoupleIndex]) {
575 good = false;
576 e += 2.0*electron_mass_c2;
577 }
578 }
579 // added secondary if it is good
580 }
581 if (good) {
582 G4Track* t = new G4Track(dp, time, track.GetPosition());
584 if (biasManager) {
585 t->SetWeight(weight * biasManager->GetWeight(i));
586 } else {
587 t->SetWeight(weight);
588 }
590
591 // define type of secondary
592 if(i < mainSecondaries) {
594 if(GetProcessSubType() == fComptonScattering && p == theGamma) {
596 }
597 } else if(i < mainSecondaries + n1) {
599 } else if(i < mainSecondaries + n1 + n2) {
601 } else {
602 if(i < num0) {
603 if(p == theGamma) {
605 } else {
607 }
608 } else {
610 }
611 }
612 /*
613 G4cout << "Secondary(post step) has weight " << t->GetWeight()
614 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
615 << GetProcessName() << " fluoID= " << fluoID
616 << " augerID= " << augerID <<G4endl;
617 */
618 } else {
619 delete dp;
620 edep += e;
621 }
622 }
623 }
625 }
626
629 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
632 }
633
634 return &fParticleChange;
635}
636
637//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
638
640 const G4String& directory,
641 G4bool ascii)
642{
643 if(!isTheMaster || part != particle) { return true; }
644 if(G4EmTableUtil::StoreTable(this, part, theLambdaTable,
645 directory, "Lambda",
646 verboseLevel, ascii) &&
647 G4EmTableUtil::StoreTable(this, part, theLambdaTablePrim,
648 directory, "LambdaPrim",
649 verboseLevel, ascii)) {
650 return true;
651 }
652 return false;
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
656
658 const G4String& dir,
659 G4bool ascii)
660{
661 if(!isTheMaster || part != particle) { return true; }
662 G4bool yes = true;
663 if(buildLambdaTable) {
664 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTable, dir,
665 "Lambda", verboseLevel,
666 ascii, splineFlag);
667 }
668 if(yes && minKinEnergyPrim < maxKinEnergy) {
669 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTablePrim, dir,
670 "LambdaPrim", verboseLevel,
671 ascii, splineFlag);
672 }
673 return yes;
674}
675
676//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677
679 const G4MaterialCutsCouple* couple)
680{
681 CurrentSetup(couple, kinEnergy);
682 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
683}
684
685//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
686
688 G4double,
690{
692 return G4VEmProcess::MeanFreePath(track);
693}
694
695//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
696
700{
702 return (currentModel) ?
703 currentModel->ComputeCrossSectionPerAtom(currentParticle, kinEnergy,
704 Z, A, cut) : 0.0;
705}
706
707//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
708
711{
712 DefineMaterial(couple);
713 G4PhysicsVector* newv = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy,
714 nLambdaBins, splineFlag);
715 return newv;
716}
717
718//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
719
721{
722 return (nullptr != currentModel) ?
723 currentModel->GetCurrentElement(currentMaterial) : nullptr;
724}
725
726//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
727
729{
730 return (nullptr != currentModel) ?
731 currentModel->GetCurrentElement(currentMaterial) : nullptr;
732}
733
734//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
735
737{
738 return (nullptr != currentModel) ?
739 currentModel->GetCurrentIsotope(GetCurrentElement()) : nullptr;
740}
741
742//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
743
745{
746 if(f > 0.0) {
747 biasFactor = f;
748 weightFlag = flag;
749 if(1 < verboseLevel) {
750 G4cout << "### SetCrossSectionBiasingFactor: for "
751 << particle->GetParticleName()
752 << " and process " << GetProcessName()
753 << " biasFactor= " << f << " weightFlag= " << flag
754 << G4endl;
755 }
756 }
757}
758
759//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
760
761void
763 G4bool flag)
764{
765 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
766 if(1 < verboseLevel) {
767 G4cout << "### ActivateForcedInteraction: for "
768 << particle->GetParticleName()
769 << " and process " << GetProcessName()
770 << " length(mm)= " << length/mm
771 << " in G4Region <" << r
772 << "> weightFlag= " << flag
773 << G4endl;
774 }
775 weightFlag = flag;
777}
778
779//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780
781void
783 G4double factor,
784 G4double energyLimit)
785{
786 if (0.0 <= factor) {
787
788 // Range cut can be applied only for e-
789 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
790 { return; }
791
793 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
794 if(1 < verboseLevel) {
795 G4cout << "### ActivateSecondaryBiasing: for "
796 << " process " << GetProcessName()
797 << " factor= " << factor
798 << " in G4Region <" << region
799 << "> energyLimit(MeV)= " << energyLimit/MeV
800 << G4endl;
801 }
802 }
803}
804
805//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
806
808{
809 if(5 < n && n < 10000000) {
810 nLambdaBins = n;
811 actBinning = true;
812 } else {
813 G4double e = (G4double)n;
814 PrintWarning("SetLambdaBinning", e);
815 }
816}
817
818//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
819
821{
822 if(1.e-3*eV < e && e < maxKinEnergy) {
823 nLambdaBins = G4lrint(nLambdaBins*G4Log(maxKinEnergy/e)
824 /G4Log(maxKinEnergy/minKinEnergy));
825 minKinEnergy = e;
826 actMinKinEnergy = true;
827 } else { PrintWarning("SetMinKinEnergy", e); }
828}
829
830//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
831
833{
834 if(minKinEnergy < e && e < 1.e+6*TeV) {
835 nLambdaBins = G4lrint(nLambdaBins*G4Log(e/minKinEnergy)
836 /G4Log(maxKinEnergy/minKinEnergy));
837 maxKinEnergy = e;
838 actMaxKinEnergy = true;
839 } else { PrintWarning("SetMaxKinEnergy", e); }
840}
841
842//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
843
845{
846 if(theParameters->MinKinEnergy() <= e &&
847 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
848 else { PrintWarning("SetMinKinEnergyPrim", e); }
849}
850
851//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
852
854{
855 return (nam == GetProcessName()) ? this : nullptr;
856}
857
858//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
859
861{
862 return theParameters->MscThetaLimit();
863}
864
865//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
866
867void G4VEmProcess::PrintWarning(G4String tit, G4double val)
868{
869 G4String ss = "G4VEmProcess::" + tit;
871 ed << "Parameter is out of range: " << val
872 << " it will have no effect!\n" << " Process "
873 << GetProcessName() << " nbins= " << theParameters->NumberOfBins()
874 << " Emin(keV)= " << theParameters->MinKinEnergy()/keV
875 << " Emax(GeV)= " << theParameters->MaxKinEnergy()/GeV;
876 G4Exception(ss, "em0044", JustWarning, ed);
877}
878
879//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
880
881void G4VEmProcess::ProcessDescription(std::ostream& out) const
882{
883 if(nullptr != particle) {
884 StreamInfo(out, *particle, true);
885 }
886}
887
888//....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)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
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 G4int Z[17]
const G4double A[17]
#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 *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: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)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4VAtomDeexcitation * AtomDeexcitation()
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:93
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
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
void SetCreatorModelID(const G4int id)
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:308
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:284
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:390
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
Definition: G4VEmModel.cc:261
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:774
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
Definition: G4VEmModel.cc:242
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
Definition: G4VEmProcess.hh:94
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:79
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
G4EmBiasingManager * biasManager
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
G4double preStepLogKinEnergy
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
G4bool isTheMaster
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
Definition: G4VProcess.hh:339
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:416
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:522
G4int verboseLevel
Definition: G4VProcess.hh:360
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62