Geant4 Cross Reference |
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 // Authors: Francesco Longo, franzlongo1969@gmail.com 27 // 28 // Code based on the hadrontherapy && radioprotection advanced example 29 30 #include "GammaRayTelPhysicsList.hh" 31 #include "GammaRayTelPhysicsListMessenger.hh" 32 #include "G4PhysListFactory.hh" 33 #include "G4VPhysicsConstructor.hh" 34 35 // Physics lists (contained inside the Geant4 distribution) 36 #include "G4EmLivermorePhysics.hh" 37 #include "G4EmLivermorePolarizedPhysics.hh" 38 #include "G4EmPenelopePhysics.hh" 39 #include "G4EmStandardPhysics_option3.hh" 40 #include "G4EmStandardPhysics_option4.hh" // to treat the new polarized process 41 42 #include "G4Decay.hh" 43 #include "G4DecayPhysics.hh" 44 #include "G4HadronDElasticPhysics.hh" 45 #include "G4HadronElasticPhysics.hh" 46 #include "G4HadronElasticPhysicsHP.hh" 47 #include "G4HadronPhysicsQGSP_BIC_HP.hh" 48 #include "G4IonBinaryCascadePhysics.hh" 49 #include "G4IonFluctuations.hh" 50 #include "G4IonParametrisedLossModel.hh" 51 #include "G4LossTableManager.hh" 52 #include "G4ProcessManager.hh" 53 #include "G4RadioactiveDecayPhysics.hh" 54 #include "G4StepLimiter.hh" 55 #include "G4SystemOfUnits.hh" 56 #include "G4UnitsTable.hh" 57 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 59 60 GammaRayTelPhysicsList::GammaRayTelPhysicsList() { 61 G4LossTableManager::Instance(); 62 63 constexpr auto DEFAULT_CUT_VALUE{100 * micrometer}; 64 defaultCutValue = DEFAULT_CUT_VALUE; 65 66 SetVerboseLevel(1); 67 68 constexpr auto ENERGY_LOWER_BOUND{250 * eV}; 69 constexpr auto ENERGY_UPPER_BOUND{1 * GeV}; 70 71 G4ProductionCutsTable::GetProductionCutsTable()->SetEnergyRange(ENERGY_LOWER_BOUND, ENERGY_UPPER_BOUND); 72 SetDefaultCutValue (defaultCutValue); 73 DumpCutValuesTable(); 74 75 helIsRegisted = false; 76 bicIsRegisted = false; 77 biciIsRegisted = false; 78 locIonIonInelasticIsRegistered = false; 79 radioactiveDecayIsRegisted = false; 80 81 pMessenger = new GammaRayTelPhysicsListMessenger(this); 82 83 SetVerboseLevel(1); 84 85 // EM physics 86 emPhysicsList = new G4EmStandardPhysics_option3(1); 87 emName = G4String("emstandard_opt3"); 88 89 // Decay physics and all particles 90 decPhysicsList = new G4DecayPhysics(); 91 } 92 93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 94 95 GammaRayTelPhysicsList::~GammaRayTelPhysicsList() { 96 delete pMessenger; 97 delete emPhysicsList; 98 delete decPhysicsList; 99 100 for (auto &hadronPhy : hadronPhys) { 101 delete hadronPhy; 102 } 103 } 104 105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 106 107 void GammaRayTelPhysicsList::AddPackage(const G4String &name) { 108 G4PhysListFactory factory; 109 auto *phys = factory.GetReferencePhysList(name); 110 111 G4int i = 0; 112 const G4VPhysicsConstructor *element = phys->GetPhysics(i); 113 auto *tmp = const_cast<G4VPhysicsConstructor*>(element); 114 115 while (element != nullptr) { 116 RegisterPhysics(tmp); 117 element = phys->GetPhysics(++i); 118 tmp = const_cast<G4VPhysicsConstructor*>(element); 119 } 120 } 121 122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 123 124 void GammaRayTelPhysicsList::ConstructParticle() { 125 decPhysicsList->ConstructParticle(); 126 } 127 128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 129 130 void GammaRayTelPhysicsList::ConstructProcess() { 131 // transportation 132 // 133 AddTransportation(); 134 135 // electromagnetic physics list 136 // 137 emPhysicsList->ConstructProcess(); 138 emConfigurator.AddModels(); 139 140 // decay physics list 141 // 142 decPhysicsList->ConstructProcess(); 143 144 // hadronic physics lists 145 for (auto &hadronPhy : hadronPhys) { 146 hadronPhy->ConstructProcess(); 147 } 148 149 // step limitation (as a full process) 150 // AddStepMax(); 151 } 152 153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 154 155 void GammaRayTelPhysicsList::AddPhysicsList(const G4String &name) { 156 if (verboseLevel > 1) { 157 G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl; 158 } 159 if (name == emName) { 160 return; 161 } 162 163 //////////////////////////////// 164 // ELECTROMAGNETIC MODELS // 165 //////////////////////////////// 166 167 if (name == "standard_opt3") { 168 emName = name; 169 delete emPhysicsList; 170 emPhysicsList = new G4EmStandardPhysics_option3(); 171 G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmStandardPhysics_option3" << G4endl; 172 } else if (name == "LowE_Livermore") { 173 emName = name; 174 delete emPhysicsList; 175 emPhysicsList = new G4EmLivermorePhysics(); 176 G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl; 177 } else if (name == "LowE_Penelope") { 178 emName = name; 179 delete emPhysicsList; 180 emPhysicsList = new G4EmPenelopePhysics(); 181 G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl; 182 } else if (name == "LowE_Polarized") { 183 emName = name; 184 delete emPhysicsList; 185 emPhysicsList = new G4EmLivermorePolarizedPhysics(); 186 G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmLivermorePhysics" << G4endl; 187 } else if (name == "standard_opt4") { 188 emName = name; 189 delete emPhysicsList; 190 emPhysicsList = new G4EmStandardPhysics_option4(); 191 G4cout << "THE FOLLOWING ELECTROMAGNETIC PHYSICS LIST HAS BEEN ACTIVATED: G4EmStandardOption_4" << G4endl; 192 193 ///////////////////////// 194 // HADRONIC MODELS // 195 ///////////////////////// 196 197 } else if (name == "elastic" && !helIsRegisted) { 198 G4cout << "THE FOLLOWING HADRONIC ELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4HadronElasticPhysics" << G4endl; 199 hadronPhys.push_back(new G4HadronElasticPhysics()); 200 helIsRegisted = true; 201 } else if (name == "DElastic" && !helIsRegisted) { 202 hadronPhys.push_back(new G4HadronDElasticPhysics()); 203 helIsRegisted = true; 204 } else if (name == "HPElastic" && !helIsRegisted) { 205 hadronPhys.push_back(new G4HadronElasticPhysicsHP()); 206 helIsRegisted = true; 207 } else if (name == "binary" && !bicIsRegisted) { 208 hadronPhys.push_back(new G4HadronPhysicsQGSP_BIC_HP()); 209 bicIsRegisted = true; 210 G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: HadronPhysicsQGSP_BIC_HP" << G4endl; 211 } else if (name == "binary_ion" && !biciIsRegisted) { 212 hadronPhys.push_back(new G4IonBinaryCascadePhysics()); 213 biciIsRegisted = true; 214 G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4IonBinaryCascadePhysics" << G4endl; 215 } else if (name == "radioactive_decay" && !radioactiveDecayIsRegisted) { 216 hadronPhys.push_back(new G4RadioactiveDecayPhysics()); 217 radioactiveDecayIsRegisted = true; 218 G4cout << "THE FOLLOWING HADRONIC INELASTIC PHYSICS LIST HAS BEEN ACTIVATED: G4RadioactiveDecayPhysics" << G4endl; 219 } else { 220 G4cout << "PhysicsList::AddPhysicsList: <" << name << "> is not defined" << G4endl; 221 } 222 } 223 224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 225 226 void GammaRayTelPhysicsList::SetCutForGamma(G4double cut) { 227 SetParticleCuts(cut, G4Gamma::Gamma()); 228 } 229 230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 231 232 void GammaRayTelPhysicsList::SetCutForElectron(G4double cut) { 233 SetParticleCuts(cut, G4Electron::Electron()); 234 } 235 236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 237 238 void GammaRayTelPhysicsList::SetCutForPositron(G4double cut) { 239 SetParticleCuts(cut, G4Positron::Positron()); 240 } 241