1 #include "G4Version.hh" 2 #if G4VERSION_NUMBER >= 1100 7 #include "G4CoulombScattering.hh" 8 #include "G4UrbanMscModel.hh" 9 #include "G4WentzelVIModel.hh" 10 #include "G4eBremsstrahlung.hh" 11 #include "G4eCoulombScatteringModel.hh" 12 #include "G4eIonisation.hh" 13 #include "G4eMultipleScattering.hh" 14 #include "G4eplusAnnihilation.hh" 16 #include "G4ElectroVDNuclearModel.hh" 17 #include "G4ElectronNuclearProcess.hh" 18 #include "G4PositronNuclearProcess.hh" 20 #include "G4ComptonScattering.hh" 21 #include "G4GammaConversion.hh" 22 #include "G4PhotoElectricEffect.hh" 24 #include "G4CascadeInterface.hh" 25 #include "G4CrossSectionDataSetRegistry.hh" 26 #include "G4ExcitedStringDecay.hh" 27 #include "G4GammaNuclearXS.hh" 28 #include "G4GammaParticipants.hh" 29 #include "G4GeneratorPrecompoundInterface.hh" 30 #include "G4HadronInelasticProcess.hh" 31 #include "G4HadronicParameters.hh" 32 #include "G4LowEGammaNuclearModel.hh" 33 #include "G4PhotoNuclearCrossSection.hh" 34 #include "G4QGSMFragmentation.hh" 35 #include "G4QGSModel.hh" 36 #include "G4TheoFSGenerator.hh" 38 #include "G4GammaGeneralProcess.hh" 39 #include "G4LossTableManager.hh" 41 #include "G4EmParameters.hh" 42 #include "G4SystemOfUnits.hh" 44 #include "G4Electron.hh" 46 #include "G4Positron.hh" 48 #include "G4RegionStore.hh" 49 #include "G4Region.hh" 52 CMSEmStandardPhysicsTrackingManager *CMSEmStandardPhysicsTrackingManager::masterTrackingManager =
nullptr;
54 CMSEmStandardPhysicsTrackingManager::CMSEmStandardPhysicsTrackingManager(
const edm::ParameterSet &
p) {
55 fRangeFactor =
p.getParameter<
double>(
"G4MscRangeFactor");
56 fGeomFactor =
p.getParameter<
double>(
"G4MscGeomFactor");
57 fSafetyFactor =
p.getParameter<
double>(
"G4MscSafetyFactor");
58 fLambdaLimit =
p.getParameter<
double>(
"G4MscLambdaLimit") * CLHEP::mm;
60 fStepLimitType = fUseSafety;
61 if (msc ==
"UseSafetyPlus") {
62 fStepLimitType = fUseSafetyPlus;
63 }
else if (msc ==
"Minimal") {
64 fStepLimitType = fMinimal;
67 G4EmParameters *param = G4EmParameters::Instance();
68 G4double highEnergyLimit = param->MscEnergyLimit();
70 const G4Region *aRegion = G4RegionStore::GetInstance()->GetRegion(
"HcalRegion",
false);
71 const G4Region *bRegion = G4RegionStore::GetInstance()->GetRegion(
"HGCalRegion",
false);
75 G4eMultipleScattering *msc =
new G4eMultipleScattering;
76 G4UrbanMscModel *msc1 =
new G4UrbanMscModel;
77 G4WentzelVIModel *msc2 =
new G4WentzelVIModel;
78 msc1->SetHighEnergyLimit(highEnergyLimit);
79 msc2->SetLowEnergyLimit(highEnergyLimit);
80 msc->SetEmModel(msc1);
81 msc->SetEmModel(msc2);
84 if (
nullptr != aRegion ||
nullptr != bRegion) {
85 G4UrbanMscModel *msc3 =
new G4UrbanMscModel();
86 msc3->SetHighEnergyLimit(highEnergyLimit);
87 msc3->SetRangeFactor(fRangeFactor);
88 msc3->SetGeomFactor(fGeomFactor);
89 msc3->SetSafetyFactor(fSafetyFactor);
90 msc3->SetLambdaLimit(fLambdaLimit);
91 msc3->SetStepLimitType(fStepLimitType);
92 msc3->SetLocked(
true);
93 if (
nullptr != aRegion) {
94 msc->AddEmModel(-1, msc3, aRegion);
96 if (
nullptr != bRegion) {
97 msc->AddEmModel(-1, msc3, bRegion);
104 electron.brems =
new G4eBremsstrahlung;
106 G4CoulombScattering *
ss =
new G4CoulombScattering;
107 G4eCoulombScatteringModel *ssm =
new G4eCoulombScatteringModel;
108 ssm->SetLowEnergyLimit(highEnergyLimit);
109 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
111 ss->SetMinKinEnergy(highEnergyLimit);
117 G4eMultipleScattering *msc =
new G4eMultipleScattering;
118 G4UrbanMscModel *msc1 =
new G4UrbanMscModel;
119 G4WentzelVIModel *msc2 =
new G4WentzelVIModel;
120 msc1->SetHighEnergyLimit(highEnergyLimit);
121 msc2->SetLowEnergyLimit(highEnergyLimit);
122 msc->SetEmModel(msc1);
123 msc->SetEmModel(msc2);
126 if (
nullptr != aRegion ||
nullptr != bRegion) {
127 G4UrbanMscModel *msc3 =
new G4UrbanMscModel();
128 msc3->SetHighEnergyLimit(highEnergyLimit);
129 msc3->SetRangeFactor(fRangeFactor);
130 msc3->SetGeomFactor(fGeomFactor);
131 msc3->SetSafetyFactor(fSafetyFactor);
132 msc3->SetLambdaLimit(fLambdaLimit);
133 msc3->SetStepLimitType(fStepLimitType);
134 msc3->SetLocked(
true);
135 if (
nullptr != aRegion) {
136 msc->AddEmModel(-1, msc3, aRegion);
138 if (
nullptr != bRegion) {
139 msc->AddEmModel(-1, msc3, bRegion);
145 positron.ioni =
new G4eIonisation;
146 positron.brems =
new G4eBremsstrahlung;
147 positron.annihilation =
new G4eplusAnnihilation;
149 G4CoulombScattering *
ss =
new G4CoulombScattering;
150 G4eCoulombScatteringModel *ssm =
new G4eCoulombScatteringModel;
151 ssm->SetLowEnergyLimit(highEnergyLimit);
152 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
154 ss->SetMinKinEnergy(highEnergyLimit);
162 gammaProc->AddEmProcess(
new G4PhotoElectricEffect);
163 gammaProc->AddEmProcess(
new G4ComptonScattering);
164 gammaProc->AddEmProcess(
new G4GammaConversion);
166 G4HadronInelasticProcess *nuc =
new G4HadronInelasticProcess(
"photonNuclear", G4Gamma::Definition());
167 auto xsreg = G4CrossSectionDataSetRegistry::Instance();
168 G4VCrossSectionDataSet *xs =
nullptr;
169 bool useGammaNuclearXS =
true;
170 if (useGammaNuclearXS) {
171 xs = xsreg->GetCrossSectionDataSet(
"GammaNuclearXS");
173 xs =
new G4GammaNuclearXS;
175 xs = xsreg->GetCrossSectionDataSet(
"PhotoNuclearXS");
177 xs =
new G4PhotoNuclearCrossSection;
181 G4QGSModel<G4GammaParticipants> *theStringModel =
new G4QGSModel<G4GammaParticipants>;
182 G4QGSMFragmentation *theFrag =
new G4QGSMFragmentation;
183 G4ExcitedStringDecay *theStringDecay =
new G4ExcitedStringDecay(theFrag);
184 theStringModel->SetFragmentationModel(theStringDecay);
186 G4GeneratorPrecompoundInterface *theCascade =
new G4GeneratorPrecompoundInterface;
188 G4TheoFSGenerator *theModel =
new G4TheoFSGenerator;
189 theModel->SetTransport(theCascade);
190 theModel->SetHighEnergyGenerator(theStringModel);
192 G4HadronicParameters *param = G4HadronicParameters::Instance();
194 G4CascadeInterface *cascade =
new G4CascadeInterface;
197 G4double gnLowEnergyLimit = 200 * CLHEP::MeV;
198 if (gnLowEnergyLimit > 0.0) {
199 G4LowEGammaNuclearModel *lemod =
new G4LowEGammaNuclearModel;
200 lemod->SetMaxEnergy(gnLowEnergyLimit);
201 nuc->RegisterMe(lemod);
202 cascade->SetMinEnergy(gnLowEnergyLimit - CLHEP::MeV);
204 cascade->SetMaxEnergy(param->GetMaxEnergyTransitionFTF_Cascade());
205 nuc->RegisterMe(cascade);
206 theModel->SetMinEnergy(param->GetMinEnergyTransitionFTF_Cascade());
207 theModel->SetMaxEnergy(param->GetMaxEnergy());
208 nuc->RegisterMe(theModel);
210 gammaProc->AddHadProcess(nuc);
211 G4LossTableManager::Instance()->SetGammaGeneralProcess(gammaProc);
216 G4ElectroVDNuclearModel *eModel =
new G4ElectroVDNuclearModel;
219 G4ElectronNuclearProcess *nuc =
new G4ElectronNuclearProcess;
220 nuc->RegisterMe(eModel);
224 G4PositronNuclearProcess *nuc =
new G4PositronNuclearProcess;
225 nuc->RegisterMe(eModel);
229 if (masterTrackingManager ==
nullptr) {
230 masterTrackingManager =
this;
232 electron.msc->SetMasterProcess(masterTrackingManager->electron.msc);
233 electron.ss->SetMasterProcess(masterTrackingManager->electron.ss);
234 electron.ioni->SetMasterProcess(masterTrackingManager->electron.ioni);
235 electron.brems->SetMasterProcess(masterTrackingManager->electron.brems);
236 electron.nuc->SetMasterProcess(masterTrackingManager->electron.nuc);
238 positron.msc->SetMasterProcess(masterTrackingManager->positron.msc);
239 positron.ss->SetMasterProcess(masterTrackingManager->positron.ss);
240 positron.ioni->SetMasterProcess(masterTrackingManager->positron.ioni);
241 positron.brems->SetMasterProcess(masterTrackingManager->positron.brems);
242 positron.annihilation->SetMasterProcess(masterTrackingManager->positron.annihilation);
243 positron.nuc->SetMasterProcess(masterTrackingManager->positron.nuc);
245 gammaProc->SetMasterProcess(masterTrackingManager->gammaProc);
249 CMSEmStandardPhysicsTrackingManager::~CMSEmStandardPhysicsTrackingManager() {
250 if (masterTrackingManager ==
this) {
251 masterTrackingManager =
nullptr;
255 void CMSEmStandardPhysicsTrackingManager::BuildPhysicsTable(
const G4ParticleDefinition &
part) {
256 if (&
part == G4Electron::Definition()) {
262 }
else if (&
part == G4Positron::Definition()) {
263 positron.msc->BuildPhysicsTable(
part);
264 positron.ioni->BuildPhysicsTable(
part);
265 positron.brems->BuildPhysicsTable(
part);
266 positron.annihilation->BuildPhysicsTable(
part);
267 positron.ss->BuildPhysicsTable(
part);
268 positron.nuc->BuildPhysicsTable(
part);
269 }
else if (&
part == G4Gamma::Definition()) {
270 gammaProc->BuildPhysicsTable(
part);
274 void CMSEmStandardPhysicsTrackingManager::PreparePhysicsTable(
const G4ParticleDefinition &
part) {
275 if (&
part == G4Electron::Definition()) {
281 }
else if (&
part == G4Positron::Definition()) {
282 positron.msc->PreparePhysicsTable(
part);
283 positron.ioni->PreparePhysicsTable(
part);
284 positron.brems->PreparePhysicsTable(
part);
285 positron.annihilation->PreparePhysicsTable(
part);
286 positron.ss->PreparePhysicsTable(
part);
287 positron.nuc->PreparePhysicsTable(
part);
288 }
else if (&
part == G4Gamma::Definition()) {
289 gammaProc->PreparePhysicsTable(
part);
293 void CMSEmStandardPhysicsTrackingManager::TrackElectron(G4Track *aTrack) {
296 ElectronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
301 electron.msc->StartTracking(aTrack);
302 electron.ioni->StartTracking(aTrack);
303 electron.brems->StartTracking(aTrack);
305 electron.nuc->StartTracking(aTrack);
307 fPreviousStepLength = 0;
321 G4double physIntLength, proposedSafety = DBL_MAX;
322 G4ForceCondition condition;
325 fProposedStep = DBL_MAX;
328 physIntLength =
electron.nuc->PostStepGPIL(
track, fPreviousStepLength, &condition);
329 if (physIntLength < fProposedStep) {
330 fProposedStep = physIntLength;
334 physIntLength =
electron.ss->PostStepGPIL(
track, fPreviousStepLength, &condition);
335 if (physIntLength < fProposedStep) {
336 fProposedStep = physIntLength;
340 physIntLength =
electron.brems->PostStepGPIL(
track, fPreviousStepLength, &condition);
341 if (physIntLength < fProposedStep) {
342 fProposedStep = physIntLength;
346 physIntLength =
electron.ioni->PostStepGPIL(
track, fPreviousStepLength, &condition);
347 if (physIntLength < fProposedStep) {
348 fProposedStep = physIntLength;
354 if (physIntLength < fProposedStep) {
355 fProposedStep = physIntLength;
361 if (physIntLength < fProposedStep) {
362 fProposedStep = physIntLength;
365 if (
selection == CandidateForSelection) {
370 return fProposedStep;
374 if (
step.GetStepLength() == fProposedStep) {
375 step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
381 G4VParticleChange *particleChange;
384 particleChange->UpdateStepForAlongStep(&
step);
385 track.SetTrackStatus(particleChange->GetTrackStatus());
386 particleChange->Clear();
389 particleChange->UpdateStepForAlongStep(&
step);
390 track.SetTrackStatus(particleChange->GetTrackStatus());
391 particleChange->Clear();
393 fPreviousStepLength =
step.GetStepLength();
400 step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
404 G4VParticleChange *particleChange =
nullptr;
425 particleChange->UpdateStepForPostStep(&
step);
428 int numSecondaries = particleChange->GetNumberOfSecondaries();
429 for (
int i = 0;
i < numSecondaries;
i++) {
430 G4Track *secondary = particleChange->GetSecondary(
i);
431 secondary->SetParentID(
track.GetTrackID());
432 secondary->SetCreatorProcess(
process);
433 secondaries.push_back(secondary);
436 track.SetTrackStatus(particleChange->GetTrackStatus());
437 particleChange->Clear();
441 CMSEmStandardPhysicsTrackingManager &fMgr;
442 G4double fPreviousStepLength;
443 G4double fProposedStep;
447 ElectronPhysics
physics(*
this);
451 void CMSEmStandardPhysicsTrackingManager::TrackPositron(G4Track *aTrack) {
454 PositronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
457 auto &positron = fMgr.positron;
459 positron.msc->StartTracking(aTrack);
460 positron.ioni->StartTracking(aTrack);
461 positron.brems->StartTracking(aTrack);
462 positron.annihilation->StartTracking(aTrack);
463 positron.ss->StartTracking(aTrack);
464 positron.nuc->StartTracking(aTrack);
466 fPreviousStepLength = 0;
469 auto &positron = fMgr.positron;
471 positron.msc->EndTracking();
472 positron.ioni->EndTracking();
473 positron.brems->EndTracking();
474 positron.annihilation->EndTracking();
475 positron.ss->EndTracking();
476 positron.nuc->EndTracking();
480 auto &positron = fMgr.positron;
481 G4double physIntLength, proposedSafety = DBL_MAX;
482 G4ForceCondition condition;
485 fProposedStep = DBL_MAX;
488 physIntLength = positron.nuc->PostStepGPIL(
track, fPreviousStepLength, &condition);
489 if (physIntLength < fProposedStep) {
490 fProposedStep = physIntLength;
494 physIntLength = positron.ss->PostStepGPIL(
track, fPreviousStepLength, &condition);
495 if (physIntLength < fProposedStep) {
496 fProposedStep = physIntLength;
500 physIntLength = positron.annihilation->PostStepGPIL(
track, fPreviousStepLength, &condition);
501 if (physIntLength < fProposedStep) {
502 fProposedStep = physIntLength;
506 physIntLength = positron.brems->PostStepGPIL(
track, fPreviousStepLength, &condition);
507 if (physIntLength < fProposedStep) {
508 fProposedStep = physIntLength;
512 physIntLength = positron.ioni->PostStepGPIL(
track, fPreviousStepLength, &condition);
513 if (physIntLength < fProposedStep) {
514 fProposedStep = physIntLength;
519 positron.ioni->AlongStepGPIL(
track, fPreviousStepLength, fProposedStep, proposedSafety, &
selection);
520 if (physIntLength < fProposedStep) {
521 fProposedStep = physIntLength;
526 positron.msc->AlongStepGPIL(
track, fPreviousStepLength, fProposedStep, proposedSafety, &
selection);
527 if (physIntLength < fProposedStep) {
528 fProposedStep = physIntLength;
531 if (
selection == CandidateForSelection) {
536 return fProposedStep;
540 if (
step.GetStepLength() == fProposedStep) {
541 step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
546 auto &positron = fMgr.positron;
547 G4VParticleChange *particleChange;
549 particleChange = positron.msc->AlongStepDoIt(
track,
step);
550 particleChange->UpdateStepForAlongStep(&
step);
551 track.SetTrackStatus(particleChange->GetTrackStatus());
552 particleChange->Clear();
554 particleChange = positron.ioni->AlongStepDoIt(
track,
step);
555 particleChange->UpdateStepForAlongStep(&
step);
556 track.SetTrackStatus(particleChange->GetTrackStatus());
557 particleChange->Clear();
559 fPreviousStepLength =
step.GetStepLength();
566 step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
568 auto &positron = fMgr.positron;
570 G4VParticleChange *particleChange =
nullptr;
575 particleChange = positron.nuc->PostStepDoIt(
track,
step);
579 particleChange = positron.ss->PostStepDoIt(
track,
step);
582 process = positron.annihilation;
583 particleChange = positron.annihilation->PostStepDoIt(
track,
step);
587 particleChange = positron.brems->PostStepDoIt(
track,
step);
591 particleChange = positron.ioni->PostStepDoIt(
track,
step);
595 particleChange->UpdateStepForPostStep(&
step);
598 int numSecondaries = particleChange->GetNumberOfSecondaries();
599 for (
int i = 0;
i < numSecondaries;
i++) {
600 G4Track *secondary = particleChange->GetSecondary(
i);
601 secondary->SetParentID(
track.GetTrackID());
602 secondary->SetCreatorProcess(
process);
603 secondaries.push_back(secondary);
606 track.SetTrackStatus(particleChange->GetTrackStatus());
607 particleChange->Clear();
613 auto &positron = fMgr.positron;
615 G4VParticleChange *particleChange = positron.annihilation->AtRestDoIt(
track,
step);
616 particleChange->UpdateStepForAtRest(&
step);
619 int numSecondaries = particleChange->GetNumberOfSecondaries();
620 for (
int i = 0;
i < numSecondaries;
i++) {
621 G4Track *secondary = particleChange->GetSecondary(
i);
622 secondary->SetParentID(
track.GetTrackID());
623 secondary->SetCreatorProcess(positron.annihilation);
624 secondaries.push_back(secondary);
627 track.SetTrackStatus(particleChange->GetTrackStatus());
628 particleChange->Clear();
632 CMSEmStandardPhysicsTrackingManager &fMgr;
633 G4double fPreviousStepLength;
634 G4double fProposedStep;
638 PositronPhysics
physics(*
this);
642 void CMSEmStandardPhysicsTrackingManager::TrackGamma(G4Track *aTrack) {
645 GammaPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
648 fMgr.gammaProc->StartTracking(aTrack);
650 fPreviousStepLength = 0;
652 void EndTracking()
override { fMgr.gammaProc->EndTracking(); }
655 G4double physIntLength;
656 G4ForceCondition condition;
658 fProposedStep = DBL_MAX;
661 physIntLength = fMgr.gammaProc->PostStepGPIL(
track, fPreviousStepLength, &condition);
662 if (physIntLength < fProposedStep) {
663 fProposedStep = physIntLength;
667 return fProposedStep;
671 if (
step.GetStepLength() == fProposedStep) {
672 step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
677 fPreviousStepLength =
step.GetStepLength();
684 step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
686 G4VProcess *
process = fMgr.gammaProc;
687 G4VParticleChange *particleChange = fMgr.gammaProc->PostStepDoIt(
track,
step);
689 particleChange->UpdateStepForPostStep(&
step);
692 int numSecondaries = particleChange->GetNumberOfSecondaries();
693 for (
int i = 0;
i < numSecondaries;
i++) {
694 G4Track *secondary = particleChange->GetSecondary(
i);
695 secondary->SetParentID(
track.GetTrackID());
696 secondary->SetCreatorProcess(
process);
697 secondaries.push_back(secondary);
700 track.SetTrackStatus(particleChange->GetTrackStatus());
701 particleChange->Clear();
705 CMSEmStandardPhysicsTrackingManager &fMgr;
706 G4double fPreviousStepLength;
707 G4double fProposedStep;
715 void CMSEmStandardPhysicsTrackingManager::HandOverOneTrack(G4Track *aTrack) {
716 const G4ParticleDefinition *
part = aTrack->GetParticleDefinition();
718 if (
part == G4Electron::Definition()) {
719 TrackElectron(aTrack);
720 }
else if (
part == G4Positron::Definition()) {
721 TrackPositron(aTrack);
722 }
else if (
part == G4Gamma::Definition()) {
726 aTrack->SetTrackStatus(fStopAndKill);
static void TrackNeutralParticle(G4Track *aTrack, PhysicsImpl &physics)
virtual void StartTracking(G4Track *)
static void TrackChargedParticle(G4Track *aTrack, PhysicsImpl &physics)
virtual void EndTracking()
virtual void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)=0
virtual bool HasAtRestProcesses()
virtual void AtRestDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)
virtual void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)=0
virtual G4double GetPhysicalInteractionLength(const G4Track &track)=0