File indexing completed on 2026-04-04 07:52:55
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 #include "PhysicsList.hh"
0030
0031 #include "DetectorConstruction.hh"
0032 #include "PhysicsListMessenger.hh"
0033
0034 #include "G4PhysicsConstructorRegistry.hh"
0035 #include "G4RunManager.hh"
0036 #include "G4SystemOfUnits.hh"
0037
0038 #include "G4ComptonScattering.hh"
0039 #include "G4DNAAttachment.hh"
0040 #include "G4DNABornExcitationModel.hh"
0041 #include "G4DNABornIonisationModel.hh"
0042 #include "G4DNAChampionElasticModel.hh"
0043 #include "G4DNAChemistryManager.hh"
0044 #include "G4DNADingfelderChargeDecreaseModel.hh"
0045 #include "G4DNADiracRMatrixExcitationModel.hh"
0046 #include "G4DNAELSEPAElasticModel.hh"
0047 #include "G4DNAElastic.hh"
0048 #include "G4DNAExcitation.hh"
0049 #include "G4DNAGenericIonsManager.hh"
0050 #include "G4DNAIonisation.hh"
0051 #include "G4DNAMeltonAttachmentModel.hh"
0052 #include "G4DNAPlasmonExcitation.hh"
0053 #include "G4DNAQuinnPlasmonExcitationModel.hh"
0054 #include "G4DNARelativisticIonisationModel.hh"
0055 #include "G4DNASancheExcitationModel.hh"
0056 #include "G4DNAVibExcitation.hh"
0057 #include "G4DummyModel.hh"
0058 #include "G4EmConfigurator.hh"
0059 #include "G4Gamma.hh"
0060 #include "G4GammaConversion.hh"
0061 #include "G4GoudsmitSaundersonMscModel.hh"
0062 #include "G4IonFluctuations.hh"
0063 #include "G4LivermoreBremsstrahlungModel.hh"
0064 #include "G4LivermoreComptonModel.hh"
0065 #include "G4LivermoreGammaConversionModel.hh"
0066 #include "G4LivermoreIonisationModel.hh"
0067 #include "G4LivermorePhotoElectricModel.hh"
0068 #include "G4LossTableManager.hh"
0069 #include "G4PenelopeBremsstrahlungModel.hh"
0070 #include "G4PenelopeIonisationModel.hh"
0071 #include "G4PhotoElectricEffect.hh"
0072 #include "G4RayleighScattering.hh"
0073 #include "G4SeltzerBergerModel.hh"
0074 #include "G4StepLimiter.hh"
0075 #include "G4UAtomicDeexcitation.hh"
0076 #include "G4UniversalFluctuation.hh"
0077 #include "G4UrbanMscModel.hh"
0078 #include "G4VEmModel.hh"
0079 #include "G4eBremsstrahlung.hh"
0080 #include "G4eIonisation.hh"
0081 #include "G4eMultipleScattering.hh"
0082 #include "G4hMultipleScattering.hh"
0083
0084 #include "G4ProcessTable.hh"
0085 #include "G4EmDNAChemistry_option3.hh"
0086 #include "G4GenericIon.hh"
0087
0088
0089 PhysicsList::PhysicsList() {
0090 fpDetector = dynamic_cast<const DetectorConstruction *>(
0091 G4RunManager::GetRunManager()->GetUserDetectorConstruction());
0092 defaultCutValue = 0.1 * CLHEP::nanometer;
0093 fPhysMessenger = std::make_unique<PhysicsListMessenger>(this);
0094 SetVerboseLevel(1);
0095 fEmDNAChemistryList = std::make_unique<G4EmDNAChemistry_option3>();
0096 G4EmParameters *param = G4EmParameters::Instance();
0097 param->SetMinEnergy(100 * eV);
0098 param->SetMaxEnergy(1 * GeV);
0099 }
0100
0101
0102
0103 void PhysicsList::ConstructParticle() {
0104 ConstructBosons();
0105 ConstructLeptons();
0106 ConstructBarions();
0107 if (fEmDNAChemistryList != nullptr) {
0108 fEmDNAChemistryList->ConstructParticle();
0109 }
0110 G4VModularPhysicsList::ConstructParticle();
0111 }
0112
0113
0114
0115 void PhysicsList::ConstructBosons() {
0116 G4Gamma::GammaDefinition();
0117 }
0118
0119
0120
0121 void PhysicsList::ConstructLeptons() {
0122 G4Electron::ElectronDefinition();
0123 G4Positron::PositronDefinition();
0124 }
0125
0126
0127
0128 void PhysicsList::ConstructBarions() {
0129 G4Proton::ProtonDefinition();
0130 G4GenericIon::GenericIonDefinition();
0131
0132 auto *genericIonsManager = G4DNAGenericIonsManager::Instance();
0133 genericIonsManager->GetIon("alpha++");
0134 genericIonsManager->GetIon("alpha+");
0135 genericIonsManager->GetIon("helium");
0136 genericIonsManager->GetIon("hydrogen");
0137 }
0138
0139
0140
0141 void PhysicsList::ConstructProcess() {
0142 ConstructEM();
0143 ConstructGeneral();
0144 if (fEmDNAChemistryList != nullptr) {
0145 fEmDNAChemistryList->ConstructProcess();
0146 }
0147
0148 G4VModularPhysicsList::ConstructProcess();
0149 }
0150
0151
0152
0153 void PhysicsList::ConstructEM() {
0154 G4bool isGNP = false;
0155 if (fphysname == "") fphysname = "DNA";
0156
0157 const G4Material *NPMaterial = fpDetector->GetNPMaterial();
0158 const G4String npname = NPMaterial->GetName();
0159
0160 if (npname == "G4_Au") isGNP = true;
0161
0162 const auto theParticleIterator = GetParticleIterator();
0163 theParticleIterator->reset();
0164
0165 while ((*theParticleIterator)()) {
0166 const G4ParticleDefinition *particle = theParticleIterator->value();
0167
0168 G4ProcessManager *pm = particle->GetProcessManager();
0169
0170 G4String particleName = particle->GetParticleName();
0171
0172 if (particleName == "e-") {
0173 const auto theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
0174 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel(), 1);
0175 pm->AddDiscreteProcess(theDNAElasticProcess);
0176
0177 const auto theDNAExcitationProcess = new G4DNAExcitation("e-_G4DNAExcitation");
0178 theDNAExcitationProcess->SetEmModel(new G4DNABornExcitationModel(), 1);
0179 pm->AddDiscreteProcess(theDNAExcitationProcess);
0180
0181 const auto theDNAIonizationProcess = new G4DNAIonisation("e-_G4DNAIonisation");
0182 theDNAIonizationProcess->SetEmModel(new G4DNABornIonisationModel(), 1);
0183 pm->AddDiscreteProcess(theDNAIonizationProcess);
0184
0185 const auto theDNAAttachmentProcess = new G4DNAAttachment("e-_G4DNAAttachment");
0186 pm->AddDiscreteProcess(theDNAAttachmentProcess);
0187
0188 const auto theDNAVibExcProcess = new G4DNAVibExcitation("e-_G4DNAVibExcitation");
0189 pm->AddDiscreteProcess(theDNAVibExcProcess);
0190
0191 const auto theDNABremProcess = new G4eBremsstrahlung("e-_G4DNABremsstrahlung");
0192 pm->AddDiscreteProcess(theDNABremProcess);
0193
0194 if (isGNP) {
0195 if (fphysname == "DNA") {
0196 const auto theDNAELSEPAElasticProcess = new G4DNAElastic("e-_G4DNAELSEPAElastic");
0197 theDNAELSEPAElasticProcess->SetEmModel(new G4DummyModel(), 1);
0198 pm->AddDiscreteProcess(theDNAELSEPAElasticProcess);
0199 const auto theDNADRMExcitationProcess =
0200 new G4DNAExcitation("e-_G4DNADRMExcitation");
0201 theDNADRMExcitationProcess->SetEmModel(new G4DummyModel(), 1);
0202 pm->AddDiscreteProcess(theDNADRMExcitationProcess);
0203 const auto theDNARelativisticIonizationProcess =
0204 new G4DNAIonisation("e-_G4DNARelativisticIonisation");
0205 theDNARelativisticIonizationProcess->SetEmModel(new G4DummyModel(), 1);
0206 pm->AddDiscreteProcess(theDNARelativisticIonizationProcess);
0207 const auto theDNAPExcitationProcess =
0208 new G4DNAPlasmonExcitation("e-_G4DNAPlasmonExcitation");
0209 theDNAPExcitationProcess->SetEmModel(new G4DummyModel(), 1);
0210 pm->AddDiscreteProcess(theDNAPExcitationProcess);
0211 const auto theDNAeBremProcess =
0212 new G4eBremsstrahlung("e-_G4DNABremsstrahlung_GNP");
0213 theDNAeBremProcess->SetEmModel(new G4DummyModel(), 1);
0214 pm->AddDiscreteProcess(theDNAeBremProcess);
0215 }
0216 if (fphysname == "Livermore") {
0217 const auto msc = new G4eMultipleScattering();
0218 msc->SetEmModel(new G4DummyModel(), 1);
0219 pm->AddDiscreteProcess(msc);
0220 const auto eIoni = new G4eIonisation("e-_G4LivermoreIoni");
0221 eIoni->SetEmModel(new G4DummyModel(), 1);
0222 pm->AddDiscreteProcess(eIoni);
0223 const auto eBrem = new G4eBremsstrahlung("e-_G4LivermoreBrem");
0224 eBrem->SetEmModel(new G4DummyModel(), 1);
0225 pm->AddDiscreteProcess(eBrem);
0226
0227 const auto steplimit = new G4StepLimiter();
0228 pm->AddDiscreteProcess(steplimit);
0229 }
0230 if (fphysname == "Penelope") {
0231 const auto msc = new G4eMultipleScattering();
0232 msc->SetEmModel(new G4DummyModel(), 1);
0233 pm->AddDiscreteProcess(msc);
0234 const auto eIoni = new G4eIonisation("e-_G4PenelopeIoni");
0235 eIoni->SetEmModel(new G4DummyModel(), 1);
0236 pm->AddDiscreteProcess(eIoni);
0237 const auto eBrem = new G4eBremsstrahlung("e-_G4PenelopeBrem");
0238 eBrem->SetEmModel(new G4DummyModel(), 1);
0239 pm->AddDiscreteProcess(eBrem);
0240 const auto steplimit = new G4StepLimiter();
0241 pm->AddDiscreteProcess(steplimit);
0242 }
0243 }
0244
0245
0246
0247
0248 } else if (particleName == "proton") {
0249 } else if (particleName == "hydrogen") {
0250 } else if (particleName == "gamma") {
0251 G4PhotoElectricEffect *thePhotoElectricEffect = new G4PhotoElectricEffect();
0252 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel(), 1);
0253 pm->AddDiscreteProcess(thePhotoElectricEffect);
0254
0255 G4ComptonScattering *theComptonScattering = new G4ComptonScattering();
0256 theComptonScattering->SetEmModel(new G4LivermoreComptonModel(), 1);
0257 pm->AddDiscreteProcess(theComptonScattering);
0258
0259 G4GammaConversion *theGammaConversion = new G4GammaConversion();
0260 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel(), 1);
0261 pm->AddDiscreteProcess(theGammaConversion);
0262
0263 G4RayleighScattering *theRayleigh = new G4RayleighScattering();
0264 pm->AddDiscreteProcess(theRayleigh);
0265 }
0266 }
0267
0268 if (isGNP) {
0269 G4EmConfigurator *em_config = G4LossTableManager::Instance()->EmConfigurator();
0270 G4VEmModel *mod;
0271 mod = new G4DNAChampionElasticModel();
0272 mod->SetActivationLowEnergyLimit(1 * GeV);
0273 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic", mod, "NP");
0274 mod = new G4DNABornExcitationModel();
0275 mod->SetActivationLowEnergyLimit(1 * GeV);
0276 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation", mod, "NP");
0277 mod = new G4DNABornIonisationModel();
0278 mod->SetActivationLowEnergyLimit(1 * GeV);
0279 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation", mod, "NP");
0280 mod = new G4DNAMeltonAttachmentModel();
0281 mod->SetActivationLowEnergyLimit(1 * GeV);
0282 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment", mod, "NP");
0283 mod = new G4DNASancheExcitationModel();
0284 mod->SetActivationLowEnergyLimit(1 * GeV);
0285 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation", mod, "NP");
0286 mod = new G4SeltzerBergerModel();
0287 mod->SetActivationLowEnergyLimit(1 * GeV);
0288 em_config->SetExtraEmModel("e-", "e-_G4DNABremsstrahlung", mod, "NP");
0289
0290 if (fphysname == "DNA") {
0291 mod = new G4DNAELSEPAElasticModel();
0292 em_config->SetExtraEmModel("e-", "e-_G4DNAELSEPAElastic", mod, "NP", 10 * eV, 1 * GeV);
0293 mod = new G4DNADiracRMatrixExcitationModel();
0294 em_config->SetExtraEmModel("e-", "e-_G4DNADRMExcitation", mod, "NP", 10 * eV, 1 * GeV);
0295 mod = new G4DNARelativisticIonisationModel();
0296 em_config->SetExtraEmModel("e-", "e-_G4DNARelativisticIonisation", mod, "NP", 10 * eV,
0297 1 * GeV);
0298 mod = new G4DNAQuinnPlasmonExcitationModel();
0299 em_config->SetExtraEmModel("e-", "e-_G4DNAPlasmonExcitation", mod, "NP", 10 * eV, 1 * GeV);
0300 mod = new G4SeltzerBergerModel();
0301 em_config->SetExtraEmModel("e-", "e-_G4DNABremsstrahlung_GNP", mod, "NP", 10 * eV, 1 * GeV);
0302 }
0303 if (fphysname == "Livermore") {
0304 mod = new G4UrbanMscModel();
0305
0306 em_config->SetExtraEmModel("e-", "msc", mod, "NP", 0 * eV, 100 * MeV,
0307 new G4UniversalFluctuation());
0308 mod = new G4LivermoreIonisationModel();
0309 em_config->SetExtraEmModel("e-", "e-_G4LivermoreIoni", mod, "NP", 0 * eV, 1.0 * MeV,
0310 new G4UniversalFluctuation());
0311 mod = new G4LivermoreBremsstrahlungModel();
0312 em_config->SetExtraEmModel("e-", "e-_G4LivermoreBrem", mod, "NP", 0 * eV, 1 * GeV,
0313 new G4UniversalFluctuation());
0314 }
0315 if (fphysname == "Penelope") {
0316 mod = new G4UrbanMscModel();
0317
0318 em_config->SetExtraEmModel("e-", "msc", mod, "NP", 0 * eV, 100 * MeV,
0319 new G4UniversalFluctuation());
0320 mod = new G4PenelopeIonisationModel();
0321 em_config->SetExtraEmModel("e-", "e-_G4PenelopeIoni", mod, "NP", 0 * eV, 1.0 * MeV,
0322 new G4UniversalFluctuation());
0323 mod = new G4PenelopeBremsstrahlungModel();
0324 em_config->SetExtraEmModel("e-", "e-_G4PenelopeBrem", mod, "NP", 0 * eV, 1 * GeV,
0325 new G4UniversalFluctuation());
0326 }
0327 }
0328
0329 G4VAtomDeexcitation *de = new G4UAtomicDeexcitation();
0330 G4LossTableManager::Instance()->SetAtomDeexcitation(de);
0331 de->SetFluo(true);
0332 de->SetPIXE(true);
0333 de->SetAuger(true);
0334 de->SetAugerCascade(true);
0335 }
0336
0337
0338
0339 void PhysicsList::SetCuts() {
0340 if (verboseLevel > 0) {
0341 G4cout << "PhysicsList::SetCuts:";
0342 G4cout << "CutLength : " << G4BestUnit(defaultCutValue, "Length") << G4endl;
0343 }
0344
0345 SetCutValue(fcutForGamma, "gamma");
0346 SetCutValue(fcutForElectron, "e-");
0347 SetCutValue(fcutForPositron, "e+");
0348
0349 if (verboseLevel > 0) {
0350 DumpCutValuesTable();
0351 }
0352 }
0353
0354 void PhysicsList::SetPhysics4NP(const G4String &name) {
0355 fphysname = name;
0356 }