File indexing completed on 2025-02-23 09:21:14
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031 #include "G4EmPenelopePhysicsMI.hh"
0032
0033 #include "G4ParticleDefinition.hh"
0034 #include "G4ParticleTable.hh"
0035 #include "G4SystemOfUnits.hh"
0036
0037
0038
0039 #include "G4ComptonScattering.hh"
0040 #include "G4GammaConversion.hh"
0041 #include "G4PenelopeComptonModel.hh"
0042 #include "G4PenelopeGammaConversionModel.hh"
0043 #include "G4PenelopePhotoElectricModel.hh"
0044 #include "G4PenelopeRayleighModel.hh"
0045 #include "G4PenelopeRayleighModelMI.hh"
0046 #include "G4PhotoElectricEffect.hh"
0047 #include "G4RayleighScattering.hh"
0048
0049
0050 #include "G4PenelopeBremsstrahlungModel.hh"
0051 #include "G4PenelopeIonisationModel.hh"
0052 #include "G4UniversalFluctuation.hh"
0053 #include "G4eBremsstrahlung.hh"
0054 #include "G4eIonisation.hh"
0055 #include "G4eMultipleScattering.hh"
0056
0057
0058 #include "G4PenelopeAnnihilationModel.hh"
0059 #include "G4eplusAnnihilation.hh"
0060
0061
0062 #include "G4MuBremsstrahlung.hh"
0063 #include "G4MuBremsstrahlungModel.hh"
0064 #include "G4MuIonisation.hh"
0065 #include "G4MuMultipleScattering.hh"
0066 #include "G4MuPairProduction.hh"
0067 #include "G4MuPairProductionModel.hh"
0068 #include "G4hBremsstrahlungModel.hh"
0069 #include "G4hPairProductionModel.hh"
0070
0071
0072 #include "G4IonParametrisedLossModel.hh"
0073 #include "G4MscStepLimitType.hh"
0074 #include "G4NuclearStopping.hh"
0075 #include "G4hBremsstrahlung.hh"
0076 #include "G4hIonisation.hh"
0077 #include "G4hMultipleScattering.hh"
0078 #include "G4hPairProduction.hh"
0079 #include "G4ionIonisation.hh"
0080
0081
0082 #include "G4CoulombScattering.hh"
0083 #include "G4GoudsmitSaundersonMscModel.hh"
0084 #include "G4LossTableManager.hh"
0085 #include "G4UAtomicDeexcitation.hh"
0086 #include "G4UrbanMscModel.hh"
0087 #include "G4VAtomDeexcitation.hh"
0088 #include "G4WentzelVIModel.hh"
0089 #include "G4eCoulombScatteringModel.hh"
0090
0091
0092 #include "G4Alpha.hh"
0093 #include "G4AntiProton.hh"
0094 #include "G4BuilderType.hh"
0095 #include "G4Deuteron.hh"
0096 #include "G4Electron.hh"
0097 #include "G4EmModelActivator.hh"
0098 #include "G4Gamma.hh"
0099 #include "G4GenericIon.hh"
0100 #include "G4He3.hh"
0101 #include "G4KaonMinus.hh"
0102 #include "G4KaonPlus.hh"
0103 #include "G4MuonMinus.hh"
0104 #include "G4MuonPlus.hh"
0105 #include "G4PhysicsListHelper.hh"
0106 #include "G4PionMinus.hh"
0107 #include "G4PionPlus.hh"
0108 #include "G4Positron.hh"
0109 #include "G4Proton.hh"
0110 #include "G4Triton.hh"
0111
0112
0113 #include "G4PhysicsConstructorFactory.hh"
0114
0115 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmPenelopePhysicsMI);
0116
0117
0118
0119 G4EmPenelopePhysicsMI::G4EmPenelopePhysicsMI(G4int ver, const G4String&, G4bool UseMIFlag)
0120 : G4VPhysicsConstructor("G4EmPenelopeMI"), fVerbose(ver), fUseMIFlag(UseMIFlag)
0121 {
0122 G4EmParameters* param = G4EmParameters::Instance();
0123 param->SetDefaults();
0124 param->SetVerbose(fVerbose);
0125 param->SetMinEnergy(100 * eV);
0126 param->SetLowestElectronEnergy(100 * eV);
0127 param->SetNumberOfBinsPerDecade(20);
0128 param->SetMscRangeFactor(0.02);
0129 param->SetMscStepLimitType(fUseDistanceToBoundary);
0130 param->SetMuHadLateralDisplacement(true);
0131 param->SetFluo(true);
0132 param->SetPIXEElectronCrossSectionModel("Penelope");
0133 SetPhysicsType(bElectromagnetic);
0134 }
0135
0136
0137
0138 G4EmPenelopePhysicsMI::~G4EmPenelopePhysicsMI() {}
0139
0140
0141
0142 void G4EmPenelopePhysicsMI::ConstructParticle()
0143 {
0144
0145 G4Gamma::Gamma();
0146
0147
0148 G4Electron::Electron();
0149 G4Positron::Positron();
0150 G4MuonPlus::MuonPlus();
0151 G4MuonMinus::MuonMinus();
0152
0153
0154 G4PionPlus::PionPlusDefinition();
0155 G4PionMinus::PionMinusDefinition();
0156 G4KaonPlus::KaonPlusDefinition();
0157 G4KaonMinus::KaonMinusDefinition();
0158
0159
0160 G4Proton::Proton();
0161 G4AntiProton::AntiProton();
0162
0163
0164 G4Deuteron::Deuteron();
0165 G4Triton::Triton();
0166 G4He3::He3();
0167 G4Alpha::Alpha();
0168 G4GenericIon::GenericIonDefinition();
0169 }
0170
0171
0172
0173 void G4EmPenelopePhysicsMI::ConstructProcess()
0174 {
0175 if (fVerbose > 1) {
0176 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
0177 }
0178
0179 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
0180
0181
0182 G4MuBremsstrahlung* mub = new G4MuBremsstrahlung();
0183 G4MuPairProduction* mup = new G4MuPairProduction();
0184 G4hBremsstrahlung* pib = new G4hBremsstrahlung();
0185 G4hPairProduction* pip = new G4hPairProduction();
0186 G4hBremsstrahlung* kb = new G4hBremsstrahlung();
0187 G4hPairProduction* kp = new G4hPairProduction();
0188 G4hBremsstrahlung* pb = new G4hBremsstrahlung();
0189 G4hPairProduction* pp = new G4hPairProduction();
0190
0191
0192 G4MuMultipleScattering* mumsc = new G4MuMultipleScattering();
0193 mumsc->SetEmModel(new G4WentzelVIModel());
0194 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
0195
0196
0197 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
0198
0199
0200 G4NuclearStopping* pnuc = new G4NuclearStopping();
0201
0202
0203
0204 G4double PenelopeHighEnergyLimit = 1.0 * GeV;
0205
0206
0207 G4ParticleTable* table = G4ParticleTable::GetParticleTable();
0208 for (const auto& particleName : fPartList.PartNames()) {
0209 G4ParticleDefinition* particle = table->FindParticle(particleName);
0210 if (!particle) {
0211 continue;
0212 }
0213 if (particleName == "gamma") {
0214
0215 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
0216 G4PenelopePhotoElectricModel* thePEPenelopeModel = new G4PenelopePhotoElectricModel();
0217 thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0218 thePhotoElectricEffect->SetEmModel(thePEPenelopeModel);
0219 ph->RegisterProcess(thePhotoElectricEffect, particle);
0220
0221
0222 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
0223 G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
0224 theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0225 theComptonScattering->SetEmModel(theComptonPenelopeModel);
0226 ph->RegisterProcess(theComptonScattering, particle);
0227
0228
0229 G4GammaConversion* theGammaConversion = new G4GammaConversion();
0230 G4PenelopeGammaConversionModel* theGCPenelopeModel = new G4PenelopeGammaConversionModel();
0231 theGammaConversion->SetEmModel(theGCPenelopeModel);
0232 ph->RegisterProcess(theGammaConversion, particle);
0233
0234
0235 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
0236 G4PenelopeRayleighModelMI* theRayleighPenelopeModel = new G4PenelopeRayleighModelMI();
0237 theRayleighPenelopeModel->SetVerbosityLevel(1);
0238 theRayleighPenelopeModel->SetMIActive(fUseMIFlag);
0239
0240 theRayleigh->SetEmModel(theRayleighPenelopeModel);
0241 ph->RegisterProcess(theRayleigh, particle);
0242 }
0243 else if (particleName == "e-") {
0244
0245 G4eMultipleScattering* msc = new G4eMultipleScattering;
0246 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
0247 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
0248 msc1->SetHighEnergyLimit(highEnergyLimit);
0249 msc2->SetLowEnergyLimit(highEnergyLimit);
0250 msc->SetEmModel(msc1);
0251 msc->SetEmModel(msc2);
0252
0253 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
0254 G4CoulombScattering* ss = new G4CoulombScattering();
0255 ss->SetEmModel(ssm);
0256 ss->SetMinKinEnergy(highEnergyLimit);
0257 ssm->SetLowEnergyLimit(highEnergyLimit);
0258 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
0259
0260
0261 G4eIonisation* eIoni = new G4eIonisation();
0262 G4PenelopeIonisationModel* theIoniPenelope = new G4PenelopeIonisationModel();
0263 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0264 eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
0265 eIoni->SetStepFunction(0.2, 100 * um);
0266
0267
0268 G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
0269 G4PenelopeBremsstrahlungModel* theBremPenelope = new G4PenelopeBremsstrahlungModel();
0270 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0271 eBrem->SetEmModel(theBremPenelope);
0272
0273
0274 ph->RegisterProcess(msc, particle);
0275 ph->RegisterProcess(eIoni, particle);
0276 ph->RegisterProcess(eBrem, particle);
0277 ph->RegisterProcess(ss, particle);
0278 }
0279 else if (particleName == "e+") {
0280
0281 G4eMultipleScattering* msc = new G4eMultipleScattering;
0282 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
0283 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
0284 msc1->SetHighEnergyLimit(highEnergyLimit);
0285 msc2->SetLowEnergyLimit(highEnergyLimit);
0286 msc->SetEmModel(msc1);
0287 msc->SetEmModel(msc2);
0288
0289 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
0290 G4CoulombScattering* ss = new G4CoulombScattering();
0291 ss->SetEmModel(ssm);
0292 ss->SetMinKinEnergy(highEnergyLimit);
0293 ssm->SetLowEnergyLimit(highEnergyLimit);
0294 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
0295
0296
0297 G4eIonisation* eIoni = new G4eIonisation();
0298 G4PenelopeIonisationModel* theIoniPenelope = new G4PenelopeIonisationModel();
0299 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0300 eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
0301 eIoni->SetStepFunction(0.2, 100 * um);
0302
0303
0304 G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
0305 G4PenelopeBremsstrahlungModel* theBremPenelope = new G4PenelopeBremsstrahlungModel();
0306 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0307 eBrem->SetEmModel(theBremPenelope);
0308
0309
0310 G4eplusAnnihilation* eAnni = new G4eplusAnnihilation();
0311 G4PenelopeAnnihilationModel* theAnnPenelope = new G4PenelopeAnnihilationModel();
0312 theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
0313 eAnni->AddEmModel(0, theAnnPenelope);
0314
0315
0316 ph->RegisterProcess(msc, particle);
0317 ph->RegisterProcess(eIoni, particle);
0318 ph->RegisterProcess(eBrem, particle);
0319 ph->RegisterProcess(eAnni, particle);
0320 ph->RegisterProcess(ss, particle);
0321 }
0322 else if (particleName == "mu+" || particleName == "mu-") {
0323 G4MuIonisation* muIoni = new G4MuIonisation();
0324 muIoni->SetStepFunction(0.2, 50 * um);
0325
0326 ph->RegisterProcess(mumsc, particle);
0327 ph->RegisterProcess(muIoni, particle);
0328 ph->RegisterProcess(mub, particle);
0329 ph->RegisterProcess(mup, particle);
0330 ph->RegisterProcess(new G4CoulombScattering(), particle);
0331 }
0332 else if (particleName == "alpha" || particleName == "He3") {
0333 G4hMultipleScattering* msc = new G4hMultipleScattering();
0334 G4ionIonisation* ionIoni = new G4ionIonisation();
0335 ionIoni->SetStepFunction(0.1, 10 * um);
0336
0337 ph->RegisterProcess(msc, particle);
0338 ph->RegisterProcess(ionIoni, particle);
0339 ph->RegisterProcess(pnuc, particle);
0340 }
0341 else if (particleName == "GenericIon") {
0342 G4ionIonisation* ionIoni = new G4ionIonisation();
0343 ionIoni->SetEmModel(new G4IonParametrisedLossModel());
0344 ionIoni->SetStepFunction(0.1, 1 * um);
0345
0346 ph->RegisterProcess(hmsc, particle);
0347 ph->RegisterProcess(ionIoni, particle);
0348 ph->RegisterProcess(pnuc, particle);
0349 }
0350 else if (particleName == "pi+" || particleName == "pi-") {
0351 G4hMultipleScattering* pimsc = new G4hMultipleScattering();
0352 G4hIonisation* hIoni = new G4hIonisation();
0353 hIoni->SetStepFunction(0.2, 50 * um);
0354
0355 ph->RegisterProcess(pimsc, particle);
0356 ph->RegisterProcess(hIoni, particle);
0357 ph->RegisterProcess(pib, particle);
0358 ph->RegisterProcess(pip, particle);
0359 }
0360 else if (particleName == "kaon+" || particleName == "kaon-") {
0361 G4hMultipleScattering* kmsc = new G4hMultipleScattering();
0362 G4hIonisation* hIoni = new G4hIonisation();
0363 hIoni->SetStepFunction(0.2, 50 * um);
0364
0365 ph->RegisterProcess(kmsc, particle);
0366 ph->RegisterProcess(hIoni, particle);
0367 ph->RegisterProcess(kb, particle);
0368 ph->RegisterProcess(kp, particle);
0369 }
0370 else if (particleName == "proton" || particleName == "anti_proton") {
0371 G4hMultipleScattering* pmsc = new G4hMultipleScattering();
0372 G4hIonisation* hIoni = new G4hIonisation();
0373 hIoni->SetStepFunction(0.2, 50 * um);
0374
0375 ph->RegisterProcess(pmsc, particle);
0376 ph->RegisterProcess(hIoni, particle);
0377 ph->RegisterProcess(pb, particle);
0378 ph->RegisterProcess(pp, particle);
0379 ph->RegisterProcess(pnuc, particle);
0380 }
0381 else if (particleName == "B+" || particleName == "B-" || particleName == "D+"
0382 || particleName == "D-" || particleName == "Ds+" || particleName == "Ds-"
0383 || particleName == "anti_He3" || particleName == "anti_alpha"
0384 || particleName == "anti_deuteron" || particleName == "anti_lambda_c+"
0385 || particleName == "anti_omega-" || particleName == "anti_sigma_c+"
0386 || particleName == "anti_sigma_c++" || particleName == "anti_sigma+"
0387 || particleName == "anti_sigma-" || particleName == "anti_triton"
0388 || particleName == "anti_xi_c+" || particleName == "anti_xi-"
0389 || particleName == "deuteron" || particleName == "lambda_c+"
0390 || particleName == "omega-" || particleName == "sigma_c+"
0391 || particleName == "sigma_c++" || particleName == "sigma+" || particleName == "sigma-"
0392 || particleName == "tau+" || particleName == "tau-" || particleName == "triton"
0393 || particleName == "xi_c+" || particleName == "xi-")
0394 {
0395 ph->RegisterProcess(hmsc, particle);
0396 ph->RegisterProcess(new G4hIonisation(), particle);
0397 }
0398 }
0399
0400
0401 pnuc->SetMaxKinEnergy(MeV);
0402
0403
0404 G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
0405 G4LossTableManager::Instance()->SetAtomDeexcitation(deexcitation);
0406
0407 G4EmModelActivator mact(GetPhysicsName());
0408 }
0409
0410