[1203] | 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 | // |
---|
[1315] | 26 | // $Id: G4EmPenelopePhysics.cc,v 1.9 2010/06/02 17:21:29 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ |
---|
[1203] | 28 | |
---|
| 29 | #include "G4EmPenelopePhysics.hh" |
---|
| 30 | |
---|
| 31 | #include "G4ParticleDefinition.hh" |
---|
| 32 | #include "G4ProcessManager.hh" |
---|
| 33 | |
---|
| 34 | // *** Processes and models |
---|
| 35 | |
---|
| 36 | // gamma |
---|
| 37 | |
---|
| 38 | #include "G4PhotoElectricEffect.hh" |
---|
| 39 | #include "G4PenelopePhotoElectricModel.hh" |
---|
| 40 | |
---|
| 41 | #include "G4ComptonScattering.hh" |
---|
| 42 | #include "G4PenelopeComptonModel.hh" |
---|
| 43 | |
---|
| 44 | #include "G4GammaConversion.hh" |
---|
| 45 | #include "G4PenelopeGammaConversionModel.hh" |
---|
| 46 | |
---|
| 47 | #include "G4RayleighScattering.hh" |
---|
| 48 | #include "G4PenelopeRayleighModel.hh" |
---|
| 49 | |
---|
| 50 | // e- and e+ |
---|
| 51 | |
---|
| 52 | #include "G4eMultipleScattering.hh" |
---|
| 53 | #include "G4UniversalFluctuation.hh" |
---|
| 54 | |
---|
| 55 | #include "G4eIonisation.hh" |
---|
| 56 | #include "G4PenelopeIonisationModel.hh" |
---|
| 57 | |
---|
| 58 | #include "G4eBremsstrahlung.hh" |
---|
| 59 | #include "G4PenelopeBremsstrahlungModel.hh" |
---|
| 60 | |
---|
| 61 | // e+ only |
---|
| 62 | |
---|
| 63 | #include "G4eplusAnnihilation.hh" |
---|
| 64 | #include "G4PenelopeAnnihilationModel.hh" |
---|
| 65 | |
---|
| 66 | // mu |
---|
| 67 | |
---|
| 68 | #include "G4MuMultipleScattering.hh" |
---|
| 69 | #include "G4MuIonisation.hh" |
---|
| 70 | #include "G4MuBremsstrahlung.hh" |
---|
| 71 | #include "G4MuPairProduction.hh" |
---|
| 72 | |
---|
| 73 | // hadrons |
---|
| 74 | |
---|
| 75 | #include "G4hMultipleScattering.hh" |
---|
| 76 | #include "G4MscStepLimitType.hh" |
---|
| 77 | |
---|
| 78 | #include "G4hBremsstrahlung.hh" |
---|
| 79 | #include "G4hPairProduction.hh" |
---|
| 80 | |
---|
| 81 | #include "G4hIonisation.hh" |
---|
| 82 | #include "G4ionIonisation.hh" |
---|
| 83 | #include "G4IonParametrisedLossModel.hh" |
---|
| 84 | #include "G4NuclearStopping.hh" |
---|
| 85 | |
---|
| 86 | // msc models |
---|
| 87 | #include "G4UrbanMscModel93.hh" |
---|
| 88 | #include "G4WentzelVIModel.hh" |
---|
[1315] | 89 | #include "G4GoudsmitSaundersonMscModel.hh" |
---|
[1203] | 90 | #include "G4CoulombScattering.hh" |
---|
| 91 | |
---|
| 92 | // |
---|
| 93 | |
---|
| 94 | #include "G4LossTableManager.hh" |
---|
| 95 | #include "G4EmProcessOptions.hh" |
---|
| 96 | |
---|
| 97 | // particles |
---|
| 98 | |
---|
| 99 | #include "G4Gamma.hh" |
---|
| 100 | #include "G4Electron.hh" |
---|
| 101 | #include "G4Positron.hh" |
---|
| 102 | #include "G4MuonPlus.hh" |
---|
| 103 | #include "G4MuonMinus.hh" |
---|
| 104 | #include "G4PionPlus.hh" |
---|
| 105 | #include "G4PionMinus.hh" |
---|
| 106 | #include "G4KaonPlus.hh" |
---|
| 107 | #include "G4KaonMinus.hh" |
---|
| 108 | #include "G4Proton.hh" |
---|
| 109 | #include "G4AntiProton.hh" |
---|
| 110 | #include "G4Deuteron.hh" |
---|
| 111 | #include "G4Triton.hh" |
---|
| 112 | #include "G4He3.hh" |
---|
| 113 | #include "G4Alpha.hh" |
---|
| 114 | #include "G4GenericIon.hh" |
---|
| 115 | |
---|
| 116 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 117 | |
---|
[1315] | 118 | G4EmPenelopePhysics::G4EmPenelopePhysics(G4int ver) |
---|
| 119 | : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver) |
---|
[1203] | 120 | { |
---|
| 121 | G4LossTableManager::Instance(); |
---|
| 122 | } |
---|
| 123 | |
---|
| 124 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 125 | |
---|
[1315] | 126 | G4EmPenelopePhysics::G4EmPenelopePhysics(G4int ver, const G4String&) |
---|
| 127 | : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver) |
---|
| 128 | { |
---|
| 129 | G4LossTableManager::Instance(); |
---|
| 130 | } |
---|
| 131 | |
---|
| 132 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 133 | |
---|
[1203] | 134 | G4EmPenelopePhysics::~G4EmPenelopePhysics() |
---|
| 135 | {} |
---|
| 136 | |
---|
| 137 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 138 | |
---|
| 139 | void G4EmPenelopePhysics::ConstructParticle() |
---|
| 140 | { |
---|
| 141 | // gamma |
---|
| 142 | G4Gamma::Gamma(); |
---|
| 143 | |
---|
| 144 | // leptons |
---|
| 145 | G4Electron::Electron(); |
---|
| 146 | G4Positron::Positron(); |
---|
| 147 | G4MuonPlus::MuonPlus(); |
---|
| 148 | G4MuonMinus::MuonMinus(); |
---|
| 149 | |
---|
| 150 | // mesons |
---|
| 151 | G4PionPlus::PionPlusDefinition(); |
---|
| 152 | G4PionMinus::PionMinusDefinition(); |
---|
| 153 | G4KaonPlus::KaonPlusDefinition(); |
---|
| 154 | G4KaonMinus::KaonMinusDefinition(); |
---|
| 155 | |
---|
| 156 | // baryons |
---|
| 157 | G4Proton::Proton(); |
---|
| 158 | G4AntiProton::AntiProton(); |
---|
| 159 | |
---|
| 160 | // ions |
---|
| 161 | G4Deuteron::Deuteron(); |
---|
| 162 | G4Triton::Triton(); |
---|
| 163 | G4He3::He3(); |
---|
| 164 | G4Alpha::Alpha(); |
---|
| 165 | G4GenericIon::GenericIonDefinition(); |
---|
| 166 | } |
---|
| 167 | |
---|
| 168 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 169 | |
---|
| 170 | void G4EmPenelopePhysics::ConstructProcess() |
---|
| 171 | { |
---|
| 172 | // Add Penelope EM Processes |
---|
| 173 | |
---|
| 174 | theParticleIterator->reset(); |
---|
| 175 | |
---|
| 176 | while( (*theParticleIterator)() ){ |
---|
| 177 | |
---|
| 178 | G4ParticleDefinition* particle = theParticleIterator->value(); |
---|
| 179 | G4ProcessManager* pmanager = particle->GetProcessManager(); |
---|
| 180 | G4String particleName = particle->GetParticleName(); |
---|
| 181 | |
---|
| 182 | if(verbose > 1) |
---|
| 183 | G4cout << "### " << GetPhysicsName() << " instantiates for " |
---|
| 184 | << particleName << G4endl; |
---|
| 185 | |
---|
| 186 | //Applicability range for Penelope models |
---|
| 187 | //for higher energies, the Standard models are used |
---|
| 188 | G4double PenelopeHighEnergyLimit = 1.0*GeV; |
---|
| 189 | |
---|
| 190 | if (particleName == "gamma") { |
---|
| 191 | |
---|
| 192 | //Photo-electric effect |
---|
| 193 | G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect(); |
---|
| 194 | G4PenelopePhotoElectricModel* thePEPenelopeModel = new |
---|
| 195 | G4PenelopePhotoElectricModel(); |
---|
| 196 | thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 197 | thePhotoElectricEffect->AddEmModel(0,thePEPenelopeModel); |
---|
| 198 | pmanager->AddDiscreteProcess(thePhotoElectricEffect); |
---|
| 199 | |
---|
| 200 | //Compton scattering |
---|
| 201 | G4ComptonScattering* theComptonScattering = new G4ComptonScattering(); |
---|
| 202 | G4PenelopeComptonModel* theComptonPenelopeModel = |
---|
| 203 | new G4PenelopeComptonModel(); |
---|
| 204 | theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 205 | theComptonScattering->AddEmModel(0,theComptonPenelopeModel); |
---|
| 206 | pmanager->AddDiscreteProcess(theComptonScattering); |
---|
| 207 | |
---|
| 208 | //Gamma conversion |
---|
| 209 | G4GammaConversion* theGammaConversion = new G4GammaConversion(); |
---|
| 210 | G4PenelopeGammaConversionModel* theGCPenelopeModel = |
---|
| 211 | new G4PenelopeGammaConversionModel(); |
---|
| 212 | theGammaConversion->AddEmModel(0,theGCPenelopeModel); |
---|
| 213 | pmanager->AddDiscreteProcess(theGammaConversion); |
---|
| 214 | |
---|
| 215 | //Rayleigh scattering |
---|
| 216 | G4RayleighScattering* theRayleigh = new G4RayleighScattering(); |
---|
| 217 | G4PenelopeRayleighModel* theRayleighPenelopeModel = |
---|
| 218 | new G4PenelopeRayleighModel(); |
---|
| 219 | theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 220 | theRayleigh->AddEmModel(0,theRayleighPenelopeModel); |
---|
| 221 | pmanager->AddDiscreteProcess(theRayleigh); |
---|
| 222 | |
---|
| 223 | } else if (particleName == "e-") { |
---|
| 224 | |
---|
| 225 | G4eMultipleScattering* msc = new G4eMultipleScattering(); |
---|
[1315] | 226 | //msc->AddEmModel(0, new G4UrbanMscModel93()); |
---|
| 227 | msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel()); |
---|
[1203] | 228 | msc->SetStepLimitType(fUseDistanceToBoundary); |
---|
| 229 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
| 230 | |
---|
| 231 | //Ionisation |
---|
| 232 | G4eIonisation* eIoni = new G4eIonisation(); |
---|
| 233 | G4PenelopeIonisationModel* theIoniPenelope = |
---|
| 234 | new G4PenelopeIonisationModel(); |
---|
| 235 | theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 236 | eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation()); |
---|
| 237 | eIoni->SetStepFunction(0.2, 100*um); // |
---|
| 238 | pmanager->AddProcess(eIoni, -1, 2, 2); |
---|
| 239 | |
---|
| 240 | //Bremsstrahlung |
---|
| 241 | G4eBremsstrahlung* eBrem = new G4eBremsstrahlung(); |
---|
| 242 | G4PenelopeBremsstrahlungModel* theBremPenelope = new |
---|
| 243 | G4PenelopeBremsstrahlungModel(); |
---|
| 244 | theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 245 | eBrem->AddEmModel(0,theBremPenelope); |
---|
| 246 | pmanager->AddProcess(eBrem, -1,-3, 3); |
---|
| 247 | |
---|
| 248 | } else if (particleName == "e+") { |
---|
| 249 | |
---|
| 250 | G4eMultipleScattering* msc = new G4eMultipleScattering(); |
---|
[1315] | 251 | //msc->AddEmModel(0, new G4UrbanMscModel93()); |
---|
| 252 | msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel()); |
---|
[1203] | 253 | msc->SetStepLimitType(fUseDistanceToBoundary); |
---|
| 254 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
| 255 | |
---|
| 256 | //Ionisation |
---|
| 257 | G4eIonisation* eIoni = new G4eIonisation(); |
---|
| 258 | G4PenelopeIonisationModel* theIoniPenelope = |
---|
| 259 | new G4PenelopeIonisationModel(); |
---|
| 260 | theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 261 | eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation()); |
---|
| 262 | eIoni->SetStepFunction(0.2, 100*um); // |
---|
| 263 | pmanager->AddProcess(eIoni, -1, 2, 2); |
---|
| 264 | |
---|
| 265 | //Bremsstrahlung |
---|
| 266 | G4eBremsstrahlung* eBrem = new G4eBremsstrahlung(); |
---|
| 267 | G4PenelopeBremsstrahlungModel* theBremPenelope = new |
---|
| 268 | G4PenelopeBremsstrahlungModel(); |
---|
| 269 | theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 270 | eBrem->AddEmModel(0,theBremPenelope); |
---|
| 271 | pmanager->AddProcess(eBrem, -1,-3, 3); |
---|
| 272 | |
---|
| 273 | //Annihilation |
---|
| 274 | G4eplusAnnihilation* eAnni = new G4eplusAnnihilation(); |
---|
| 275 | G4PenelopeAnnihilationModel* theAnnPenelope = new |
---|
| 276 | G4PenelopeAnnihilationModel(); |
---|
| 277 | theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit); |
---|
| 278 | eAnni->AddEmModel(0,theAnnPenelope); |
---|
| 279 | pmanager->AddProcess(eAnni,0,-1, 4); |
---|
| 280 | |
---|
| 281 | } else if (particleName == "mu+" || |
---|
| 282 | particleName == "mu-" ) { |
---|
| 283 | |
---|
| 284 | // Identical to G4EmStandardPhysics_option3 |
---|
| 285 | |
---|
| 286 | G4MuMultipleScattering* msc = new G4MuMultipleScattering(); |
---|
| 287 | msc->AddEmModel(0, new G4WentzelVIModel()); |
---|
| 288 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
| 289 | |
---|
| 290 | G4MuIonisation* muIoni = new G4MuIonisation(); |
---|
| 291 | muIoni->SetStepFunction(0.2, 50*um); |
---|
| 292 | |
---|
| 293 | pmanager->AddProcess(muIoni, -1, 2, 2); |
---|
| 294 | pmanager->AddProcess(new G4MuBremsstrahlung, -1,-3, 3); |
---|
| 295 | pmanager->AddProcess(new G4MuPairProduction, -1,-4, 4); |
---|
| 296 | pmanager->AddDiscreteProcess(new G4CoulombScattering()); |
---|
| 297 | |
---|
| 298 | } else if (particleName == "GenericIon") { |
---|
| 299 | |
---|
| 300 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
| 301 | |
---|
| 302 | G4ionIonisation* ionIoni = new G4ionIonisation(); |
---|
| 303 | ionIoni->SetEmModel(new G4IonParametrisedLossModel()); |
---|
| 304 | ionIoni->SetStepFunction(0.1, 10*um); |
---|
| 305 | pmanager->AddProcess(ionIoni, -1, 2, 2); |
---|
| 306 | pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1); |
---|
| 307 | |
---|
| 308 | } else if (particleName == "alpha" || |
---|
| 309 | particleName == "He3" ) { |
---|
| 310 | |
---|
| 311 | // Identical to G4EmStandardPhysics_option3 |
---|
| 312 | |
---|
| 313 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
| 314 | |
---|
| 315 | G4ionIonisation* ionIoni = new G4ionIonisation(); |
---|
| 316 | ionIoni->SetStepFunction(0.1, 20*um); |
---|
| 317 | pmanager->AddProcess(ionIoni, -1, 2, 2); |
---|
| 318 | pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1); |
---|
| 319 | |
---|
| 320 | } else if (particleName == "pi+" || |
---|
| 321 | particleName == "pi-" || |
---|
| 322 | particleName == "kaon+" || |
---|
| 323 | particleName == "kaon-" || |
---|
| 324 | particleName == "proton" ) { |
---|
| 325 | |
---|
| 326 | // Identical to G4EmStandardPhysics_option3 |
---|
| 327 | |
---|
| 328 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
| 329 | |
---|
| 330 | G4hIonisation* hIoni = new G4hIonisation(); |
---|
| 331 | hIoni->SetStepFunction(0.2, 50*um); |
---|
| 332 | |
---|
| 333 | pmanager->AddProcess(hIoni, -1, 2, 2); |
---|
| 334 | pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3); |
---|
| 335 | pmanager->AddProcess(new G4hPairProduction, -1,-4, 4); |
---|
| 336 | |
---|
| 337 | } else if (particleName == "B+" || |
---|
| 338 | particleName == "B-" || |
---|
| 339 | particleName == "D+" || |
---|
| 340 | particleName == "D-" || |
---|
| 341 | particleName == "Ds+" || |
---|
| 342 | particleName == "Ds-" || |
---|
| 343 | particleName == "anti_lambda_c+" || |
---|
| 344 | particleName == "anti_omega-" || |
---|
| 345 | particleName == "anti_proton" || |
---|
| 346 | particleName == "anti_sigma_c+" || |
---|
| 347 | particleName == "anti_sigma_c++" || |
---|
| 348 | particleName == "anti_sigma+" || |
---|
| 349 | particleName == "anti_sigma-" || |
---|
| 350 | particleName == "anti_xi_c+" || |
---|
| 351 | particleName == "anti_xi-" || |
---|
| 352 | particleName == "deuteron" || |
---|
| 353 | particleName == "lambda_c+" || |
---|
| 354 | particleName == "omega-" || |
---|
| 355 | particleName == "sigma_c+" || |
---|
| 356 | particleName == "sigma_c++" || |
---|
| 357 | particleName == "sigma+" || |
---|
| 358 | particleName == "sigma-" || |
---|
| 359 | particleName == "tau+" || |
---|
| 360 | particleName == "tau-" || |
---|
| 361 | particleName == "triton" || |
---|
| 362 | particleName == "xi_c+" || |
---|
| 363 | particleName == "xi-" ) { |
---|
| 364 | |
---|
| 365 | // Identical to G4EmStandardPhysics_option3 |
---|
| 366 | |
---|
| 367 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
| 368 | pmanager->AddProcess(new G4hIonisation, -1, 2, 2); |
---|
| 369 | |
---|
| 370 | } |
---|
| 371 | } |
---|
| 372 | |
---|
| 373 | // Em options |
---|
| 374 | // |
---|
| 375 | G4EmProcessOptions opt; |
---|
| 376 | opt.SetVerbose(verbose); |
---|
| 377 | |
---|
| 378 | // Multiple Coulomb scattering |
---|
| 379 | // |
---|
| 380 | //opt.SetMscStepLimitation(fUseDistanceToBoundary); |
---|
| 381 | //opt.SetMscRangeFactor(0.02); |
---|
| 382 | |
---|
| 383 | // Physics tables |
---|
| 384 | // |
---|
| 385 | |
---|
| 386 | opt.SetMinEnergy(100*eV); |
---|
| 387 | opt.SetMaxEnergy(10*TeV); |
---|
| 388 | opt.SetDEDXBinning(220); |
---|
| 389 | opt.SetLambdaBinning(220); |
---|
| 390 | |
---|
| 391 | //opt.SetSplineFlag(true); |
---|
| 392 | opt.SetPolarAngleLimit(0.2); |
---|
| 393 | |
---|
| 394 | // Ionization |
---|
| 395 | // |
---|
| 396 | //opt.SetSubCutoff(true); |
---|
| 397 | } |
---|
| 398 | |
---|
| 399 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|