Geant4-11
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"
63#include "G4Region.hh"
64#include "G4Gamma.hh"
65#include "G4Electron.hh"
66#include "G4Positron.hh"
68#include "G4EmBiasingManager.hh"
69#include "G4EmParameters.hh"
70#include "G4EmProcessSubType.hh"
72#include "G4DNAModelSubType.hh"
73#include "G4GenericIon.hh"
74#include "G4Log.hh"
75#include <iostream>
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
81{
84
85 // Size of tables
88
89 // default lambda factor
91
92 // particle types
96
99 secParticles.reserve(5);
100
103 lManager->Register(this);
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110
112{
113 /*
114 if(1 < verboseLevel) {
115 G4cout << "G4VEmProcess destruct " << GetProcessName()
116 << " " << this << " " << theLambdaTable <<G4endl;
117 }
118 */
119 if(isTheMaster) {
120 delete theData;
122 }
123 delete modelManager;
124 delete biasManager;
125 lManager->DeRegister(this);
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
129
131{
132 currentCouple = nullptr;
133 preStepLambda = 0.0;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
139 const G4Material*)
140{
141 return 0.0;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145
147 const G4Region* region)
148{
149 if(nullptr == ptr) { return; }
150 G4VEmFluctuationModel* fm = nullptr;
151 modelManager->AddEmModel(order, ptr, fm, region);
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
158{
159 if(nullptr == ptr) { return; }
160 if(!emModels.empty()) {
161 for(auto & em : emModels) { if(em == ptr) { return; } }
162 }
163 emModels.push_back(ptr);
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167
169{
170 return modelManager->GetModel(idx, ver);
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
176{
178 if(nullptr == particle) { SetParticle(&part); }
179
180 if(part.GetParticleType() == "nucleus" &&
181 part.GetParticleSubType() == "generic") {
182
184 if(pname != "deuteron" && pname != "triton" &&
185 pname != "alpha" && pname != "He3" &&
186 pname != "alpha+" && pname != "helium" &&
187 pname != "hydrogen") {
188
190 isIon = true;
191 }
192 }
193
194 if(1 < verboseLevel) {
195 G4cout << "G4VEmProcess::PreparePhysicsTable() for "
196 << GetProcessName()
197 << " and particle " << part.GetParticleName()
198 << " local particle " << particle->GetParticleName()
199 << G4endl;
200 }
201
202 if(particle != &part) { return; }
203
205
206 Clear();
208
210 const G4ProductionCutsTable* theCoupleTable=
212
213 // initialisation of the process
216
217 if(isTheMaster) {
219 if(nullptr == theData) { theData = new G4EmDataHandler(2); }
220 if(fEmOnePeak == fXSType) {
221 if(nullptr == theEnergyOfCrossSectionMax) {
222 theEnergyOfCrossSectionMax = new std::vector<G4double>;
223 }
224 size_t n = theCoupleTable->GetTableSize();
226 }
227 } else {
229 }
234
235 // integral option may be disabled
237
238 // prepare tables
242 }
243 // high energy table
247 }
249
250 // initialisation of models
252 for(G4int i=0; i<numberOfModels; ++i) {
254 if(nullptr == mod) { continue; }
255 if(nullptr == currentModel) { currentModel = mod; }
258 if(mod->HighEnergyLimit() > maxKinEnergy) {
260 }
261 SetEmModel(mod);
263 }
264
265 if(nullptr != lManager->AtomDeexcitation()) {
267 }
268 fLambdaEnergy = 0.0;
269
270 theCuts =
275
276 // forced biasing
277 if(biasManager) {
279 biasFlag = false;
280 }
281
282 // defined ID of secondary particles
283 G4int stype = GetProcessSubType();
284 if(stype == fAnnihilation) {
287 } else if(stype == fGammaConversion) {
289 mainSecondaries = 2;
290 } else if(stype == fPhotoElectricEffect) {
292 } else if(stype == fComptonScattering) {
294 } else if(stype >= fLowEnergyElastic) {
296 }
297 if(1 < verboseLevel) {
298 G4cout << "### G4VEmProcess::PreparePhysicsTable() done for "
299 << GetProcessName()
300 << " and particle " << part.GetParticleName()
301 << " baseMat=" << baseMat << G4endl;
302 }
303}
304
305//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
306
308{
309 if(nullptr == masterProc) {
310 if(isTheMaster) { masterProc = this; }
311 else { masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());}
312 }
313
314 G4String num = part.GetParticleName();
315 if(1 < verboseLevel) {
316 G4cout << "### G4VEmProcess::BuildPhysicsTable() for "
317 << GetProcessName()
318 << " and particle " << num
319 << " buildLambdaTable= " << buildLambdaTable
320 << " isTheMaster= " << isTheMaster
321 << " " << masterProc
322 << G4endl;
323 }
324
325 if(particle == &part) {
326
327 // worker initialisation
328 if(!isTheMaster) {
331 if(fXSType == fEmOnePeak) {
333 }
335
336 // local initialisation of models
337 G4bool printing = true;
338 for(G4int i=0; i<numberOfModels; ++i) {
339 G4VEmModel* mod = GetModelByIndex(i, printing);
340 G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing);
341 //G4cout << i << ". " << mod << " " << mod0 << " "
342 // << particle->GetParticleName() << G4endl;
344 mod->InitialiseLocal(particle, mod0);
345 }
346 // master thread
347 } else {
350 }
351 if(fXSType == fEmOnePeak) {
355 }
356 }
357 }
358 // protection against double printout
359 if(theParameters->IsPrintLocked()) { return; }
360
361 // explicitly defined printout by particle name
362 if(1 < verboseLevel ||
363 (0 < verboseLevel && (num == "gamma" || num == "e-" ||
364 num == "e+" || num == "mu+" ||
365 num == "mu-" || num == "proton"||
366 num == "pi+" || num == "pi-" ||
367 num == "kaon+" || num == "kaon-" ||
368 num == "alpha" || num == "anti_proton" ||
369 num == "GenericIon"|| num == "alpha++" ||
370 num == "alpha+" || num == "helium" ||
371 num == "hydrogen")))
372 {
373 StreamInfo(G4cout, part);
374 }
375
376 if(1 < verboseLevel) {
377 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
378 << GetProcessName()
379 << " and particle " << num
380 << " baseMat=" << baseMat
381 << G4endl;
382 }
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386
388{
389 if(1 < verboseLevel) {
390 G4cout << "G4EmProcess::BuildLambdaTable() for process "
391 << GetProcessName() << " and particle "
392 << particle->GetParticleName() << " " << this
393 << G4endl;
394 }
395
396 // Access to materials
397 const G4ProductionCutsTable* theCoupleTable=
399 size_t numOfCouples = theCoupleTable->GetTableSize();
400
402
403 G4PhysicsLogVector* aVector = nullptr;
404 G4PhysicsLogVector* aVectorPrim = nullptr;
405 G4PhysicsLogVector* bVectorPrim = nullptr;
406
408 G4int nbin =
409 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
410 scale = G4Log(scale);
411 if(actBinning) { nbin = std::max(nbin, nLambdaBins); }
413
414 for(size_t i=0; i<numOfCouples; ++i) {
415
416 if (bld->GetFlag(i)) {
417
418 // create physics vector and fill it
419 const G4MaterialCutsCouple* couple =
420 theCoupleTable->GetMaterialCutsCouple(i);
421
422 // build main table
423 if(buildLambdaTable) {
424 delete (*theLambdaTable)[i];
425
426 // if start from zero then change the scale
427 G4double emin = minKinEnergy;
428 G4bool startNull = false;
429 if(startFromNull) {
431 if(e >= emin) {
432 emin = e;
433 startNull = true;
434 }
435 }
436 G4double emax = emax1;
437 if(emax <= emin) { emax = 2*emin; }
438 G4int bin = G4lrint(nbin*G4Log(emax/emin)/scale);
439 if(bin < 3) { bin = 3; }
440 aVector = new G4PhysicsLogVector(emin, emax, bin, splineFlag);
441 modelManager->FillLambdaVector(aVector, couple, startNull);
442 if(splineFlag) { aVector->FillSecondDerivatives(); }
444 }
445 // build high energy table
447 delete (*theLambdaTablePrim)[i];
448
449 // start not from zero and always use spline
450 if(!bVectorPrim) {
452 if(bin < 3) { bin = 3; }
453 aVectorPrim =
455 bVectorPrim = aVectorPrim;
456 } else {
457 aVectorPrim = new G4PhysicsLogVector(*bVectorPrim);
458 }
459 modelManager->FillLambdaVector(aVectorPrim, couple, false,
461 aVectorPrim->FillSecondDerivatives();
463 aVectorPrim);
464 }
465 }
466 }
467
468 if(1 < verboseLevel) {
469 G4cout << "Lambda table is built for "
471 << G4endl;
472 }
473}
474
475//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
476
477void G4VEmProcess::StreamInfo(std::ostream& out,
478 const G4ParticleDefinition& part, G4bool rst) const
479{
480 G4String indent = (rst ? " " : "");
481 out << std::setprecision(6);
482 out << G4endl << indent << GetProcessName() << ": ";
483 if (!rst) {
484 out << " for " << part.GetParticleName();
485 }
486 if(fXSType != fEmNoIntegral) { out << " XStype:" << fXSType; }
487 if(applyCuts) { out << " applyCuts:1 "; }
488 out << " SubType=" << GetProcessSubType();
489 if(biasFactor != 1.0) { out << " BiasingFactor= " << biasFactor; }
490 out << " BuildTable=" << buildLambdaTable << G4endl;
491 if(buildLambdaTable) {
492 if(particle == &part) {
493 size_t length = theLambdaTable->length();
494 for(size_t i=0; i<length; ++i) {
495 G4PhysicsVector* v = (*theLambdaTable)[i];
496 if(v) {
497 out << " Lambda table from ";
498 G4double emin = v->Energy(0);
500 G4int nbin = v->GetVectorLength() - 1;
501 if(emin > minKinEnergy) { out << "threshold "; }
502 else { out << G4BestUnit(emin,"Energy"); }
503 out << " to "
504 << G4BestUnit(emax,"Energy")
505 << ", " << G4lrint(nbin/std::log10(emax/emin))
506 << " bins/decade, spline: "
507 << splineFlag << G4endl;
508 break;
509 }
510 }
511 } else {
512 out << " Used Lambda table of "
514 }
515 }
517 if(particle == &part) {
518 size_t length = theLambdaTablePrim->length();
519 for(size_t i=0; i<length; ++i) {
520 G4PhysicsVector* v = (*theLambdaTablePrim)[i];
521 if(v) {
522 out << " LambdaPrime table from "
523 << G4BestUnit(v->Energy(0),"Energy")
524 << " to "
525 << G4BestUnit(v->GetMaxEnergy(),"Energy")
526 << " in " << v->GetVectorLength()-1
527 << " bins " << G4endl;
528 break;
529 }
530 }
531 } else {
532 out << " Used LambdaPrime table of "
534 }
535 }
538
539 if(verboseLevel > 2 && buildLambdaTable) {
540 out << " LambdaTable address= " << theLambdaTable << G4endl;
541 if(theLambdaTable && particle == &part) {
542 out << (*theLambdaTable) << G4endl;
543 }
544 }
545}
546
547//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
548
550{
551 // reset parameters for the new track
555
557
558 // forced biasing only for primary particles
559 if(biasManager) {
560 if(0 == track->GetParentID()) {
561 // primary particle
562 biasFlag = true;
564 }
565 }
566}
567
568//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
569
571 const G4Track& track,
572 G4double previousStepSize,
574{
576 G4double x = DBL_MAX;
577
581 const G4double scaledEnergy = preStepKinEnergy*massRatio;
582 SelectModel(scaledEnergy, currentCoupleIndex);
583 /*
584 G4cout << "PostStepGetPhysicalInteractionLength: idx= " << currentCoupleIndex
585 << " couple: " << currentCouple << G4endl;
586 */
587 if(!currentModel->IsActive(scaledEnergy)) {
590 return x;
591 }
592
593 // forced biasing only for primary particles
594 if(biasManager) {
595 if(0 == track.GetParentID()) {
596 if(biasFlag &&
598 return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
599 }
600 }
601 }
602
603 // compute mean free path
605
606 // zero cross section
607 if(preStepLambda <= 0.0) {
610
611 } else {
612
613 // non-zero cross section
615
616 // beggining of tracking (or just after DoIt of this process)
619
620 } else if(currentInteractionLength < DBL_MAX) {
621
623 previousStepSize/currentInteractionLength;
626 }
627
628 // new mean free path and step limit for the next step
631 }
632 return x;
633}
634
635//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
636
638{
639 if(fXSType == fEmNoIntegral) {
641
642 } else if(fXSType == fEmIncreasing) {
643 if(e/lambdaFactor < mfpKinEnergy) {
644 mfpKinEnergy = e;
646 }
647
648 } else if(fXSType == fEmDecreasing) {
649 if(e < mfpKinEnergy) {
650 const G4double e1 = e*lambdaFactor;
653 }
654
655 } else if(fXSType == fEmOnePeak) {
656 const G4double epeak = (*theEnergyOfCrossSectionMax)[currentCoupleIndex];
657 if(e <= epeak) {
658 if(e/lambdaFactor < mfpKinEnergy) {
659 mfpKinEnergy = e;
661 }
662 } else if(e < mfpKinEnergy) {
663 const G4double e1 = std::max(epeak, e*lambdaFactor);
666 }
667
668 } else {
670 }
671}
672
673//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
674
676 const G4Step& step)
677{
678 // In all cases clear number of interaction lengths
681
683
684 // Do not make anything if particle is stopped, the annihilation then
685 // should be performed by the AtRestDoIt!
686 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
687
688 const G4double finalT = track.GetKineticEnergy();
689
690 // forced process - should happen only once per track
691 if(biasFlag) {
693 biasFlag = false;
694 }
695 }
696
697 // check active and select model
698 const G4double scaledEnergy = finalT*massRatio;
699 SelectModel(scaledEnergy, currentCoupleIndex);
700 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
701
702 // Integral approach
703 if (fXSType != fEmNoIntegral) {
704 const G4double logFinalT = track.GetDynamicParticle()->GetLogKineticEnergy();
705 const G4double lx = std::max(GetCurrentLambda(finalT, logFinalT), 0.0);
706 const G4double lg = preStepLambda;
707 if(finalT < mfpKinEnergy) {
708 mfpKinEnergy = finalT;
709 preStepLambda = lx;
710 }
711#ifdef G4VERBOSE
712 if(lg < lx && 1 < verboseLevel) {
713 G4cout << "WARNING: for " << currentParticle->GetParticleName()
714 << " and " << GetProcessName()
715 << " E(MeV)= " << finalT/MeV
716 << " preLambda= " << lg << " < " << lx << " (postLambda) "
717 << G4endl;
718 }
719#endif
720 if(lg*G4UniformRand() >= lx) {
721 return &fParticleChange;
722 }
723 }
724
725 // define new weight for primary and secondaries
727 if(weightFlag) {
728 weight /= biasFactor;
730 }
731
732#ifdef G4VERBOSE
733 if(1 < verboseLevel) {
734 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
735 << finalT/MeV
736 << " MeV; model= (" << currentModel->LowEnergyLimit()
737 << ", " << currentModel->HighEnergyLimit() << ")"
738 << G4endl;
739 }
740#endif
741
742 // sample secondaries
743 secParticles.clear();
746 track.GetDynamicParticle(),
748
749 G4int num0 = secParticles.size();
750
751 // splitting or Russian roulette
752 if(biasManager) {
754 G4double eloss = 0.0;
756 secParticles, track, currentModel, &fParticleChange, eloss,
758 step.GetPostStepPoint()->GetSafety());
759 if(eloss > 0.0) {
762 }
763 }
764 }
765
766 // save secondaries
767 G4int num = secParticles.size();
768 if(num > 0) {
769
772 G4double time = track.GetGlobalTime();
773
774 G4int n1(0), n2(0);
775 if(num > mainSecondaries) {
777 }
778
779 for (G4int i=0; i<num; ++i) {
781 if (nullptr != dp) {
783 G4double e = dp->GetKineticEnergy();
784 G4bool good = true;
785 if(applyCuts) {
786 if (p == theGamma) {
787 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
788
789 } else if (p == theElectron) {
790 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
791
792 } else if (p == thePositron) {
795 good = false;
796 e += 2.0*electron_mass_c2;
797 }
798 }
799 // added secondary if it is good
800 }
801 if (good) {
802 G4Track* t = new G4Track(dp, time, track.GetPosition());
804 if (biasManager) {
805 t->SetWeight(weight * biasManager->GetWeight(i));
806 } else {
807 t->SetWeight(weight);
808 }
810
811 // define type of secondary
812 if(i < mainSecondaries) {
816 }
817 } else if(i < mainSecondaries + n1) {
819 } else if(i < mainSecondaries + n1 + n2) {
821 } else {
822 if(i < num0) {
823 if(p == theGamma) {
825 } else {
827 }
828 } else {
830 }
831 }
832 /*
833 G4cout << "Secondary(post step) has weight " << t->GetWeight()
834 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
835 << GetProcessName() << " fluoID= " << fluoID
836 << " augerID= " << augerID <<G4endl;
837 */
838 } else {
839 delete dp;
840 edep += e;
841 }
842 }
843 }
845 }
846
852 }
853
854 return &fParticleChange;
855}
856
857//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
858
860 const G4String& directory,
861 G4bool ascii)
862{
863 G4bool yes = true;
864 if(!isTheMaster) { return yes; }
865
866 if ( theLambdaTable && part == particle) {
867 const G4String& nam =
868 GetPhysicsTableFileName(part,directory,"Lambda",ascii);
869 yes = theLambdaTable->StorePhysicsTable(nam,ascii);
870
871 if ( yes ) {
872 if(0 < verboseLevel) G4cout << "Stored: " << nam << G4endl;
873 } else {
874 G4cout << "Fail to store Physics Table for "
876 << " and process " << GetProcessName()
877 << " in the directory <" << directory
878 << "> " << G4endl;
879 }
880 }
881 if ( theLambdaTablePrim && part == particle) {
882 const G4String& name =
883 GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
885
886 if ( yes ) {
887 if(0 < verboseLevel) {
888 G4cout << "Physics table prim is stored for "
890 << " and process " << GetProcessName()
891 << " in the directory <" << directory
892 << "> " << G4endl;
893 }
894 } else {
895 G4cout << "Fail to store Physics Table Prim for "
897 << " and process " << GetProcessName()
898 << " in the directory <" << directory
899 << "> " << G4endl;
900 }
901 }
902 return yes;
903}
904
905//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
906
908 const G4String& directory,
909 G4bool ascii)
910{
911 if(1 < verboseLevel) {
912 G4cout << "G4VEmProcess::RetrievePhysicsTable() for "
913 << part->GetParticleName() << " and process "
914 << GetProcessName() << G4endl;
915 }
916 G4bool yes = true;
917
919 || particle != part) { return yes; }
920
921 const G4String particleName = part->GetParticleName();
922
923 if(buildLambdaTable) {
924 const G4String& filename =
925 GetPhysicsTableFileName(part,directory,"Lambda",ascii);
927 filename,ascii,
928 splineFlag);
929 if ( yes ) {
930 if (0 < verboseLevel) {
931 G4cout << "Lambda table for " << particleName
932 << " is Retrieved from <"
933 << filename << ">"
934 << G4endl;
935 }
936 if(splineFlag) {
937 for(auto & v : *theLambdaTable) {
938 if(nullptr != v) { v->FillSecondDerivatives(); }
939 }
940 }
941
942 } else {
943 if (1 < verboseLevel) {
944 G4cout << "Lambda table for " << particleName << " in file <"
945 << filename << "> is not exist"
946 << G4endl;
947 }
948 }
949 }
951 const G4String& filename =
952 GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
954 filename,ascii,true);
955 if ( yes ) {
956 if (0 < verboseLevel) {
957 G4cout << "Lambda table prim for " << particleName
958 << " is Retrieved from <"
959 << filename << ">"
960 << G4endl;
961 }
962 for(auto & v : *theLambdaTablePrim) {
963 if(nullptr != v) { v->FillSecondDerivatives(); }
964 }
965 } else {
966 if (1 < verboseLevel) {
967 G4cout << "Lambda table prim for " << particleName << " in file <"
968 << filename << "> is not exist"
969 << G4endl;
970 }
971 }
972 }
973 return yes;
974}
975
976//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
977
979 const G4MaterialCutsCouple* couple,
980 G4double)
981{
982 G4double cross = RecalculateLambda(kinEnergy, couple);
983 return std::max(cross, 0.0);
984}
985
986//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
987
989 G4double,
991{
993 return G4VEmProcess::MeanFreePath(track);
994}
995
996//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
997
999{
1000 const G4double kinEnergy = track.GetKineticEnergy();
1001 CurrentSetup(track.GetMaterialCutsCouple(), kinEnergy);
1002 const G4double xs = GetCurrentLambda(kinEnergy,
1004 return (0.0 < xs) ? 1.0/xs : DBL_MAX;
1005}
1006
1007//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1008
1009G4double
1012{
1013 SelectModel(kinEnergy, currentCoupleIndex);
1014 return (currentModel) ?
1016 Z, A, cut) : 0.0;
1017}
1018
1019//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1020
1021std::vector<G4double>* G4VEmProcess::FindLambdaMax()
1022{
1023 if(1 < verboseLevel) {
1024 G4cout << "### G4VEmProcess::FindLambdaMax: "
1026 << " and process " << GetProcessName() << " " << G4endl;
1027 }
1028 std::vector<G4double>* ptr = nullptr;
1029 if(fXSType != fEmOnePeak) { return ptr; }
1030
1031 const G4ProductionCutsTable* theCoupleTable=
1033 size_t n = theCoupleTable->GetTableSize();
1034 ptr = new std::vector<G4double>;
1035 ptr->resize(n, DBL_MAX);
1036
1037 G4bool isPeak = false;
1038 const G4double g4log10 = G4Log(10.);
1039 const G4double scale = theParameters->NumberOfBinsPerDecade()/g4log10;
1040
1041 for(size_t i=0; i<n; ++i) {
1042 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
1044 G4double emax = std::max(maxKinEnergy, emin + emin);
1045 G4double ee = G4Log(emax/emin);
1046 G4int nbin = G4lrint(ee*scale);
1047 if(nbin < 4) { nbin = 4; }
1048 G4double x = G4Exp(ee/nbin);
1049 G4double sm = 0.0;
1050 G4double em = emin;
1051 G4double e = emin;
1052 for(G4int j=0; j<=nbin; ++j) {
1053 G4double sig = RecalculateLambda(e, couple);
1054 //G4cout << j << " E=" << e << " Lambda=" << sig << G4endl;
1055 if(sig >= sm) {
1056 em = e;
1057 sm = sig;
1058 e *= x;
1059 } else {
1060 isPeak = true;
1061 (*ptr)[i] = em;
1062 break;
1063 }
1064 }
1065 if(1 < verboseLevel) {
1066 G4cout << " " << i << ". Epeak(GeV)=" << em/GeV
1067 << " SigmaMax(1/mm)=" << sm
1068 << " Emin(GeV)=" << emin/GeV << " Emax(GeV)=" << emax/GeV
1069 << " " << couple->GetMaterial()->GetName() << G4endl;
1070 }
1071 }
1072 // there is no peak for any material
1073 if(!isPeak) {
1074 delete ptr;
1075 ptr = nullptr;
1076 }
1077 return ptr;
1078}
1079
1080//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1081
1082void G4VEmProcess::SetEnergyOfCrossSectionMax(std::vector<G4double>* ptr)
1083{
1084 if(nullptr == ptr) {
1086 } else {
1088 }
1089}
1090
1091//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1092
1095{
1096 DefineMaterial(couple);
1099 return newv;
1100}
1101
1102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1103
1105{
1106 return (nullptr != currentModel) ? currentModel->GetCurrentElement() : nullptr;
1107}
1108
1109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1110
1112{
1113 if(f > 0.0) {
1114 biasFactor = f;
1115 weightFlag = flag;
1116 if(1 < verboseLevel) {
1117 G4cout << "### SetCrossSectionBiasingFactor: for "
1119 << " and process " << GetProcessName()
1120 << " biasFactor= " << f << " weightFlag= " << flag
1121 << G4endl;
1122 }
1123 }
1124}
1125
1126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1127
1128void
1130 G4bool flag)
1131{
1132 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
1133 if(1 < verboseLevel) {
1134 G4cout << "### ActivateForcedInteraction: for "
1136 << " and process " << GetProcessName()
1137 << " length(mm)= " << length/mm
1138 << " in G4Region <" << r
1139 << "> weightFlag= " << flag
1140 << G4endl;
1141 }
1142 weightFlag = flag;
1144}
1145
1146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1147
1148void
1150 G4double factor,
1151 G4double energyLimit)
1152{
1153 if (0.0 <= factor) {
1154
1155 // Range cut can be applied only for e-
1156 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1157 { return; }
1158
1160 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1161 if(1 < verboseLevel) {
1162 G4cout << "### ActivateSecondaryBiasing: for "
1163 << " process " << GetProcessName()
1164 << " factor= " << factor
1165 << " in G4Region <" << region
1166 << "> energyLimit(MeV)= " << energyLimit/MeV
1167 << G4endl;
1168 }
1169 }
1170}
1171
1172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1173
1175{
1176 if(5 < n && n < 10000000) {
1177 nLambdaBins = n;
1178 actBinning = true;
1179 } else {
1180 G4double e = (G4double)n;
1181 PrintWarning("SetLambdaBinning", e);
1182 }
1183}
1184
1185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1186
1188{
1189 if(1.e-3*eV < e && e < maxKinEnergy) {
1192 minKinEnergy = e;
1193 actMinKinEnergy = true;
1194 } else { PrintWarning("SetMinKinEnergy", e); }
1195}
1196
1197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1198
1200{
1201 if(minKinEnergy < e && e < 1.e+6*TeV) {
1204 maxKinEnergy = e;
1205 actMaxKinEnergy = true;
1206 } else { PrintWarning("SetMaxKinEnergy", e); }
1207}
1208
1209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1210
1212{
1213 if(theParameters->MinKinEnergy() <= e &&
1214 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
1215 else { PrintWarning("SetMinKinEnergyPrim", e); }
1216}
1217
1218//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1219
1221{
1222 return (nam == GetProcessName()) ? this : nullptr;
1223}
1224
1225//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1226
1227G4double
1229{
1230 CurrentSetup(couple, kinEnergy);
1231 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
1232}
1233
1234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1235
1237{
1238 return theParameters->MscThetaLimit();
1239}
1240
1241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1242
1244{
1245 G4String ss = "G4VEmProcess::" + tit;
1247 ed << "Parameter is out of range: " << val
1248 << " it will have no effect!\n" << " Process "
1249 << GetProcessName() << " nbins= " << theParameters->NumberOfBins()
1250 << " Emin(keV)= " << theParameters->MinKinEnergy()/keV
1251 << " Emax(GeV)= " << theParameters->MaxKinEnergy()/GeV;
1252 G4Exception(ss, "em0044", JustWarning, ed);
1253}
1254
1255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1256
1257void G4VEmProcess::ProcessDescription(std::ostream& out) const
1258{
1259 if(particle) {
1260 StreamInfo(out, *particle, true);
1261 }
1262}
1263
1264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
static const G4double e1[44]
static const G4double emax
@ fDNAUnknownModel
@ fGammaConversion
@ fComptonScattering
@ fAnnihilation
@ fPhotoElectricEffect
@ fIsCrossSectionPrim
@ 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:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ProcessType
@ idxG4ElectronCut
@ idxG4GammaCut
@ idxG4PositronCut
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double TeV
Definition: G4SIunits.hh:204
#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
void SetFluoFlag(G4bool val)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double, G4int verb)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
G4bool IsPrintLocked() const
static G4EmParameters * Instance()
G4int NumberOfBins() const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MscThetaLimit() const
void DefineRegParamForEM(G4VEmProcess *) const
G4int Verbose() const
G4int WorkerVerbose() const
G4bool ApplyCuts() const
G4double MaxKinEnergy() const
G4bool Integral() const
G4double LambdaFactor() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
G4bool GetFlag(size_t idx)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:173
void InitializeForPostStep(const G4Track &)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii, G4bool spline)
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
std::size_t length() const
G4double GetMaxEnergy() const
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetSafety() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
void SetCreatorModelID(const G4int id)
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:365
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:802
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:739
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:505
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:795
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:753
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:204
const std::vector< G4double > * theCutsElectron
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4double RecalculateLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4double fLambdaEnergy
G4double maxKinEnergy
void ComputeIntegralLambda(G4double kinEnergy, G4double logKinEnergy)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double MeanFreePath(const G4Track &track)
G4bool applyCuts
G4int mainSecondaries
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
virtual void StreamProcessInfo(std::ostream &) const
Definition: G4VEmProcess.hh:94
G4double massRatio
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:79
G4LossTableManager * lManager
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
G4EmBiasingManager * biasManager
const G4ParticleDefinition * secondaryParticle
G4bool weightFlag
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
G4double preStepLogKinEnergy
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
const std::vector< G4double > * theCutsGamma
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4VEmModel * currentModel
void SetEnergyOfCrossSectionMax(std::vector< G4double > *)
G4bool splineFlag
G4PhysicsTable * LambdaTable() const
std::vector< G4double > * theEnergyOfCrossSectionMax
G4bool actMaxKinEnergy
void SetMinKinEnergy(G4double e)
G4int numberOfModels
const G4ParticleDefinition * currentParticle
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
std::vector< G4VEmModel * > emModels
void StartTracking(G4Track *) override
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4bool actMinKinEnergy
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)
const std::vector< G4double > * theDensityFactor
void PrintWarning(G4String tit, G4double val)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
std::vector< G4double > * EnergyOfCrossSectionMax() const
G4bool UseBaseMaterial() const
G4CrossSectionType fXSType
void ProcessDescription(std::ostream &outFile) const override
G4double GetCurrentLambda(G4double kinEnergy)
G4PhysicsTable * theLambdaTablePrim
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4bool actBinning
G4double lambdaFactor
virtual G4VEmProcess * GetEmProcess(const G4String &name)
G4double MaxKinEnergy() const
G4double logLambdaFactor
G4bool startFromNull
const G4ParticleDefinition * thePositron
G4bool isTheMaster
G4double biasFactor
std::vector< G4DynamicParticle * > secParticles
std::vector< G4double > * FindLambdaMax()
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
const G4ParticleDefinition * theElectron
G4EmParameters * theParameters
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * theGamma
G4double minKinEnergyPrim
G4bool buildLambdaTable
G4PhysicsTable * LambdaTablePrim() const
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
const G4VEmProcess * masterProc
const std::vector< G4double > * theCuts
size_t currentCoupleIndex
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4ParticleChangeForGamma fParticleChange
G4PhysicsTable * theLambdaTable
const std::vector< G4int > * theDensityIdx
const std::vector< G4double > * theCutsPositron
G4double minKinEnergy
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * particle
void SetMinKinEnergyPrim(G4double e)
void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMaxKinEnergy(G4double e)
G4EmModelManager * modelManager
G4EmDataHandler * theData
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:335
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
static constexpr double TeV
static constexpr double keV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const char * name(G4int ptype)
string pname
Definition: eplot.py:33
float electron_mass_c2
Definition: hepunit.py:273
float proton_mass_c2
Definition: hepunit.py:274
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62