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 // 27 //--------------------------------------------------------------------------- 28 // 29 // ClassName: G4EmStandardPhysicsWVI 30 // 31 // Author: V.Ivanchenko 09.11.2005 32 // 33 // Modified: 34 // 35 // 21.04.2008 V.Ivanchenko add long-lived D and B mesons 36 // 37 //---------------------------------------------------------------------------- 38 // 39 40 #include "G4EmStandardPhysicsWVI.hh" 41 #include "G4SystemOfUnits.hh" 42 #include "G4ParticleDefinition.hh" 43 #include "G4EmParameters.hh" 44 #include "G4EmBuilder.hh" 45 #include "G4LossTableManager.hh" 46 47 #include "G4ComptonScattering.hh" 48 #include "G4GammaConversion.hh" 49 #include "G4PhotoElectricEffect.hh" 50 #include "G4RayleighScattering.hh" 51 52 #include "G4KleinNishinaModel.hh" 53 #include "G4LivermorePhotoElectricModel.hh" 54 #include "G4eMultipleScattering.hh" 55 #include "G4hMultipleScattering.hh" 56 #include "G4CoulombScattering.hh" 57 #include "G4WentzelVIModel.hh" 58 #include "G4WentzelVIRelModel.hh" 59 #include "G4UrbanMscModel.hh" 60 #include "G4hCoulombScatteringModel.hh" 61 #include "G4eCoulombScatteringModel.hh" 62 63 #include "G4eIonisation.hh" 64 #include "G4eBremsstrahlung.hh" 65 #include "G4eplusAnnihilation.hh" 66 #include "G4UAtomicDeexcitation.hh" 67 68 #include "G4hIonisation.hh" 69 #include "G4ionIonisation.hh" 70 #include "G4BetheHeitler5DModel.hh" 71 #include "G4AtimaEnergyLossModel.hh" 72 #include "G4AtimaFluctuations.hh" 73 #include "G4IonParametrisedLossModel.hh" 74 #include "G4LindhardSorensenIonModel.hh" 75 #include "G4BraggIonModel.hh" 76 #include "G4IonFluctuations.hh" 77 #include "G4NuclearStopping.hh" 78 #include "G4eplusTo2or3GammaModel.hh" 79 80 #include "G4Gamma.hh" 81 #include "G4Electron.hh" 82 #include "G4Positron.hh" 83 #include "G4GenericIon.hh" 84 85 #include "G4PhysicsListHelper.hh" 86 #include "G4BuilderType.hh" 87 #include "G4EmModelActivator.hh" 88 89 // factory 90 #include "G4PhysicsConstructorFactory.hh" 91 // 92 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmStandardPhysicsWVI); 93 94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 95 96 G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI(G4int ver) 97 : G4VPhysicsConstructor("G4EmStandardWVI") 98 { 99 G4EmParameters* param = G4EmParameters::Instance(); 100 param->SetDefaults(); 101 param->SetVerbose(ver); 102 param->SetMinEnergy(10*CLHEP::eV); 103 param->SetLowestElectronEnergy(100*CLHEP::eV); 104 param->SetNumberOfBinsPerDecade(20); 105 param->ActivateAngularGeneratorForIonisation(true); 106 param->SetStepFunction(0.2, 100*CLHEP::um); 107 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um); 108 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um); 109 param->SetStepFunctionIons(0.1, 1*CLHEP::um); 110 param->SetUseMottCorrection(true); 111 param->SetMuHadLateralDisplacement(true); 112 param->SetUseICRU90Data(true); 113 param->SetMscThetaLimit(0.15); 114 param->SetFluo(true); 115 param->SetMaxNIELEnergy(1*CLHEP::MeV); 116 SetPhysicsType(bElectromagnetic); 117 } 118 119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 120 121 G4EmStandardPhysicsWVI::~G4EmStandardPhysicsWVI() = default; 122 123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 124 125 void G4EmStandardPhysicsWVI::ConstructParticle() 126 { 127 // minimal set of particles for EM physics 128 G4EmBuilder::ConstructMinimalEmSet(); 129 } 130 131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 132 133 void G4EmStandardPhysicsWVI::ConstructProcess() 134 { 135 if(verboseLevel > 1) { 136 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl; 137 } 138 G4EmBuilder::PrepareEMPhysics(); 139 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 140 G4EmParameters* param = G4EmParameters::Instance(); 141 142 // common processes 143 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc"); 144 145 // nuclear stopping is enabled if th eenergy limit above zero 146 G4double nielEnergyLimit = param->MaxNIELEnergy(); 147 G4NuclearStopping* pnuc = nullptr; 148 if(nielEnergyLimit > 0.0) { 149 pnuc = new G4NuclearStopping(); 150 pnuc->SetMaxKinEnergy(nielEnergyLimit); 151 } 152 153 // high energy limit for e+- scattering models 154 G4double highEnergyLimit = 1*CLHEP::MeV; 155 156 // Add gamma EM processes 157 G4ParticleDefinition* particle = G4Gamma::Gamma(); 158 159 G4PhotoElectricEffect* pee = new G4PhotoElectricEffect(); 160 pee->SetEmModel(new G4LivermorePhotoElectricModel()); 161 162 G4ComptonScattering* cs = new G4ComptonScattering; 163 cs->SetEmModel(new G4KleinNishinaModel()); 164 165 G4GammaConversion* gc = new G4GammaConversion(); 166 if(param->EnablePolarisation()) { 167 gc->SetEmModel(new G4BetheHeitler5DModel()); 168 } 169 170 ph->RegisterProcess(pee, particle); 171 ph->RegisterProcess(cs, particle); 172 ph->RegisterProcess(gc, particle); 173 ph->RegisterProcess(new G4RayleighScattering(), particle); 174 175 // e- 176 particle = G4Electron::Electron(); 177 178 G4eMultipleScattering* msc = new G4eMultipleScattering; 179 G4UrbanMscModel* msc1 = new G4UrbanMscModel(); 180 G4WentzelVIModel* msc2 = new G4WentzelVIModel(); 181 msc1->SetHighEnergyLimit(highEnergyLimit); 182 msc2->SetLowEnergyLimit(highEnergyLimit); 183 msc->SetEmModel(msc1); 184 msc->SetEmModel(msc2); 185 186 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 187 G4CoulombScattering* ss = new G4CoulombScattering(); 188 ss->SetEmModel(ssm); 189 ss->SetMinKinEnergy(highEnergyLimit); 190 ssm->SetLowEnergyLimit(highEnergyLimit); 191 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 192 193 ph->RegisterProcess(msc, particle); 194 ph->RegisterProcess(new G4eIonisation(), particle); 195 ph->RegisterProcess(new G4eBremsstrahlung(), particle); 196 ph->RegisterProcess(ss, particle); 197 198 // e+ 199 particle = G4Positron::Positron(); 200 201 msc = new G4eMultipleScattering; 202 msc1 = new G4UrbanMscModel(); 203 msc2 = new G4WentzelVIModel(); 204 msc1->SetHighEnergyLimit(highEnergyLimit); 205 msc2->SetLowEnergyLimit(highEnergyLimit); 206 msc->SetEmModel(msc1); 207 msc->SetEmModel(msc2); 208 209 ssm = new G4eCoulombScatteringModel(); 210 ss = new G4CoulombScattering(); 211 ss->SetEmModel(ssm); 212 ss->SetMinKinEnergy(highEnergyLimit); 213 ssm->SetLowEnergyLimit(highEnergyLimit); 214 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 215 216 G4eplusAnnihilation* ann = new G4eplusAnnihilation(); 217 ann->SetEmModel(new G4eplusTo2or3GammaModel()); 218 219 ph->RegisterProcess(msc, particle); 220 ph->RegisterProcess(new G4eIonisation(), particle); 221 ph->RegisterProcess(new G4eBremsstrahlung(), particle); 222 ph->RegisterProcess(ann, particle); 223 ph->RegisterProcess(ss, particle); 224 225 // generic ion 226 particle = G4GenericIon::GenericIon(); 227 G4ionIonisation* ionIoni = new G4ionIonisation(); 228 auto fluc = new G4AtimaFluctuations(); 229 ionIoni->SetFluctModel(fluc); 230 ionIoni->SetEmModel(new G4AtimaEnergyLossModel()); 231 ph->RegisterProcess(hmsc, particle); 232 ph->RegisterProcess(ionIoni, particle); 233 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); } 234 235 // muons, hadrons, ions 236 G4EmBuilder::ConstructCharged(hmsc, pnuc); 237 238 // extra configuration 239 G4EmModelActivator mact(GetPhysicsName()); 240 } 241 242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 243