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/TestEm7/src/PhysListEmStandardNR.cc 27 /// \brief Implementation of the PhysListEmStandardNR class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "PhysListEmStandardNR.hh" 34 35 #include "G4BuilderType.hh" 36 #include "G4ComptonScattering.hh" 37 #include "G4CoulombScattering.hh" 38 #include "G4EmParameters.hh" 39 #include "G4GammaConversion.hh" 40 #include "G4Generator2BS.hh" 41 #include "G4IonParametrisedLossModel.hh" 42 #include "G4KleinNishinaModel.hh" 43 #include "G4LivermorePhotoElectricModel.hh" 44 #include "G4LossTableManager.hh" 45 #include "G4LowEPComptonModel.hh" 46 #include "G4MuBremsstrahlung.hh" 47 #include "G4MuIonisation.hh" 48 #include "G4MuMultipleScattering.hh" 49 #include "G4MuPairProduction.hh" 50 #include "G4PEEffectFluoModel.hh" 51 #include "G4ParticleDefinition.hh" 52 #include "G4PenelopeGammaConversionModel.hh" 53 #include "G4PenelopeIonisationModel.hh" 54 #include "G4PhotoElectricEffect.hh" 55 #include "G4PhysicsListHelper.hh" 56 #include "G4RayleighScattering.hh" 57 #include "G4ScreenedNuclearRecoil.hh" 58 #include "G4SeltzerBergerModel.hh" 59 #include "G4SystemOfUnits.hh" 60 #include "G4UAtomicDeexcitation.hh" 61 #include "G4UniversalFluctuation.hh" 62 #include "G4UrbanMscModel.hh" 63 #include "G4eBremsstrahlung.hh" 64 #include "G4eCoulombScatteringModel.hh" 65 #include "G4eIonisation.hh" 66 #include "G4eMultipleScattering.hh" 67 #include "G4eplusAnnihilation.hh" 68 #include "G4hIonisation.hh" 69 #include "G4hMultipleScattering.hh" 70 #include "G4ionIonisation.hh" 71 72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 73 74 PhysListEmStandardNR::PhysListEmStandardNR(const G4String& name) : G4VPhysicsConstructor(name) 75 { 76 G4EmParameters* param = G4EmParameters::Instance(); 77 param->SetDefaults(); 78 param->SetStepFunction(0.2, 100 * um); 79 param->SetStepFunctionMuHad(0.2, 50 * um); 80 param->SetStepFunctionLightIons(0.1, 10 * um); 81 param->SetStepFunctionIons(0.1, 1 * um); 82 param->SetFluo(true); 83 SetPhysicsType(bElectromagnetic); 84 } 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 88 PhysListEmStandardNR::~PhysListEmStandardNR() {} 89 90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 91 92 void PhysListEmStandardNR::ConstructProcess() 93 { 94 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 95 96 // muon & hadron bremsstrahlung and pair production 97 G4MuBremsstrahlung* mub = new G4MuBremsstrahlung(); 98 G4MuPairProduction* mup = new G4MuPairProduction(); 99 100 G4ScreenedNuclearRecoil* nucr = new G4ScreenedNuclearRecoil(); 101 G4double energyLimit = 100. * MeV; 102 nucr->SetMaxEnergyForScattering(energyLimit); 103 G4eCoulombScatteringModel* csm = new G4eCoulombScatteringModel(); 104 csm->SetActivationLowEnergyLimit(energyLimit); 105 106 auto particleIterator = GetParticleIterator(); 107 particleIterator->reset(); 108 while ((*particleIterator)()) { 109 G4ParticleDefinition* particle = particleIterator->value(); 110 G4String particleName = particle->GetParticleName(); 111 112 if (particleName == "gamma") { 113 // Compton scattering 114 G4ComptonScattering* cs = new G4ComptonScattering; 115 cs->SetEmModel(new G4KleinNishinaModel()); 116 ph->RegisterProcess(cs, particle); 117 118 // Photoelectric 119 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect(); 120 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel(); 121 theLivermorePEModel->SetHighEnergyLimit(10 * GeV); 122 pe->SetEmModel(theLivermorePEModel); 123 ph->RegisterProcess(pe, particle); 124 125 // Gamma conversion 126 G4GammaConversion* gc = new G4GammaConversion(); 127 G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel(); 128 thePenelopeGCModel->SetHighEnergyLimit(1 * GeV); 129 gc->SetEmModel(thePenelopeGCModel); 130 ph->RegisterProcess(gc, particle); 131 132 // Rayleigh scattering 133 ph->RegisterProcess(new G4RayleighScattering(), particle); 134 } 135 else if (particleName == "e-") { 136 // ionisation 137 G4eIonisation* eIoni = new G4eIonisation(); 138 139 // bremsstrahlung 140 G4eBremsstrahlung* eBrem = new G4eBremsstrahlung(); 141 142 ph->RegisterProcess(new G4eMultipleScattering(), particle); 143 ph->RegisterProcess(eIoni, particle); 144 ph->RegisterProcess(eBrem, particle); 145 } 146 else if (particleName == "e+") { 147 // ionisation 148 G4eIonisation* eIoni = new G4eIonisation(); 149 150 // bremsstrahlung 151 G4eBremsstrahlung* eBrem = new G4eBremsstrahlung(); 152 153 ph->RegisterProcess(new G4eMultipleScattering(), particle); 154 ph->RegisterProcess(eIoni, particle); 155 ph->RegisterProcess(eBrem, particle); 156 157 // annihilation at rest and in flight 158 ph->RegisterProcess(new G4eplusAnnihilation(), particle); 159 } 160 else if (particleName == "mu+" || particleName == "mu-") { 161 G4MuIonisation* muIoni = new G4MuIonisation(); 162 163 ph->RegisterProcess(muIoni, particle); 164 ph->RegisterProcess(mub, particle); 165 ph->RegisterProcess(mup, particle); 166 ph->RegisterProcess(new G4CoulombScattering(), particle); 167 } 168 else if (particleName == "alpha" || particleName == "He3") { 169 G4hMultipleScattering* msc = new G4hMultipleScattering(); 170 G4UrbanMscModel* model = new G4UrbanMscModel(); 171 model->SetActivationLowEnergyLimit(energyLimit); 172 msc->SetEmModel(model); 173 ph->RegisterProcess(msc, particle); 174 175 G4ionIonisation* ionIoni = new G4ionIonisation(); 176 ph->RegisterProcess(ionIoni, particle); 177 178 ph->RegisterProcess(nucr, particle); 179 } 180 else if (particleName == "GenericIon") { 181 G4hMultipleScattering* msc = new G4hMultipleScattering(); 182 G4UrbanMscModel* model = new G4UrbanMscModel(); 183 model->SetActivationLowEnergyLimit(energyLimit); 184 msc->SetEmModel(model); 185 ph->RegisterProcess(msc, particle); 186 187 G4ionIonisation* ionIoni = new G4ionIonisation(); 188 ionIoni->SetEmModel(new G4IonParametrisedLossModel()); 189 ph->RegisterProcess(ionIoni, particle); 190 191 ph->RegisterProcess(nucr, particle); 192 } 193 else if (particleName == "proton" || particleName == "deuteron" || particleName == "triton") { 194 G4hMultipleScattering* msc = new G4hMultipleScattering(); 195 G4UrbanMscModel* model = new G4UrbanMscModel(); 196 model->SetActivationLowEnergyLimit(energyLimit); 197 msc->SetEmModel(model); 198 ph->RegisterProcess(msc, particle); 199 200 G4hIonisation* hIoni = new G4hIonisation(); 201 ph->RegisterProcess(hIoni, particle); 202 ph->RegisterProcess(nucr, particle); 203 } 204 else if ((!particle->IsShortLived()) && (particle->GetPDGCharge() != 0.0) 205 && (particle->GetParticleName() != "chargedgeantino")) 206 { 207 // all others charged particles except geantino 208 209 ph->RegisterProcess(new G4hMultipleScattering(), particle); 210 ph->RegisterProcess(new G4hIonisation(), particle); 211 } 212 } 213 214 // Deexcitation 215 G4VAtomDeexcitation* de = new G4UAtomicDeexcitation(); 216 G4LossTableManager::Instance()->SetAtomDeexcitation(de); 217 } 218 219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 220