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 /// \file electromagnetic/TestEm1/src/PhysicsList.cc 27 /// \brief Implementation of the PhysicsList class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "PhysicsList.hh" 34 35 #include "DetectorConstruction.hh" 36 #include "PhysListEmStandard.hh" 37 #include "PhysicsListMessenger.hh" 38 39 #include "G4EmLivermorePhysics.hh" 40 #include "G4EmLowEPPhysics.hh" 41 #include "G4EmParameters.hh" 42 #include "G4EmPenelopePhysics.hh" 43 #include "G4EmStandardPhysics.hh" 44 #include "G4EmStandardPhysicsGS.hh" 45 #include "G4EmStandardPhysicsSS.hh" 46 #include "G4EmStandardPhysicsWVI.hh" 47 #include "G4EmStandardPhysics_option1.hh" 48 #include "G4EmStandardPhysics_option2.hh" 49 #include "G4EmStandardPhysics_option3.hh" 50 #include "G4EmStandardPhysics_option4.hh" 51 #include "G4LossTableManager.hh" 52 #include "G4SystemOfUnits.hh" 53 #include "G4UnitsTable.hh" 54 55 // particles 56 57 #include "StepMax.hh" 58 59 #include "G4BaryonConstructor.hh" 60 #include "G4BosonConstructor.hh" 61 #include "G4Decay.hh" 62 #include "G4Electron.hh" 63 #include "G4Gamma.hh" 64 #include "G4GenericIon.hh" 65 #include "G4IonConstructor.hh" 66 #include "G4LeptonConstructor.hh" 67 #include "G4Material.hh" 68 #include "G4MesonConstructor.hh" 69 #include "G4NuclideTable.hh" 70 #include "G4PhysicsListHelper.hh" 71 #include "G4ProcessManager.hh" 72 #include "G4Proton.hh" 73 #include "G4RadioactiveDecay.hh" 74 #include "G4ShortLivedConstructor.hh" 75 76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 77 78 PhysicsList::PhysicsList(DetectorConstruction* det) : fDet(det) 79 { 80 fMessenger = new PhysicsListMessenger(this); 81 SetVerboseLevel(1); 82 83 // EM physics 84 AddPhysicsList("emstandard_opt3"); 85 86 // fix lower limit for cut 87 G4ProductionCutsTable::GetProductionCutsTable()->SetEnergyRange(10 * eV, 1 * GeV); 88 SetDefaultCutValue(1 * mm); 89 } 90 91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 92 93 PhysicsList::~PhysicsList() 94 { 95 delete fMessenger; 96 } 97 98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 99 100 void PhysicsList::ConstructParticle() 101 { 102 G4BosonConstructor pBosonConstructor; 103 pBosonConstructor.ConstructParticle(); 104 105 G4LeptonConstructor pLeptonConstructor; 106 pLeptonConstructor.ConstructParticle(); 107 108 G4MesonConstructor pMesonConstructor; 109 pMesonConstructor.ConstructParticle(); 110 111 G4BaryonConstructor pBaryonConstructor; 112 pBaryonConstructor.ConstructParticle(); 113 114 G4IonConstructor pIonConstructor; 115 pIonConstructor.ConstructParticle(); 116 117 G4ShortLivedConstructor pShortLivedConstructor; 118 pShortLivedConstructor.ConstructParticle(); 119 } 120 121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 122 123 void PhysicsList::ConstructProcess() 124 { 125 // Transportation 126 // 127 AddTransportation(); 128 129 // Electromagnetic physics list 130 // 131 fEmPhysicsList->ConstructProcess(); 132 133 // Decay Process 134 // 135 AddDecay(); 136 137 // Decay Process 138 // 139 AddRadioactiveDecay(); 140 141 // step limitation (as a full process) 142 // 143 AddStepMax(); 144 145 // example of Get process 146 auto process = GetProcess("RadioactiveDecay"); 147 if (process != nullptr) { 148 G4cout << "\n GetProcess : " << process->GetProcessName() << G4endl; 149 } 150 } 151 152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 153 154 void PhysicsList::AddPhysicsList(const G4String& name) 155 { 156 if (verboseLevel > 0) { 157 G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl; 158 } 159 160 if (name == fEmName) return; 161 162 if (name == "local") { 163 fEmName = name; 164 delete fEmPhysicsList; 165 fEmPhysicsList = new PhysListEmStandard(name); 166 } 167 else if (name == "emstandard_opt0") { 168 fEmName = name; 169 delete fEmPhysicsList; 170 fEmPhysicsList = new G4EmStandardPhysics(); 171 } 172 else if (name == "emstandard_opt1") { 173 fEmName = name; 174 delete fEmPhysicsList; 175 fEmPhysicsList = new G4EmStandardPhysics_option1(); 176 } 177 else if (name == "emstandard_opt2") { 178 fEmName = name; 179 delete fEmPhysicsList; 180 fEmPhysicsList = new G4EmStandardPhysics_option2(); 181 } 182 else if (name == "emstandard_opt3") { 183 fEmName = name; 184 delete fEmPhysicsList; 185 fEmPhysicsList = new G4EmStandardPhysics_option3(); 186 } 187 else if (name == "emstandard_opt4") { 188 fEmName = name; 189 delete fEmPhysicsList; 190 fEmPhysicsList = new G4EmStandardPhysics_option4(); 191 } 192 else if (name == "emstandardSS") { 193 fEmName = name; 194 delete fEmPhysicsList; 195 fEmPhysicsList = new G4EmStandardPhysicsSS(); 196 } 197 else if (name == "emstandardGS") { 198 fEmName = name; 199 delete fEmPhysicsList; 200 fEmPhysicsList = new G4EmStandardPhysicsGS(); 201 } 202 else if (name == "emstandardWVI") { 203 fEmName = name; 204 delete fEmPhysicsList; 205 fEmPhysicsList = new G4EmStandardPhysicsWVI(); 206 } 207 else if (name == "emlivermore") { 208 fEmName = name; 209 delete fEmPhysicsList; 210 fEmPhysicsList = new G4EmLivermorePhysics(); 211 } 212 else if (name == "empenelope") { 213 fEmName = name; 214 delete fEmPhysicsList; 215 fEmPhysicsList = new G4EmPenelopePhysics(); 216 } 217 else if (name == "emlowenergy") { 218 fEmName = name; 219 delete fEmPhysicsList; 220 fEmPhysicsList = new G4EmLowEPPhysics(); 221 } 222 else { 223 G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" 224 << " is not defined" << G4endl; 225 } 226 227 // Em options 228 // 229 G4EmParameters::Instance()->SetBuildCSDARange(true); 230 G4EmParameters::Instance()->SetGeneralProcessActive(false); 231 } 232 233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 234 235 void PhysicsList::AddDecay() 236 { 237 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 238 239 // Decay Process 240 // 241 G4Decay* fDecayProcess = new G4Decay(); 242 243 auto particleIterator = GetParticleIterator(); 244 particleIterator->reset(); 245 while ((*particleIterator)()) { 246 G4ParticleDefinition* particle = particleIterator->value(); 247 if (fDecayProcess->IsApplicable(*particle) && !particle->IsShortLived()) 248 ph->RegisterProcess(fDecayProcess, particle); 249 } 250 } 251 252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 253 254 void PhysicsList::AddRadioactiveDecay() 255 { 256 G4RadioactiveDecay* radioactiveDecay = new G4RadioactiveDecay(); 257 258 G4bool armFlag = false; 259 radioactiveDecay->SetARM(armFlag); // Atomic Rearangement 260 261 // atomic de-excitation module 262 if (armFlag) { 263 G4EmParameters::Instance()->SetAuger(true); 264 G4EmParameters::Instance()->SetDeexcitationIgnoreCut(true); 265 } 266 267 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 268 ph->RegisterProcess(radioactiveDecay, G4GenericIon::GenericIon()); 269 270 // mandatory for G4NuclideTable 271 // 272 const G4double meanLife = 1 * picosecond, halfLife = meanLife * std::log(2); 273 G4NuclideTable::GetInstance()->SetThresholdOfHalfLife(halfLife); 274 } 275 276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 277 278 void PhysicsList::AddStepMax() 279 { 280 // Step limitation seen as a process 281 StepMax* stepMaxProcess = new StepMax(); 282 283 auto particleIterator = GetParticleIterator(); 284 particleIterator->reset(); 285 while ((*particleIterator)()) { 286 G4ParticleDefinition* particle = particleIterator->value(); 287 G4ProcessManager* pmanager = particle->GetProcessManager(); 288 289 if (stepMaxProcess->IsApplicable(*particle) && !particle->IsShortLived()) 290 pmanager->AddDiscreteProcess(stepMaxProcess); 291 } 292 } 293 294 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 295 296 void PhysicsList::GetRange(G4double val) 297 { 298 G4LogicalVolume* lBox = fDet->GetWorld()->GetLogicalVolume(); 299 const G4MaterialCutsCouple* couple = lBox->GetMaterialCutsCouple(); 300 const G4Material* currMat = lBox->GetMaterial(); 301 302 G4ParticleDefinition* part; 303 G4double cut; 304 part = G4Electron::Electron(); 305 cut = G4LossTableManager::Instance()->GetRange(part, val, couple); 306 G4cout << "material : " << currMat->GetName() << G4endl; 307 G4cout << "particle : " << part->GetParticleName() << G4endl; 308 G4cout << "energy : " << G4BestUnit(val, "Energy") << G4endl; 309 G4cout << "range : " << G4BestUnit(cut, "Length") << G4endl; 310 } 311 312 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 313 314 G4VProcess* PhysicsList::GetProcess(const G4String& processName) const 315 { 316 G4ParticleDefinition* particle = G4GenericIon::GenericIon(); 317 G4ProcessVector* procList = particle->GetProcessManager()->GetProcessList(); 318 G4int nbProc = particle->GetProcessManager()->GetProcessListLength(); 319 for (G4int k = 0; k < nbProc; k++) { 320 G4VProcess* process = (*procList)[k]; 321 if (process->GetProcessName() == processName) return process; 322 } 323 return nullptr; 324 } 325 326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 327